ITM Grid Service Library: Fortran 90
|
00001 program itm_grid_example8_separatrix 00002 00003 00004 use itm_types ! ITM standard Fortran data types, e.g. real(R8) 00005 use itm_constants 00006 use euITM_schemas ! ITM data structure/CPO definitions 00007 use euITM_routines ! ITM UAL routines 00008 use itm_grid_common ! ITM grid service library: constant definitions like COORDTYPE_R, COORDTYPE_Z 00009 use itm_grid_structured ! ITM grid service libary, routines for handling structured grids 00010 00011 implicit none 00012 00013 ! parameters 00014 integer, parameter :: NPOINTR = 6 00015 integer, parameter :: NPOINTZ = 5 00016 00017 integer, parameter :: SHOTNUM = 17151 00018 integer, parameter :: RUNNUM = 20407 00019 00020 ! variables 00021 type (type_edge) :: edgecpo 00022 integer :: idx 00023 integer :: ir, iz, i, j, k 00024 real(R8) :: cellData(NPOINTR - 1, NPOINTZ - 1) 00025 real(R8) :: nodeData(NPOINTR, NPOINTZ) 00026 00027 integer :: sepId,ompId,n_sep_objs,n_omp_objs,match 00028 type(GridObject) :: obj,face,node 00029 type(GridObject),dimension(:),allocatable :: nodes 00030 integer,dimension(:),allocatable :: sepNodes, ompNodes 00031 00032 real(R8) :: ne, Te, Ti 00033 00034 ! === 1. Set up CPO === 00035 00036 write (*,*) "Example 8: reading from shot ", SHOTNUM, ", run ", RUNNUM 00037 call euitm_open_env( 'euitm', SHOTNUM, RUNNUM, idx, 'coster', 'aug', '4.09a') 00038 call euitm_get_slice( idx, "edge", edgecpo, 2.5d0, UAL_PREVIOUS_SAMPLE ) 00039 00040 ! === 00041 00042 ! figure out subgrid id for separatrix 00043 sepId = gridFindSubgridByName(edgecpo%grid, "Separatrix") 00044 ompId = gridFindSubgridByName(edgecpo%grid, "Outer midplane") 00045 00046 !print *,sepId,ompId 00047 00048 n_sep_objs = gridSubGridSize(edgecpo%grid%subgrids(sepId)) 00049 n_omp_objs = gridSubGridSize(edgecpo%grid%subgrids(ompId)) 00050 00051 allocate(sepNodes(2*n_sep_objs),ompNodes(n_omp_objs)) 00052 sepNodes = -9999 00053 ompNodes = -9999 00054 00055 k = 0 00056 do i = 1,n_sep_objs 00057 face = subGridGetObject(edgecpo%grid%subgrids(sepId),i) 00058 !call gridWriteGridObject(face) 00059 call getComposingObjects(edgecpo%grid,face,nodes) 00060 !call writeObjectList(nodes) 00061 do j = 1,size(nodes) 00062 !call gridWriteGridObject(nodes(j)) 00063 !if (k .eq. 1) then 00064 ! sepNodes(k) = nodes(j)%ind 00065 ! k = k + 1 00066 !else 00067 if (all(sepNodes .ne. nodes(j)%ind(1),1)) then 00068 k = k + 1 00069 !print *,k,nodes(j)%ind(1) 00070 sepNodes(k) = nodes(j)%ind(1) 00071 end if 00072 !end if 00073 end do 00074 deallocate(nodes) 00075 end do 00076 00077 !print *, sepNodes(:k) 00078 00079 do i = 1,n_omp_objs 00080 node = subGridGetObject(edgecpo%grid%subgrids(ompId),i) 00081 !call gridWriteGridObject(node) 00082 if (any(sepNodes(:k) .eq. node%ind(1))) then 00083 match = i !node%ind(1) 00084 !print *, "Node ", match, " at position ", i, " is in outer midplane and separatrix" 00085 end if 00086 end do 00087 00088 do i = 1,size(edgecpo%fluid%ne%value) 00089 !print *,ompId, edgecpo%fluid%ne%value(i)%subgrid 00090 if (edgecpo%fluid%ne%value(i)%subgrid .eq. ompId) then 00091 ne = edgecpo%fluid%ne%value(i)%scalar(match) 00092 te = edgecpo%fluid%te%value(i)%scalar(match) 00093 ti = edgecpo%fluid%ti(1)%value(i)%scalar(match) 00094 print *, ne, te, ti 00095 else 00096 !print *, "Cannot find the outer midplane values in the edge cpo." 00097 end if 00098 end do 00099 00100 end program itm_grid_example8_separatrix