grid.F90 coverage: 63.33 %func 59.33 %block
1) module Grid_module
2)
3) use Grid_Structured_module
4) use Grid_Unstructured_module
5) use Grid_Unstructured_Aux_module
6) use Grid_Unstructured_Polyhedra_module
7) use Connection_module
8)
9) use PFLOTRAN_Constants_module
10)
11) implicit none
12)
13) private
14)
15) #include "petsc/finclude/petscsys.h"
16) #include "petsc/finclude/petscvec.h"
17) #include "petsc/finclude/petscvec.h90"
18) #include "petsc/finclude/petscis.h"
19) #include "petsc/finclude/petscis.h90"
20) #include "petsc/finclude/petscmat.h"
21) #include "petsc/finclude/petscmat.h90"
22)
23) type, public :: grid_type
24)
25) character(len=MAXWORDLENGTH) :: ctype
26) PetscInt :: itype ! type of grid (e.g. structured_grid, implicit_unstructured_grid, etc.)
27)
28) PetscInt :: nmax ! Total number of nodes in global domain
29) PetscInt :: nlmax ! Total number of non-ghosted nodes in local domain.
30) PetscInt :: ngmax ! Number of ghosted & non-ghosted nodes in local domain.
31) PetscInt :: global_offset ! Offset of first cell on process in petsc ordering
32) PetscInt :: nlmax_faces ! Total number of non-ghosted faces in local domain.
33) PetscInt :: ngmax_faces ! Number of ghosted & non-ghosted faces in local domain.
34) PetscInt :: nmax_faces ! Number of ghosted & non-ghosted faces in local domain.
35) PetscInt :: global_cell_offset, global_faces_offset ! offsets for LP formulation
36)
37) ! Below, we define several arrays used for mapping between different
38) ! types of array indices. Our terminology is as follows:
39) !
40) ! 'Local' indices are used to access arrays containing values that are
41) ! entirely local to the MPI process -- these arrays contain no "ghost"
42) ! entries used to hold copies of values that are owned by neighboring
43) ! processes.
44) !
45) ! 'Ghosted local' (or simply 'ghost') indices are used to access arrays
46) ! that contain additional entries that hold copies of values that are
47) ! owned by neighboring processes. (These entries are filled in by
48) ! DMGlobalToLocalBegin/End() in the structured grid case.)
49) !
50) ! Entries of a vector created with DMCreateGlobalVector() should be
51) ! indexed using 'local' indices. The array returned from a call to
52) ! VecGetArrayF90() on such a vector consists of local entries only and
53) ! NO ghost points.
54) !
55) ! Entries of a vector created with DMCreateLocalVector() should be
56) ! indexed using 'ghosted local' indices. The array returned from a call
57) ! to VecGetArrayF90() on such a vector contains the truly3 local entries
58) ! as well as ghost points.
59) !
60) ! The index mapping arrays are the following:
61) ! nL2G : not collective, local processor: local => ghosted local
62) ! nG2L : not collective, local processor: ghosted local => local
63) ! nG2A : not collective, ghosted local => natural
64)
65) PetscInt, pointer :: nL2G(:), nG2L(:)
66) PetscInt, pointer :: nG2A(:)
67)
68) PetscReal, pointer :: x(:), y(:), z(:) ! coordinates of ghosted grid cells
69)
70) PetscReal :: x_min_global, y_min_global, z_min_global
71) PetscReal :: x_max_global, y_max_global, z_max_global
72) PetscReal :: x_min_local, y_min_local, z_min_local
73) PetscReal :: x_max_local, y_max_local, z_max_local
74)
75) PetscInt, pointer :: hash(:,:,:)
76) PetscInt :: num_hash_bins
77)
78) type(grid_structured_type), pointer :: structured_grid
79) type(grid_unstructured_type), pointer :: unstructured_grid
80)
81) type(connection_set_list_type), pointer :: internal_connection_set_list
82) type(connection_set_list_type), pointer :: boundary_connection_set_list
83)
84) end type grid_type
85)
86) type, public :: face_type
87) type(connection_set_type), pointer :: conn_set_ptr
88) PetscInt :: id
89) end type face_type
90)
91) public :: GridCreate, &
92) GridDestroy, &
93) GridComputeInternalConnect, &
94) GridMapIndices, &
95) GridComputeSpacing, &
96) GridComputeCoordinates, &
97) GridComputeVolumes, &
98) GridComputeAreas, &
99) GridLocalizeRegions, &
100) GridLocalizeRegionsFromCellIDsUGrid, &
101) GridPopulateConnection, &
102) GridCopyIntegerArrayToVec, &
103) GridCopyRealArrayToVec, &
104) GridCopyVecToIntegerArray, &
105) GridCopyVecToRealArray, &
106) GridCreateNaturalToGhostedHash, &
107) GridDestroyHashTable, &
108) GridGetLocalGhostedIdFromHash, &
109) GridIndexToCellID, &
110) GridGetGhostedNeighbors, &
111) GridGetGhostedNeighborsWithCorners, &
112) GridMapCellsInPolVol, &
113) GridGetLocalIDFromCoordinate
114)
115) contains
116)
117) ! ************************************************************************** !
118)
119) function GridCreate()
120) !
121) ! Creates a structured or unstructured grid
122) !
123) ! Author: Glenn Hammond
124) ! Date: 10/23/07
125) !
126)
127) implicit none
128)
129) type(grid_type), pointer :: GridCreate
130)
131) type(grid_type), pointer :: grid
132)
133) allocate(grid)
134) grid%ctype = ''
135) grid%itype = 0
136)
137) nullify(grid%structured_grid)
138) nullify(grid%unstructured_grid)
139)
140) nullify(grid%internal_connection_set_list)
141)
142) nullify(grid%nL2G)
143) nullify(grid%nG2L)
144) nullify(grid%nG2A)
145)
146) nullify(grid%x)
147) nullify(grid%y)
148) nullify(grid%z)
149)
150) grid%x_min_global = 1.d20
151) grid%x_max_global = -1.d20
152) grid%y_min_global = 1.d20
153) grid%y_max_global = -1.d20
154) grid%z_min_global = 1.d20
155) grid%z_max_global = -1.d20
156)
157) grid%x_min_local = 1.d20
158) grid%x_max_local = -1.d20
159) grid%y_min_local = 1.d20
160) grid%y_max_local = -1.d20
161) grid%z_min_local = 1.d20
162) grid%z_max_local = -1.d20
163)
164) grid%nmax = 0
165) grid%nlmax = 0
166) grid%ngmax = 0
167) grid%global_offset = 0
168)
169) nullify(grid%hash)
170) grid%num_hash_bins = 1000
171)
172) GridCreate => grid
173)
174) end function GridCreate
175)
176) ! ************************************************************************** !
177)
178) subroutine GridComputeInternalConnect(grid,option,ugdm)
179) !
180) ! computes internal connectivity of a grid
181) ! sp modified December 2010
182) !
183) ! Author: Glenn Hammond
184) ! Date: 10/17/07
185) !
186)
187) use Connection_module
188) use Option_module
189) use Grid_Unstructured_Explicit_module
190) use Grid_Unstructured_Polyhedra_module
191)
192) implicit none
193)
194) PetscInt ierr
195)
196) type(grid_type) :: grid
197) type(option_type) :: option
198) type(ugdm_type), optional :: ugdm
199)
200) type(connection_set_type), pointer :: connection_set, connection_bound_set
201) type(connection_set_type), pointer :: connection_set_2
202) nullify(connection_set); nullify(connection_bound_set)
203)
204) select case(grid%itype)
205) case(STRUCTURED_GRID)
206) connection_set => &
207) StructGridComputeInternConnect( grid%structured_grid, grid%x, grid%y, &
208) grid%z, option)
209) case(IMPLICIT_UNSTRUCTURED_GRID)
210) connection_set => &
211) UGridComputeInternConnect(grid%unstructured_grid,grid%x,grid%y, &
212) grid%z,option)
213) case(EXPLICIT_UNSTRUCTURED_GRID)
214) connection_set => &
215) UGridExplicitSetInternConnect(grid%unstructured_grid%explicit_grid, &
216) option)
217) case(POLYHEDRA_UNSTRUCTURED_GRID)
218) connection_set => &
219) UGridPolyhedraComputeInternConnect(grid%unstructured_grid, &
220) grid%x, grid%y, grid%z, &
221) option)
222) call UGridPolyhedraComputeOutputInfo(grid%unstructured_grid, grid%nL2G, &
223) grid%nG2L, grid%nG2A, option)
224) end select
225)
226) allocate(grid%internal_connection_set_list)
227) call ConnectionInitList(grid%internal_connection_set_list)
228) call ConnectionAddToList(connection_set,grid%internal_connection_set_list)
229)
230)
231) select case(grid%itype)
232) case(IMPLICIT_UNSTRUCTURED_GRID)
233) ! connection_bound_set => &
234) ! UGridComputeBoundConnect(grid%unstructured_grid,option)
235) end select
236)
237) end subroutine GridComputeInternalConnect
238)
239) ! ************************************************************************** !
240)
241) subroutine GridPopulateConnection(grid,connection,iface,iconn,cell_id_local, &
242) option)
243) !
244) ! computes connectivity coupler to a grid
245) !
246) ! Author: Glenn Hammond
247) ! Date: 11/09/07
248) !
249)
250) use Connection_module
251) use Grid_Structured_module
252) use Option_module
253)
254) implicit none
255)
256) type(grid_type) :: grid
257) type(connection_set_type) :: connection
258) PetscInt :: iface
259) PetscInt :: iconn
260) PetscInt :: cell_id_local
261) type(option_type) :: option
262)
263) PetscInt :: cell_id_ghosted
264)
265) cell_id_ghosted = grid%nL2G(cell_id_local)
266) ! Use ghosted index to access dx, dy, dz because we have
267) ! already done a global-to-local scatter for computing the
268) ! interior node connections.
269)
270) select case(grid%itype)
271) case(STRUCTURED_GRID)
272) call StructGridPopulateConnection(grid%x,grid%structured_grid,connection, &
273) iface,iconn,cell_id_ghosted,option)
274) case(IMPLICIT_UNSTRUCTURED_GRID)
275) call UGridPopulateConnection(grid%unstructured_grid,connection,iface,&
276) iconn,cell_id_ghosted,option)
277) case(POLYHEDRA_UNSTRUCTURED_GRID)
278) call UGridPolyhedraPopulateConnection(grid%unstructured_grid,connection,iface, &
279) iconn,cell_id_ghosted,option)
280) end select
281)
282) end subroutine GridPopulateConnection
283)
284) ! ************************************************************************** !
285)
286) subroutine GridMapIndices(grid, dm_ptr, sgrid_stencil_type,option)
287) !
288) ! maps global, local and natural indices of cells
289) ! to each other
290) !
291) ! Author: Glenn Hammond
292) ! Date: 10/24/07
293) !
294)
295) use Option_module
296) use DM_Kludge_module
297)
298) implicit none
299)
300) type(grid_type) :: grid
301) type(dm_ptr_type) :: dm_ptr
302) PetscEnum :: sgrid_stencil_type
303) type(option_type) :: option
304)
305) PetscInt, allocatable :: int_tmp(:)
306) ! PetscInt, pointer :: int_tmp(:)
307) PetscInt :: n
308) PetscOffset :: i_da
309)
310) select case(grid%itype)
311) case(STRUCTURED_GRID)
312) call StructGridMapIndices(grid%structured_grid,sgrid_stencil_type, &
313) grid%nG2L,grid%nL2G,grid%nG2A, &
314) option)
315) case(IMPLICIT_UNSTRUCTURED_GRID,EXPLICIT_UNSTRUCTURED_GRID, &
316) POLYHEDRA_UNSTRUCTURED_GRID)
317) call UGridMapIndices(grid%unstructured_grid, &
318) dm_ptr%ugdm, &
319) grid%nG2L,grid%nL2G,grid%nG2A,option)
320) end select
321)
322)
323) end subroutine GridMapIndices
324)
325) ! ************************************************************************** !
326)
327) subroutine GridComputeSpacing(grid,origin_global,option)
328) !
329) ! Computes grid spacing (only for structured grid
330) !
331) ! Author: Glenn Hammond
332) ! Date: 10/26/07
333) !
334) use Option_module
335)
336) implicit none
337)
338) type(grid_type) :: grid
339) PetscReal :: origin_global(3)
340) type(option_type) :: option
341)
342) select case(grid%itype)
343) case(STRUCTURED_GRID)
344) call StructGridComputeSpacing(grid%structured_grid,origin_global,option)
345) case(IMPLICIT_UNSTRUCTURED_GRID)
346) end select
347)
348) end subroutine GridComputeSpacing
349)
350) ! ************************************************************************** !
351)
352) subroutine GridComputeCoordinates(grid,origin_global,option,ugdm)
353) !
354) ! Computes x,y,z coordinates of grid cells
355) !
356) ! Author: Glenn Hammond
357) ! Date: 10/24/07
358) !
359)
360) use Option_module
361) use Grid_Unstructured_Explicit_module
362) use Grid_Unstructured_Polyhedra_module
363)
364) implicit none
365)
366) type(grid_type) :: grid
367) PetscReal :: origin_global(3)
368) type(option_type) :: option
369) type(ugdm_type), optional :: ugdm ! sp
370) PetscInt :: icell
371)
372) PetscErrorCode :: ierr
373)
374) allocate(grid%x(grid%ngmax))
375) grid%x = 0.d0
376) allocate(grid%y(grid%ngmax))
377) grid%y = 0.d0
378) allocate(grid%z(grid%ngmax))
379) grid%z = 0.d0
380)
381) select case(grid%itype)
382) case(STRUCTURED_GRID)
383) call StructGridComputeCoord(grid%structured_grid,option, &
384) origin_global, &
385) grid%x,grid%y,grid%z, &
386) grid%x_min_local,grid%x_max_local, &
387) grid%y_min_local,grid%y_max_local, &
388) grid%z_min_local,grid%z_max_local)
389) case(IMPLICIT_UNSTRUCTURED_GRID)
390) call UGridComputeCoord(grid%unstructured_grid,option, &
391) grid%x,grid%y,grid%z, &
392) grid%x_min_local,grid%x_max_local, &
393) grid%y_min_local,grid%y_max_local, &
394) grid%z_min_local,grid%z_max_local)
395) case(EXPLICIT_UNSTRUCTURED_GRID)
396) call UGridExplicitSetCellCentroids(grid%unstructured_grid% &
397) explicit_grid, &
398) grid%x,grid%y,grid%z, &
399) grid%x_min_local,grid%x_max_local, &
400) grid%y_min_local,grid%y_max_local, &
401) grid%z_min_local,grid%z_max_local)
402) case(POLYHEDRA_UNSTRUCTURED_GRID)
403) call UGridPolyhedraSetCellCentroids(grid%unstructured_grid%polyhedra_grid, &
404) grid%x,grid%y,grid%z, &
405) grid%x_min_local,grid%x_max_local, &
406) grid%y_min_local,grid%y_max_local, &
407) grid%z_min_local,grid%z_max_local,option)
408)
409) end select
410)
411) if (associated(grid%structured_grid)) then
412) ! compute global max/min from the local max/in
413) call MPI_Allreduce(grid%x_min_local,grid%x_min_global,ONE_INTEGER_MPI, &
414) MPI_DOUBLE_PRECISION,MPI_MIN,option%mycomm,ierr)
415) call MPI_Allreduce(grid%y_min_local,grid%y_min_global,ONE_INTEGER_MPI, &
416) MPI_DOUBLE_PRECISION,MPI_MIN,option%mycomm,ierr)
417) call MPI_Allreduce(grid%z_min_local,grid%z_min_global,ONE_INTEGER_MPI, &
418) MPI_DOUBLE_PRECISION,MPI_MIN,option%mycomm,ierr)
419) call MPI_Allreduce(grid%x_max_local,grid%x_max_global,ONE_INTEGER_MPI, &
420) MPI_DOUBLE_PRECISION,MPI_MAX,option%mycomm,ierr)
421) call MPI_Allreduce(grid%y_max_local,grid%y_max_global,ONE_INTEGER_MPI, &
422) MPI_DOUBLE_PRECISION,MPI_MAX,option%mycomm,ierr)
423) call MPI_Allreduce(grid%z_max_local,grid%z_max_global,ONE_INTEGER_MPI, &
424) MPI_DOUBLE_PRECISION,MPI_MAX,option%mycomm,ierr)
425) endif
426)
427) if (associated(grid%unstructured_grid)) then
428) ! compute global max/min from the local max/in
429) call MPI_Allreduce(grid%x_min_local,grid%x_min_global,ONE_INTEGER_MPI, &
430) MPI_DOUBLE_PRECISION,MPI_MIN,option%mycomm,ierr)
431) call MPI_Allreduce(grid%y_min_local,grid%y_min_global,ONE_INTEGER_MPI, &
432) MPI_DOUBLE_PRECISION,MPI_MIN,option%mycomm,ierr)
433) call MPI_Allreduce(grid%z_min_local,grid%z_min_global,ONE_INTEGER_MPI, &
434) MPI_DOUBLE_PRECISION,MPI_MIN,option%mycomm,ierr)
435) call MPI_Allreduce(grid%x_max_local,grid%x_max_global,ONE_INTEGER_MPI, &
436) MPI_DOUBLE_PRECISION,MPI_MAX,option%mycomm,ierr)
437) call MPI_Allreduce(grid%y_max_local,grid%y_max_global,ONE_INTEGER_MPI, &
438) MPI_DOUBLE_PRECISION,MPI_MAX,option%mycomm,ierr)
439) call MPI_Allreduce(grid%z_max_local,grid%z_max_global,ONE_INTEGER_MPI, &
440) MPI_DOUBLE_PRECISION,MPI_MAX,option%mycomm,ierr)
441) !endif
442) endif
443)
444) end subroutine GridComputeCoordinates
445)
446) ! ************************************************************************** !
447)
448) subroutine GridComputeVolumes(grid,volume,option)
449) !
450) ! Computes the volumes of cells in structured grid
451) !
452) ! Author: Glenn Hammond
453) ! Date: 10/25/07
454) !
455)
456) use Option_module
457) use Grid_Unstructured_Explicit_module
458) use Grid_Unstructured_Polyhedra_module
459)
460) implicit none
461)
462) #include "petsc/finclude/petscvec.h"
463) #include "petsc/finclude/petscvec.h90"
464)
465) type(grid_type) :: grid
466) type(option_type) :: option
467) Vec :: volume
468)
469) select case(grid%itype)
470) case(STRUCTURED_GRID)
471) call StructGridComputeVolumes(grid%x,grid%structured_grid,option, &
472) grid%nL2G,volume)
473) case(IMPLICIT_UNSTRUCTURED_GRID)
474) call UGridComputeVolumes(grid%unstructured_grid,option,volume)
475) call UGridComputeQuality(grid%unstructured_grid,option)
476) case(EXPLICIT_UNSTRUCTURED_GRID)
477) call UGridExplicitComputeVolumes(grid%unstructured_grid, &
478) option,volume)
479) case(POLYHEDRA_UNSTRUCTURED_GRID)
480) call UGridPolyhedraComputeVolumes(grid%unstructured_grid,option,volume)
481) end select
482)
483) end subroutine GridComputeVolumes
484)
485) ! ************************************************************************** !
486)
487) subroutine GridComputeAreas(grid,area,option)
488) !
489) ! Computes the areas for 2D-mesh
490) !
491) ! Author: Gautam Bisht
492) ! Date: 03/07/2012
493) !
494)
495) use Option_module
496)
497) implicit none
498)
499) #include "petsc/finclude/petscvec.h"
500) #include "petsc/finclude/petscvec.h90"
501)
502) type(grid_type) :: grid
503) type(option_type) :: option
504) Vec :: area
505)
506) select case(grid%itype)
507) case(IMPLICIT_UNSTRUCTURED_GRID)
508) call UGridComputeAreas(grid%unstructured_grid,option,area)
509) call UGridComputeQuality(grid%unstructured_grid,option)
510) case default
511) option%io_buffer = 'ERROR: GridComputeAreas only implemented for Unstructured grid'
512) call printErrMsg(option)
513) end select
514)
515) end subroutine GridComputeAreas
516)
517) ! ************************************************************************** !
518)
519) subroutine GridLocalizeRegions(grid,region_list,option)
520) !
521) ! Resticts regions to cells local to processor
522) !
523) ! Author: Glenn Hammond
524) ! Date: 10/29/07
525) !
526)
527) use Option_module
528) use Region_module
529)
530) implicit none
531)
532) type(region_list_type), pointer :: region_list
533) type(grid_type), pointer :: grid
534) type(option_type) :: option
535)
536) type(region_type), pointer :: region
537) character(len=MAXSTRINGLENGTH) :: string
538) PetscInt, allocatable :: temp_int_array(:)
539) PetscInt :: i, j, k, count, local_count, ghosted_id, local_id
540) PetscInt :: i_min, i_max, j_min, j_max, k_min, k_max
541) PetscReal :: x_min, x_max, y_min, y_max, z_min, z_max
542) PetscReal, parameter :: pert = 1.d-8, tol = 1.d-20
543) PetscReal :: x_shift, y_shift, z_shift
544) PetscReal :: del_x, del_y, del_z
545) PetscInt :: iflag, global_cell_count
546) PetscBool :: same_point
547) PetscErrorCode :: ierr
548)
549) iflag = 0
550) region => region_list%first
551) do
552) if (.not.associated(region)) exit
553)
554) select case(region%def_type)
555) case (DEFINED_BY_BLOCK)
556) call GridLocalizeRegionFromBlock(grid,region,option)
557) case (DEFINED_BY_CARTESIAN_BOUNDARY)
558) call GridLocalizeRegionFromCartBound(grid,region,option)
559) case (DEFINED_BY_COORD)
560) call GridLocalizeRegionFromCoordinates(grid,region,option)
561) case (DEFINED_BY_CELL_IDS)
562) select case(grid%itype)
563) ! case(STRUCTURED_GRID)
564) ! The region is localized in InitCommonReadRegionFiles->
565) ! HDF5ReadRegionFromFile->HDF5MapLocalToNaturalIndices
566) case(IMPLICIT_UNSTRUCTURED_GRID)
567) if (region%hdf5_ugrid_kludge) then
568) call GridLocalizeRegionsFromCellIDsUGrid(grid,region,option)
569) endif
570) case(EXPLICIT_UNSTRUCTURED_GRID)
571) call GridLocalizeRegionsFromCellIDsUGrid(grid,region,option)
572) ! case(STRUCTURED_GRID)
573) ! The region is localized in
574) end select
575) case (DEFINED_BY_CELL_AND_FACE_IDS)
576) select case(grid%itype)
577) case (STRUCTURED_GRID)
578) ! Do nothing since the region was localized during the reading process
579) case default
580) option%io_buffer = 'GridLocalizeRegions() must tbe extended ' // &
581) 'for unstructured region DEFINED_BY_CELL_AND_FACE_IDS'
582) call printErrMsg(option)
583) end select
584) case (DEFINED_BY_VERTEX_IDS)
585) option%io_buffer = 'GridLocalizeRegions() must tbe extended ' // &
586) 'for unstructured region DEFINED_BY_VERTEX_IDS'
587) call printErrMsg(option)
588) case (DEFINED_BY_SIDESET_UGRID)
589) call UGridMapSideSet(grid%unstructured_grid, &
590) region%sideset%face_vertices, &
591) region%sideset%nfaces,region%name, &
592) option,region%cell_ids,region%faces)
593) region%num_cells = size(region%cell_ids)
594) case (DEFINED_BY_FACE_UGRID_EXP)
595) call GridLocalizeExplicitFaceset(grid%unstructured_grid,region, &
596) option)
597) case (DEFINED_BY_POLY_BOUNDARY_FACE)
598) call UGridMapBoundFacesInPolVol(grid%unstructured_grid, &
599) region%polygonal_volume, &
600) region%name,option, &
601) region%cell_ids,region%faces)
602) region%num_cells = size(region%cell_ids)
603) case (DEFINED_BY_POLY_CELL_CENTER)
604) call GridMapCellsInPolVol(grid, &
605) region%polygonal_volume, &
606) region%name,option, &
607) region%cell_ids)
608) region%num_cells = size(region%cell_ids)
609) case default
610) option%io_buffer = 'GridLocalizeRegions: Region definition not recognized'
611) call printErrMsg(option)
612) end select
613)
614) if (region%num_cells == 0 .and. associated(region%cell_ids)) then
615) deallocate(region%cell_ids)
616) nullify(region%cell_ids)
617) endif
618)
619) if (region%num_cells == 0 .and. associated(region%faces)) then
620) deallocate(region%faces)
621) nullify(region%faces)
622) endif
623)
624) ! check to ensure that there is at least one grid cell in each region
625) call MPI_Allreduce(region%num_cells,global_cell_count,ONE_INTEGER_MPI, &
626) MPI_INTEGER, MPI_SUM, option%mycomm,ierr)
627) if (global_cell_count == 0) then
628) option%io_buffer = 'No cells assigned to REGION "' // &
629) trim(region%name) // '".'
630) call printErrMsg(option)
631) endif
632) region => region%next
633)
634) enddo
635)
636) #if 0
637) !
638) ! GB: Older formulation. Need to remove it.
639) !
640) iflag = 0
641) region => region_list%first
642) do
643)
644) if (.not.associated(region)) exit
645)
646) if (.not.(associated(region%cell_ids) .or. &
647) associated(region%sideset) .or. &
648) associated(region%polygonal_volume) .or. &
649) associated(region%vertex_ids))) then
650) ! i, j, k block
651) if (region%i1 > 0 .and. region%i2 > 0 .and. &
652) region%j1 > 0 .and. region%j2 > 0 .and. &
653) region%k1 > 0 .and. region%k2 > 0) then
654)
655) call GridLocalizeRegionFromBlock(grid,region,option)
656)
657) else if (associated(region%coordinates)) then
658) call GridLocalizeRegionFromCoordinates(grid,region,option)
659) endif
660) else if (associated(region%sideset)) then
661) call UGridMapSideSet(grid%unstructured_grid, &
662) region%sideset%face_vertices, &
663) region%sideset%nfaces,region%name, &
664) option,region%cell_ids,region%faces)
665) region%num_cells = size(region%cell_ids)
666) else if (associated(region%polygonal_volume)) then
667) select case(region%def_type)
668) case(DEFINED_BY_POLY_BOUND_UGRID)
669) call UGridMapBoundFacesInPolVol(grid%unstructured_grid, &
670) region%polygonal_volume, &
671) region%name,option, &
672) region%cell_ids,region%faces)
673) case(DEFINED_BY_POLY_VOLUME)
674) call GridMapCellsInPolVol(grid,region%polygonal_volume, &
675) region%name,option,region%cell_ids)
676) end select
677) region%num_cells = size(region%cell_ids)
678) else if (associated(region%cell_ids)) then
679) select case(grid%itype)
680) case(IMPLICIT_UNSTRUCTURED_GRID)
681) call GridLocalizeRegionsFromCellIDsUGrid(grid,region,option)
682)
683) case(STRUCTURED_GRID)
684) !sp following was commented out
685) !sp remove?
686) !geh: Not for now. The below maps a list of natural ids to local. This is now done
687) ! in hdf5.F90 for lists of cells from hdf5 files. If we elect to support ascii
688) ! files, we will need this functionality here.
689) #if 0
690) allocate(temp_int_array(region%num_cells))
691) do count=1,region%num_cells
692) i = mod(region%cell_ids(count),grid%structured_grid%nx) - &
693) grid%structured_grid%lxs
694) j = mod((region%cell_ids(count)-1)/grid%structured_grid%nx, &
695) grid%structured_grid%ny)+1 - &
696) grid%structured_grid%lys
697) k = ((region%cell_ids(count)-1)/grid%structured_grid%nxy)+1 - &
698) grid%structured_grid%lzs
699) if (i > 0 .and. i <= grid%structured_grid%nlx .and. &
700) j > 0 .and. j <= grid%structured_grid%nly .and. &
701) k > 0 .and. k <= grid%structured_grid%nlz) then
702) temp_int_array(local_count) = &
703) i + (j-1)*grid%structured_grid%nlx + &
704) (k-1)*grid%structured_grid%nlxy
705) local_count = local_count + 1
706) endif
707) enddo
708) #endif
709) case(EXPLICIT_UNSTRUCTURED_GRID)
710) call GridLocalizeExplicitFaceset(grid%unstructured_grid,region, &
711) option)
712) case default
713) option%io_buffer = 'GridLocalizeRegions: define region by list ' // &
714) 'of cells not implemented: ' // trim(region%name)
715) call printErrMsg(option)
716) end select
717) !sp end
718) endif
719)
720) if (region%num_cells == 0 .and. associated(region%cell_ids)) then
721) deallocate(region%cell_ids)
722) nullify(region%cell_ids)
723) endif
724) if (region%num_cells == 0 .and. associated(region%faces)) then
725) deallocate(region%faces)
726) nullify(region%faces)
727) endif
728) region => region%next
729)
730) enddo
731) #endif
732)
733) ! Assign
734)
735) end subroutine GridLocalizeRegions
736)
737) ! ************************************************************************** !
738)
739) subroutine GridLocalizeRegionsFromCellIDsUGrid(grid, region, option)
740) !
741) ! Resticts regions to cells local
742) ! to processor for unstrucutred mesh when the region is defined by a
743) ! list of cell ids (in natural order)
744) !
745) ! Author: Gautam Bisht
746) ! Date: 5/30/2011
747) !
748)
749) use Option_module
750) use Region_module
751)
752) implicit none
753)
754) #include "petsc/finclude/petsclog.h"
755) #include "petsc/finclude/petscviewer.h"
756) #include "petsc/finclude/petscvec.h"
757) #include "petsc/finclude/petscvec.h90"
758) #include "petsc/finclude/petscis.h"
759) #include "petsc/finclude/petscis.h90"
760) #include "petsc/finclude/petscmat.h"
761)
762) type(grid_type) :: grid
763) type(region_type) :: region
764) type(option_type) :: option
765)
766) ! local
767) type(grid_unstructured_type),pointer :: ugrid
768) Vec :: vec_cell_ids,vec_cell_ids_loc
769) Vec :: vec_face_ids,vec_face_ids_loc
770) IS :: is_from, is_to
771) Mat :: adj, adj_t,adj_d, adj_o
772) VecScatter :: vec_scat
773) PetscErrorCode :: ierr
774) PetscViewer :: viewer
775) PetscInt :: ii,jj,kk,count
776) PetscInt :: istart, iend
777) PetscInt :: ghosted_id,local_id,natural_id
778) PetscInt,pointer :: tmp_int_array(:),tmp_int_array2(:)
779) PetscScalar,pointer :: v_loc_p(:),v_loc2_p(:)
780) PetscScalar,pointer :: tmp_scl_array(:)
781)
782)
783) PetscInt, pointer :: ia_p(:), ja_p(:)
784) PetscInt :: n,rstart,rend,icol(1)
785) PetscInt :: index
786) PetscInt :: vertex_id
787) PetscOffset :: iia,jja,aaa,iicol
788) PetscBool :: done,found
789) PetscScalar :: aa(1)
790) PetscInt :: cell_id_max_local
791) PetscInt :: cell_id_max_global
792) ! PetscScalar, pointer :: aa(:)
793) ! Would like to use the above, but I have to fix MatGetArrayF90() first. --RTM
794)
795) ugrid => grid%unstructured_grid
796)
797) if (associated(region%cell_ids)) then
798)
799) call VecCreateMPI(option%mycomm, ugrid%nlmax, PETSC_DECIDE, &
800) vec_cell_ids, ierr);CHKERRQ(ierr)
801) call VecCreateMPI(option%mycomm, ugrid%nlmax, PETSC_DECIDE, &
802) vec_cell_ids_loc, ierr);CHKERRQ(ierr)
803)
804) call VecZeroEntries(vec_cell_ids, ierr);CHKERRQ(ierr)
805)
806) allocate(tmp_int_array(region%num_cells))
807) allocate(tmp_scl_array(region%num_cells))
808)
809) count = 0
810) cell_id_max_local = -1
811) do ii = 1, region%num_cells
812) count = count + 1
813) tmp_int_array(count) = region%cell_ids(ii) - 1
814) tmp_scl_array(count) = 1.d0
815) cell_id_max_local = max(cell_id_max_local, region%cell_ids(ii))
816) enddo
817)
818) call MPI_Allreduce(cell_id_max_local, cell_id_max_global, ONE_INTEGER_MPI, &
819) MPI_INTEGER, MPI_MAX, option%mycomm,ierr)
820) if (cell_id_max_global > grid%nmax) then
821) option%io_buffer = 'The following region includes a cell-id that is greater than ' // &
822) 'number of control volumes present in the grid: ' // trim(region%name)
823) call printErrMsg(option)
824) endif
825)
826) call VecSetValues(vec_cell_ids, region%num_cells, tmp_int_array, &
827) tmp_scl_array, ADD_VALUES, ierr);CHKERRQ(ierr)
828)
829) deallocate(tmp_int_array)
830) deallocate(tmp_scl_array)
831)
832) call VecAssemblyBegin(vec_cell_ids, ierr);CHKERRQ(ierr)
833) call VecAssemblyEnd(vec_cell_ids, ierr);CHKERRQ(ierr)
834)
835) allocate(tmp_int_array(ugrid%nlmax))
836) count = 0
837) do ghosted_id = 1, ugrid%ngmax
838) local_id = grid%nG2L(ghosted_id)
839) if (local_id < 1) cycle
840) count = count + 1
841) natural_id = grid%nG2A(ghosted_id)
842) tmp_int_array(count) = natural_id
843) enddo
844)
845) tmp_int_array = tmp_int_array - 1
846) call ISCreateBlock(option%mycomm, 1, ugrid%nlmax, &
847) tmp_int_array, PETSC_COPY_VALUES, is_from, &
848) ierr);CHKERRQ(ierr)
849)
850) call VecGetOwnershipRange(vec_cell_ids_loc,istart,iend,ierr);CHKERRQ(ierr)
851) do ii=1,ugrid%nlmax
852) tmp_int_array(ii) = ii + istart
853) enddo
854)
855) tmp_int_array = tmp_int_array - 1
856) call ISCreateBlock(option%mycomm, 1, ugrid%nlmax, &
857) tmp_int_array, PETSC_COPY_VALUES, is_to, &
858) ierr);CHKERRQ(ierr)
859) deallocate(tmp_int_array)
860)
861) call VecScatterCreate(vec_cell_ids,is_from,vec_cell_ids_loc,is_to, &
862) vec_scat, ierr);CHKERRQ(ierr)
863) call ISDestroy(is_from, ierr);CHKERRQ(ierr)
864) call ISDestroy(is_to, ierr);CHKERRQ(ierr)
865)
866) call VecScatterBegin(vec_scat, vec_cell_ids, vec_cell_ids_loc, &
867) INSERT_VALUES, SCATTER_FORWARD, ierr);CHKERRQ(ierr)
868) call VecScatterEnd(vec_scat, vec_cell_ids, vec_cell_ids_loc, &
869) INSERT_VALUES, SCATTER_FORWARD, ierr);CHKERRQ(ierr)
870) call VecScatterDestroy(vec_scat, ierr);CHKERRQ(ierr)
871)
872) call VecGetArrayF90(vec_cell_ids_loc, v_loc_p, ierr);CHKERRQ(ierr)
873) count = 0
874) do ii=1, ugrid%nlmax
875) if (v_loc_p(ii) == 1) count = count + 1
876) enddo
877)
878) region%num_cells = count
879) if (count > 0) then
880) allocate(tmp_int_array(count))
881) count = 0
882) do ii =1, ugrid%nlmax
883) if (v_loc_p(ii) == 1) then
884) count = count + 1
885) tmp_int_array(count) = ii
886) endif
887) enddo
888)
889) deallocate(region%cell_ids)
890) allocate(region%cell_ids(region%num_cells))
891) region%cell_ids = tmp_int_array
892) deallocate(tmp_int_array)
893) else
894) deallocate(region%cell_ids)
895) allocate(region%cell_ids(region%num_cells))
896) endif
897)
898) call VecRestoreArrayF90(vec_cell_ids_loc,v_loc_p,ierr);CHKERRQ(ierr)
899)
900) call VecDestroy(vec_cell_ids,ierr);CHKERRQ(ierr)
901) call VecDestroy(vec_cell_ids_loc,ierr);CHKERRQ(ierr)
902)
903) endif
904)
905) end subroutine GridLocalizeRegionsFromCellIDsUGrid
906)
907) ! ************************************************************************** !
908)
909) subroutine GridLocalizeExplicitFaceset(ugrid,region,option)
910) !
911) ! GridLocalizeExplicitFaceset
912) !
913) ! Author: Glenn Hammond
914) ! Date: 10/10/12
915) !
916)
917) use Region_module
918) use Option_module
919)
920) implicit none
921)
922) type(grid_unstructured_type) :: ugrid
923) type(region_type) :: region
924) type(option_type) :: option
925) Vec :: volume
926)
927) type(unstructured_explicit_type), pointer :: explicit_grid
928) type(region_explicit_face_type), pointer :: faceset
929)
930) character(len=MAXSTRINGLENGTH) :: string
931) PetscInt :: icell, count
932) PetscInt, allocatable :: int_array(:)
933) PetscReal, allocatable :: real_array_2d(:,:)
934) PetscErrorCode :: ierr
935)
936) explicit_grid => ugrid%explicit_grid
937) faceset => region%explicit_faceset
938)
939) ! convert ids to petsc
940) region%cell_ids = region%cell_ids - 1
941) call AOApplicationToPetsc(ugrid%ao_natural_to_petsc,size(region%cell_ids), &
942) region%cell_ids,ierr);CHKERRQ(ierr)
943) region%cell_ids = region%cell_ids + 1
944)
945) ! if petsc ids are below global_offset or above global_offset + nlmax,
946) ! they are off processor; otherwise, local
947)
948) allocate(int_array(size(region%cell_ids)))
949) ! negate off processor ids
950) int_array = UNINITIALIZED_INTEGER
951) ! only if faceset exists
952) if (associated(faceset)) then
953) allocate(real_array_2d(4,size(region%cell_ids)))
954) real_array_2d = UNINITIALIZED_DOUBLE
955) endif
956) count = 0
957) do icell = 1, size(region%cell_ids)
958) if (region%cell_ids(icell) > ugrid%global_offset .and. &
959) region%cell_ids(icell) <= ugrid%global_offset + ugrid%nlmax) then
960) count = count + 1
961) ! local cell id
962) int_array(count) = region%cell_ids(icell) - ugrid%global_offset
963) if (associated(faceset)) then
964) real_array_2d(1,count) = faceset%face_centroids(icell)%x
965) real_array_2d(2,count) = faceset%face_centroids(icell)%y
966) real_array_2d(3,count) = faceset%face_centroids(icell)%z
967) real_array_2d(4,count) = faceset%face_areas(icell)
968) endif
969) endif
970) enddo
971)
972) deallocate(region%cell_ids)
973) allocate(region%cell_ids(count))
974) region%cell_ids = int_array(1:count)
975) deallocate(int_array)
976)
977) if (associated(faceset)) then
978) deallocate(faceset%face_centroids)
979) deallocate(faceset%face_areas)
980) allocate(faceset%face_centroids(count))
981) allocate(faceset%face_areas(count))
982)
983) do icell = 1, count
984) faceset%face_centroids(icell)%x = real_array_2d(1,icell)
985) faceset%face_centroids(icell)%y = real_array_2d(2,icell)
986) faceset%face_centroids(icell)%z = real_array_2d(3,icell)
987) faceset%face_areas(icell) = real_array_2d(4,icell)
988) enddo
989) deallocate(real_array_2d)
990) endif
991)
992)
993) region%num_cells = count
994)
995) if (region%num_cells == 0) then
996) deallocate(region%cell_ids)
997) nullify(region%cell_ids)
998) if (associated(faceset)) then
999) deallocate(faceset%face_centroids)
1000) nullify(faceset%face_centroids)
1001) deallocate(faceset%face_areas)
1002) nullify(faceset%face_areas)
1003) ! note that have to use full reference
1004) deallocate(region%explicit_faceset)
1005) nullify(region%explicit_faceset)
1006) endif
1007) endif
1008)
1009) #if UGRID_DEBUG
1010) if (region%num_cells > 0) then
1011) write(string,*) option%myrank
1012) string = 'region_faceset_' // trim(region%name) // trim(adjustl(string)) // '.out'
1013) open(unit=86,file=trim(string))
1014) if (associated(faceset)) then
1015) do icell = 1, region%num_cells
1016) write(86,'(i5,4f7.3)') region%cell_ids(icell), &
1017) faceset%face_centroids(icell)%x, &
1018) faceset%face_centroids(icell)%y, &
1019) faceset%face_centroids(icell)%z, &
1020) faceset%face_areas(icell)
1021) enddo
1022) else
1023) do icell = 1, region%num_cells
1024) write(86,'(i5)') region%cell_ids(icell)
1025) enddo
1026) endif
1027) close(86)
1028) endif
1029) #endif
1030)
1031) end subroutine GridLocalizeExplicitFaceset
1032)
1033) ! ************************************************************************** !
1034)
1035) subroutine GridCopyIntegerArrayToVec(grid, array,vector,num_values)
1036) !
1037) ! Copies values from an integer array into a
1038) ! PETSc Vec
1039) !
1040) ! Author: Glenn Hammond
1041) ! Date: 12/18/07
1042) !
1043)
1044) implicit none
1045)
1046) #include "petsc/finclude/petscvec.h"
1047) #include "petsc/finclude/petscvec.h90"
1048)
1049) type(grid_type) :: grid
1050) PetscInt :: array(:)
1051) Vec :: vector
1052) PetscInt :: num_values
1053)
1054) PetscReal, pointer :: vec_ptr(:)
1055) PetscErrorCode :: ierr
1056)
1057) call VecGetArrayF90( vector,vec_ptr,ierr);CHKERRQ(ierr)
1058) vec_ptr(1:num_values) = array(1:num_values)
1059) call VecRestoreArrayF90( vector,vec_ptr,ierr);CHKERRQ(ierr)
1060)
1061) end subroutine GridCopyIntegerArrayToVec
1062)
1063) ! ************************************************************************** !
1064)
1065) subroutine GridCopyRealArrayToVec(grid,array,vector,num_values)
1066) !
1067) ! Copies values from an integer array into a
1068) ! PETSc Vec
1069) !
1070) ! Author: Glenn Hammond
1071) ! Date: 12/18/07
1072) !
1073)
1074) implicit none
1075)
1076) #include "petsc/finclude/petscvec.h"
1077) #include "petsc/finclude/petscvec.h90"
1078)
1079) type(grid_type) :: grid
1080) PetscReal :: array(:)
1081) Vec :: vector
1082) PetscInt :: num_values
1083)
1084) PetscReal, pointer :: vec_ptr(:)
1085) PetscErrorCode :: ierr
1086)
1087) call VecGetArrayF90(vector,vec_ptr,ierr);CHKERRQ(ierr)
1088) vec_ptr(1:num_values) = array(1:num_values)
1089) call VecRestoreArrayF90(vector,vec_ptr,ierr);CHKERRQ(ierr)
1090)
1091) end subroutine GridCopyRealArrayToVec
1092)
1093) ! ************************************************************************** !
1094)
1095) subroutine GridCopyVecToIntegerArray(grid,array,vector,num_values)
1096) !
1097) ! Copies values from a PETSc Vec to an
1098) ! integer array
1099) !
1100) ! Author: Glenn Hammond
1101) ! Date: 12/18/07
1102) !
1103)
1104) implicit none
1105)
1106) #include "petsc/finclude/petscvec.h"
1107) #include "petsc/finclude/petscvec.h90"
1108)
1109) type(grid_type) :: grid
1110) PetscInt :: array(:)
1111) Vec :: vector
1112) PetscInt :: num_values
1113)
1114) PetscInt :: i
1115) PetscReal, pointer :: vec_ptr(:)
1116) PetscErrorCode :: ierr
1117)
1118) call VecGetArrayF90(vector,vec_ptr,ierr);CHKERRQ(ierr)
1119) do i=1,num_values
1120) if (vec_ptr(i) > 0.d0) then
1121) array(i) = int(vec_ptr(i)+1.d-4)
1122) else
1123) array(i) = int(vec_ptr(i)-1.d-4)
1124) endif
1125) enddo
1126) call VecRestoreArrayF90(vector,vec_ptr,ierr);CHKERRQ(ierr)
1127)
1128) end subroutine GridCopyVecToIntegerArray
1129)
1130) ! ************************************************************************** !
1131)
1132) subroutine GridCopyVecToRealArray(grid,array,vector,num_values)
1133) !
1134) ! Copies values from a PETSc Vec to an integer
1135) ! array
1136) !
1137) ! Author: Glenn Hammond
1138) ! Date: 12/18/07
1139) !
1140)
1141) implicit none
1142)
1143) #include "petsc/finclude/petscvec.h"
1144) #include "petsc/finclude/petscvec.h90"
1145)
1146) type(grid_type) :: grid
1147) PetscReal :: array(:)
1148) Vec :: vector
1149) PetscInt :: num_values
1150)
1151) PetscReal, pointer :: vec_ptr(:)
1152) PetscErrorCode :: ierr
1153)
1154) call VecGetArrayF90(vector,vec_ptr,ierr);CHKERRQ(ierr)
1155) array(1:num_values) = vec_ptr(1:num_values)
1156) call VecRestoreArrayF90(vector,vec_ptr,ierr);CHKERRQ(ierr)
1157)
1158) end subroutine GridCopyVecToRealArray
1159)
1160) ! ************************************************************************** !
1161)
1162) subroutine GridCreateNaturalToGhostedHash(grid,option)
1163) !
1164) ! Creates a hash table for looking up the
1165) ! local ghosted id of a natural id, if it
1166) ! exists
1167) !
1168) ! Author: Glenn Hammond
1169) ! Date: 03/07/07
1170) !
1171)
1172) use Option_module
1173) use Logging_module
1174)
1175) implicit none
1176)
1177) type(grid_type) :: grid
1178) type(option_type) :: option
1179)
1180) character(len=MAXSTRINGLENGTH) :: string
1181) PetscInt :: local_ghosted_id, natural_id
1182) PetscInt :: num_in_hash, num_ids_per_hash, hash_id, id, ierr, hash_id_2
1183) PetscInt :: max_num_ids_per_hash
1184) PetscInt, pointer :: hash(:,:,:), temp_hash(:,:,:)
1185)
1186) if (associated(grid%hash)) return
1187)
1188) call PetscLogEventBegin(logging%event_hash_create,ierr);CHKERRQ(ierr)
1189)
1190) max_num_ids_per_hash = 0
1191) ! initial guess of 10% of ids per hash
1192) ! must be at least 5 so that reallocation (*1.2) works below
1193) num_ids_per_hash = max(grid%nlmax/(grid%num_hash_bins/10),5)
1194)
1195) allocate(hash(2,0:num_ids_per_hash,grid%num_hash_bins))
1196) hash(:,:,:) = 0
1197)
1198)
1199) do local_ghosted_id = 1, grid%ngmax
1200) natural_id = grid%nG2A(local_ghosted_id) !nG2A is 1-based
1201) hash_id = mod(natural_id,grid%num_hash_bins)+1
1202) num_in_hash = hash(1,0,hash_id)
1203) num_in_hash = num_in_hash+1
1204) if (num_in_hash > max_num_ids_per_hash) max_num_ids_per_hash = num_in_hash
1205) ! if a hash runs out of space reallocate
1206) if (num_in_hash > num_ids_per_hash) then
1207) allocate(temp_hash(2,0:num_ids_per_hash,grid%num_hash_bins))
1208) ! copy old hash
1209) temp_hash(1:2,0:num_ids_per_hash,1:grid%num_hash_bins) = &
1210) hash(1:2,0:num_ids_per_hash,1:grid%num_hash_bins)
1211) deallocate(hash)
1212) ! recompute hash 20% larger
1213) num_ids_per_hash = int(dble(num_ids_per_hash)*1.2)
1214) allocate(hash(2,0:num_ids_per_hash,grid%num_hash_bins))
1215) ! copy old to new
1216) do hash_id_2 = 1, grid%num_hash_bins
1217) do id = 1, temp_hash(1,0,hash_id_2)
1218) hash(1:2,id,hash_id_2) = temp_hash(1:2,id,hash_id_2)
1219) enddo
1220) hash(1,0,hash_id_2) = temp_hash(1,0,hash_id_2)
1221) enddo
1222) deallocate(temp_hash)
1223) endif
1224) hash(1,0,hash_id) = num_in_hash
1225) hash(1,num_in_hash,hash_id) = natural_id
1226) hash(2,num_in_hash,hash_id) = local_ghosted_id
1227) enddo
1228)
1229) grid%hash => hash
1230)
1231) ! call GridPrintHashTable(grid)
1232) call MPI_Allreduce(max_num_ids_per_hash,num_in_hash,ONE_INTEGER_MPI, &
1233) MPIU_INTEGER,MPI_MAX,option%mycomm,ierr)
1234) write(option%io_buffer,'("max_num_ids_per_hash: ",i5)') num_in_hash
1235) call printMsg(option)
1236)
1237) call PetscLogEventEnd(logging%event_hash_create,ierr);CHKERRQ(ierr)
1238)
1239) end subroutine GridCreateNaturalToGhostedHash
1240)
1241) ! ************************************************************************** !
1242)
1243) PetscInt function GridGetLocalIdFromNaturalId(grid,natural_id)
1244) !
1245) ! GetLocalIdFromNaturalId: Returns the local id corresponding to a natural
1246) ! id or 0, if the natural id is off-processor
1247) ! WARNING: Extremely inefficient for large jobs
1248) !
1249) ! Author: Glenn Hammond
1250) ! Date: 03/07/07
1251) !
1252)
1253) implicit none
1254)
1255) type(grid_type) :: grid
1256)
1257) PetscInt :: natural_id, local_id
1258)
1259) do local_id = 1, grid%nlmax
1260) if (natural_id == grid%nG2A(grid%nL2G(local_id))) then
1261) GridGetLocalIdFromNaturalId = local_id
1262) return
1263) endif
1264) enddo
1265) GridGetLocalIdFromNaturalId = 0
1266)
1267) end function GridGetLocalIdFromNaturalId
1268)
1269) ! ************************************************************************** !
1270)
1271) PetscInt function GridGetLocalGhostedIdFromNatId(grid,natural_id)
1272) !
1273) ! Returns the local ghosted id corresponding
1274) ! to a natural id or 0, if the natural id
1275) ! is off-processor
1276) ! WARNING: Extremely inefficient for large jobs
1277) !
1278) ! Author: Glenn Hammond
1279) ! Date: 03/07/07
1280) !
1281)
1282) implicit none
1283)
1284) type(grid_type) :: grid
1285) PetscInt :: natural_id
1286)
1287) PetscInt :: local_ghosted_id
1288)
1289) do local_ghosted_id = 1, grid%ngmax
1290) !geh: nG2A is 1-based
1291) if (natural_id == grid%nG2A(local_ghosted_id)) then
1292) GridGetLocalGhostedIdFromNatId = local_ghosted_id
1293) return
1294) endif
1295) enddo
1296) GridGetLocalGhostedIdFromNatId = 0
1297)
1298) end function GridGetLocalGhostedIdFromNatId
1299)
1300) ! ************************************************************************** !
1301)
1302) PetscInt function GridGetLocalGhostedIdFromHash(grid,natural_id)
1303) !
1304) ! Returns the local ghosted id of a natural
1305) ! id, if it exists. Otherwise 0 is returned
1306) !
1307) ! Author: Glenn Hammond
1308) ! Date: 03/07/07
1309) !
1310)
1311) implicit none
1312)
1313) type(grid_type) :: grid
1314) PetscInt :: natural_id
1315)
1316) PetscInt :: hash_id, id
1317)
1318) GridGetLocalGhostedIdFromHash = 0
1319) hash_id = mod(natural_id,grid%num_hash_bins)+1
1320) do id = 1, grid%hash(1,0,hash_id)
1321) if (grid%hash(1,id,hash_id) == natural_id) then
1322) GridGetLocalGhostedIdFromHash = grid%hash(2,id,hash_id)
1323) return
1324) endif
1325) enddo
1326)
1327) end function GridGetLocalGhostedIdFromHash
1328)
1329) ! ************************************************************************** !
1330)
1331) subroutine GridDestroyHashTable(grid)
1332) !
1333) ! Deallocates the hash table
1334) !
1335) ! Author: Glenn Hammond
1336) ! Date: 03/07/07
1337) !
1338) use Utility_module, only : DeallocateArray
1339)
1340) implicit none
1341)
1342) type(grid_type), pointer :: grid
1343)
1344) call DeallocateArray(grid%hash)
1345)
1346) nullify(grid%hash)
1347) grid%num_hash_bins = 100
1348)
1349) end subroutine GridDestroyHashTable
1350)
1351) ! ************************************************************************** !
1352)
1353) subroutine GridPrintHashTable(grid)
1354) !
1355) ! UnstructGridPrintHashTable: Prints the hashtable for viewing
1356) !
1357) ! Author: Glenn Hammond
1358) ! Date: 03/09/07
1359) !
1360)
1361) implicit none
1362)
1363) type(grid_type) :: grid
1364)
1365) PetscInt :: ihash, id, fid
1366)
1367) fid = 87
1368) open(fid,file='hashtable.dat',action='write')
1369) do ihash=1,grid%num_hash_bins
1370) write(fid,'(a4,i3,a,i5,a2,x,200(i6,x))') 'Hash',ihash,'(', &
1371) grid%hash(1,0,ihash), &
1372) '):', &
1373) (grid%hash(1,id,ihash),id=1,grid%hash(1,0,ihash))
1374) enddo
1375) close(fid)
1376)
1377) end subroutine GridPrintHashTable
1378)
1379) ! ************************************************************************** !
1380)
1381) subroutine GridGetGhostedNeighbors(grid,ghosted_id,stencil_type, &
1382) stencil_width_i,stencil_width_j, &
1383) stencil_width_k,x_count,y_count, &
1384) z_count, &
1385) ghosted_neighbors,option)
1386) !
1387) ! GridGetNeighbors: Returns an array of neighboring cells
1388) !
1389) ! Author: Glenn Hammond
1390) ! Date: 01/28/11
1391) !
1392)
1393) use Option_module
1394)
1395) implicit none
1396)
1397) type(grid_type) :: grid
1398) type(option_type) :: option
1399) PetscInt :: ghosted_id
1400) PetscEnum :: stencil_type
1401) PetscInt :: stencil_width_i
1402) PetscInt :: stencil_width_j
1403) PetscInt :: stencil_width_k
1404) PetscInt :: x_count
1405) PetscInt :: y_count
1406) PetscInt :: z_count
1407) PetscInt :: ghosted_neighbors(*)
1408)
1409) select case(grid%itype)
1410) case(STRUCTURED_GRID)
1411) call StructGridGetGhostedNeighbors(grid%structured_grid, &
1412) ghosted_id,stencil_type, &
1413) stencil_width_i, &
1414) stencil_width_j,stencil_width_k, &
1415) x_count,y_count,z_count, &
1416) ghosted_neighbors,option)
1417) case(IMPLICIT_UNSTRUCTURED_GRID,EXPLICIT_UNSTRUCTURED_GRID)
1418) option%io_buffer = 'GridGetNeighbors not currently supported for ' // &
1419) 'unstructured grids.'
1420) call printErrMsg(option)
1421) end select
1422)
1423) end subroutine GridGetGhostedNeighbors
1424)
1425) ! ************************************************************************** !
1426)
1427) subroutine GridGetGhostedNeighborsWithCorners(grid,ghosted_id,stencil_type, &
1428) stencil_width_i,stencil_width_j, &
1429) stencil_width_k,icount, &
1430) ghosted_neighbors,option)
1431) !
1432) ! GridGetNeighborsWithCorners: Returns an array of neighboring cells along with corner
1433) ! cells
1434) !
1435) ! Author: Satish Karra, LANL
1436) ! Date: 02/19/12
1437) !
1438)
1439) use Option_module
1440)
1441) implicit none
1442)
1443) type(grid_type) :: grid
1444) type(option_type) :: option
1445) PetscInt :: ghosted_id
1446) PetscEnum :: stencil_type
1447) PetscInt :: stencil_width_i
1448) PetscInt :: stencil_width_j
1449) PetscInt :: stencil_width_k
1450) PetscInt :: icount
1451) PetscInt :: ghosted_neighbors(*)
1452)
1453) select case(grid%itype)
1454) case(STRUCTURED_GRID)
1455) call StructGridGetGhostedNeighborsCorners(grid%structured_grid, &
1456) ghosted_id,stencil_type, &
1457) stencil_width_i, &
1458) stencil_width_j, &
1459) stencil_width_k, &
1460) icount, &
1461) ghosted_neighbors,option)
1462) case(IMPLICIT_UNSTRUCTURED_GRID,EXPLICIT_UNSTRUCTURED_GRID)
1463) option%io_buffer = 'GridGetNeighbors not currently supported for ' // &
1464) 'unstructured grids.'
1465) call printErrMsg(option)
1466) end select
1467)
1468) end subroutine GridGetGhostedNeighborsWithCorners
1469)
1470) ! ************************************************************************** !
1471)
1472) subroutine GridDestroy(grid)
1473) !
1474) ! Deallocates a grid
1475) !
1476) ! Author: Glenn Hammond
1477) ! Date: 11/01/07
1478) !
1479) use Utility_module, only : DeallocateArray
1480)
1481) implicit none
1482)
1483) type(grid_type), pointer :: grid
1484) PetscErrorCode :: ierr
1485) PetscInt :: ghosted_id
1486)
1487) if (.not.associated(grid)) return
1488)
1489) call DeallocateArray(grid%nL2G)
1490) call DeallocateArray(grid%nG2L)
1491) call DeallocateArray(grid%nG2A)
1492)
1493) call DeallocateArray(grid%x)
1494) call DeallocateArray(grid%y)
1495) call DeallocateArray(grid%z)
1496)
1497) if (associated(grid%hash)) call GridDestroyHashTable(grid)
1498)
1499) call UGridDestroy(grid%unstructured_grid)
1500) call StructGridDestroy(grid%structured_grid)
1501)
1502) call ConnectionDestroyList(grid%internal_connection_set_list)
1503)
1504) if (associated(grid)) deallocate(grid)
1505) nullify(grid)
1506)
1507) end subroutine GridDestroy
1508)
1509) ! ************************************************************************** !
1510)
1511) function GridIndexToCellID(vec,index,grid,vec_type)
1512) !
1513) ! Returns the local grid cell id of a Petsc Vec index
1514) !
1515) ! Author: Glenn Hammond
1516) ! Date: 01/07/09
1517) !
1518)
1519) implicit none
1520)
1521) Vec :: vec
1522) PetscInt :: index
1523) type(grid_type) :: grid
1524) PetscInt :: vec_type
1525)
1526) PetscInt :: GridIndexToCellID
1527)
1528) PetscInt :: low
1529) PetscInt :: high
1530) PetscInt :: ndof
1531) PetscInt :: cell_id
1532) PetscErrorCode :: ierr
1533)
1534)
1535) cell_id = -1
1536) call VecGetOwnershipRange(vec,low,high,ierr);CHKERRQ(ierr)
1537) call VecGetBlockSize(vec,ndof,ierr);CHKERRQ(ierr)
1538) if (index >= low .and. index < high) then
1539) cell_id = (index-low)/ndof+1
1540) if (vec_type == GLOBAL) then
1541) cell_id = grid%nG2A(grid%nL2G(cell_id))
1542) else if (vec_type == LOCAL) then
1543) cell_id = grid%nG2A(cell_id) !nG2A is 1-based
1544) endif
1545) endif
1546)
1547) call MPI_Allreduce(cell_id,GridIndexToCellID,ONE_INTEGER_MPI,MPIU_INTEGER, &
1548) MPI_MAX,PETSC_COMM_WORLD,ierr)
1549)
1550) end function GridIndexToCellID
1551)
1552) ! ************************************************************************** !
1553)
1554) subroutine GridLocalizeRegionFromBlock(grid,region,option)
1555) !
1556) ! This routine resticts regions to cells local to processor when the region
1557) ! was defined using a BLOCK from inputfile.
1558) !
1559) ! Author: Gautam Bisht, LBNL
1560) ! Date: 09/04/12
1561) !
1562)
1563) use Option_module
1564) use Region_module
1565)
1566) implicit none
1567)
1568) type(region_type), pointer :: region
1569) type(grid_type), pointer :: grid
1570) type(option_type) :: option
1571)
1572) character(len=MAXSTRINGLENGTH) :: string
1573) PetscInt, allocatable :: temp_int_array(:)
1574) PetscInt :: i, j, k, count, local_count, ghosted_id, local_id
1575) PetscInt :: i_min, i_max, j_min, j_max, k_min, k_max
1576) PetscReal :: x_min, x_max, y_min, y_max, z_min, z_max
1577) PetscReal, parameter :: pert = 1.d-8, tol = 1.d-20
1578) PetscReal :: x_shift, y_shift, z_shift
1579) PetscReal :: del_x, del_y, del_z
1580) PetscInt :: iflag
1581) PetscBool :: same_point
1582) PetscErrorCode :: ierr
1583)
1584) if (grid%itype /= STRUCTURED_GRID) then
1585) option%io_buffer='Region definition using BLOCK is only supported for ' //&
1586) ' structured grids'
1587) call printErrMsg(option)
1588) endif
1589)
1590) ! convert indexing from global (entire domain) to local processor
1591) region%i1 = region%i1 - grid%structured_grid%lxs
1592) region%i2 = region%i2 - grid%structured_grid%lxs
1593) region%j1 = region%j1 - grid%structured_grid%lys
1594) region%j2 = region%j2 - grid%structured_grid%lys
1595) region%k1 = region%k1 - grid%structured_grid%lzs
1596) region%k2 = region%k2 - grid%structured_grid%lzs
1597)
1598) ! clip region to within local processor domain
1599) region%i1 = max(region%i1,1)
1600) region%i2 = min(region%i2,grid%structured_grid%nlx)
1601) region%j1 = max(region%j1,1)
1602) region%j2 = min(region%j2,grid%structured_grid%nly)
1603) region%k1 = max(region%k1,1)
1604) region%k2 = min(region%k2,grid%structured_grid%nlz)
1605)
1606) count = 0
1607) if (region%i1 <= region%i2 .and. &
1608) region%j1 <= region%j2 .and. &
1609) region%k1 <= region%k2) then
1610) region%num_cells = (region%i2-region%i1+1)* &
1611) (region%j2-region%j1+1)* &
1612) (region%k2-region%k1+1)
1613) allocate(region%cell_ids(region%num_cells))
1614) if (region%iface /= 0) then
1615) allocate(region%faces(region%num_cells))
1616) region%faces = region%iface
1617) endif
1618) region%cell_ids = 0
1619) do k=region%k1,region%k2
1620) do j=region%j1,region%j2
1621) do i=region%i1,region%i2
1622) count = count + 1
1623) region%cell_ids(count) = &
1624) i + (j-1)*grid%structured_grid%nlx + &
1625) (k-1)*grid%structured_grid%nlxy
1626) enddo
1627) enddo
1628) enddo
1629) ! if (region%num_cells > 0) then
1630) ! region%coordinates(1)%x = grid%x(region%cell_ids(ONE_INTEGER))
1631) ! region%coordinates(1)%y = grid%y(region%cell_ids(ONE_INTEGER))
1632) ! region%coordinates(1)%z = grid%z(region%cell_ids(ONE_INTEGER))
1633) ! endif
1634) else
1635) region%num_cells = 0
1636) endif
1637)
1638) if (count /= region%num_cells) then
1639) option%io_buffer = 'Mismatch in number of cells in block region'
1640) call printErrMsg(option)
1641) endif
1642)
1643) end subroutine GridLocalizeRegionFromBlock
1644)
1645) ! ************************************************************************** !
1646)
1647) subroutine GridLocalizeRegionFromCartBound(grid,region,option)
1648) !
1649) ! This routine resticts regions to cells local to processor when the region
1650) ! was defined using a BLOCK from inputfile.
1651) !
1652) ! Author: Gautam Bisht, LBNL
1653) ! Date: 09/04/12
1654) !
1655)
1656) use Option_module
1657) use Region_module
1658)
1659) implicit none
1660)
1661) type(region_type), pointer :: region
1662) type(grid_type), pointer :: grid
1663) type(option_type) :: option
1664)
1665) character(len=MAXSTRINGLENGTH) :: string
1666) PetscInt, allocatable :: temp_int_array(:)
1667) PetscInt :: i, j, k, count, local_count, ghosted_id, local_id
1668) PetscInt :: i_min, i_max, j_min, j_max, k_min, k_max
1669) PetscReal :: x_min, x_max, y_min, y_max, z_min, z_max
1670) PetscReal, parameter :: pert = 1.d-8, tol = 1.d-20
1671) PetscReal :: x_shift, y_shift, z_shift
1672) PetscReal :: del_x, del_y, del_z
1673) PetscInt :: iflag
1674) PetscBool :: same_point
1675) PetscErrorCode :: ierr
1676)
1677) if (grid%itype /= STRUCTURED_GRID) then
1678) option%io_buffer='Region definition using CARTESIAN_BOUNDARY is ' // &
1679) 'only supported for structured grids.'
1680) call printErrMsg(option)
1681) endif
1682)
1683) region%i1 = 1
1684) region%i2 = grid%structured_grid%nx
1685) region%j1 = 1
1686) region%j2 = grid%structured_grid%ny
1687) region%k1 = 1
1688) region%k2 = grid%structured_grid%nz
1689)
1690) select case(region%iface)
1691) case(WEST_FACE)
1692) region%i2 = region%i1
1693) case(EAST_FACE)
1694) region%i1 = region%i2
1695) case(SOUTH_FACE)
1696) region%j2 = region%j1
1697) case(NORTH_FACE)
1698) region%j1 = region%j2
1699) case(BOTTOM_FACE)
1700) region%k2 = region%k1
1701) case(TOP_FACE)
1702) region%k1 = region%k2
1703) end select
1704)
1705) call GridLocalizeRegionFromBlock(grid,region,option)
1706)
1707) end subroutine GridLocalizeRegionFromCartBound
1708)
1709) ! ************************************************************************** !
1710)
1711) subroutine GridLocalizeRegionFromCoordinates(grid,region,option)
1712) !
1713) ! This routine resticts regions to cells local to processor when the region
1714) ! was defined using COORDINATES from inputfile.
1715) !
1716) ! Author: Gautam Bisht, LBNL
1717) ! Date: 09/04/12
1718) !
1719)
1720) use Option_module
1721) use Region_module
1722)
1723) implicit none
1724)
1725) type(region_type), pointer :: region
1726) type(grid_type), pointer :: grid
1727) type(option_type) :: option
1728)
1729) character(len=MAXSTRINGLENGTH) :: string
1730) PetscInt, allocatable :: temp_int_array(:)
1731) PetscInt :: i, j, k, count, local_count, ghosted_id, local_id
1732) PetscInt :: i_min, i_max, j_min, j_max, k_min, k_max
1733) PetscReal :: x_min, x_max, y_min, y_max, z_min, z_max
1734) PetscReal, parameter :: pert = 1.d-8, tol = 1.d-20
1735) PetscReal :: x_shift, y_shift, z_shift
1736) PetscReal :: del_x, del_y, del_z
1737) PetscInt :: iflag
1738) PetscBool :: same_point
1739) PetscErrorCode :: ierr
1740)
1741) iflag = 0
1742) if (size(region%coordinates) > TWO_INTEGER) then
1743) option%io_buffer = 'GridLocalizeRegions: more than 2 coordinates' // &
1744) ' not supported in region object'
1745) call printErrMsg(option)
1746) endif
1747)
1748) same_point = PETSC_FALSE
1749) ! if two coordinates, determine whether they are the same point
1750) if (size(region%coordinates) == TWO_INTEGER) then
1751) if (dabs(region%coordinates(ONE_INTEGER)%x - &
1752) region%coordinates(TWO_INTEGER)%x) < tol .and. &
1753) dabs(region%coordinates(ONE_INTEGER)%y - &
1754) region%coordinates(TWO_INTEGER)%y) < tol .and. &
1755) dabs(region%coordinates(ONE_INTEGER)%z - &
1756) region%coordinates(TWO_INTEGER)%z) < tol) then
1757) same_point = PETSC_TRUE
1758) endif
1759) endif
1760)
1761) ! treat two identical coordinates the same as a single coordinate
1762) if (size(region%coordinates) == ONE_INTEGER .or. same_point) then
1763) call GridGetLocalIDFromCoordinate(grid,region%coordinates(ONE_INTEGER), &
1764) option,local_id)
1765) if (INITIALIZED(local_id)) then
1766) region%num_cells = 1
1767) allocate(region%cell_ids(region%num_cells))
1768) region%cell_ids(1) = local_id
1769) select case(grid%itype)
1770) case(STRUCTURED_GRID)
1771) if (region%iface /= 0) then
1772) allocate(region%faces(region%num_cells))
1773) region%faces = region%iface
1774) endif
1775) end select
1776) else
1777) region%num_cells = 0
1778) endif
1779) ! the next test as designed will only work on a uniform grid
1780) call MPI_Allreduce(region%num_cells,count, &
1781) ONE_INTEGER_MPI,MPIU_INTEGER,MPI_SUM, &
1782) option%mycomm,ierr)
1783) if (count == 0) then
1784) write(option%io_buffer,*) 'Region: (coord)', &
1785) region%coordinates(ONE_INTEGER)%x, &
1786) region%coordinates(ONE_INTEGER)%y, &
1787) region%coordinates(ONE_INTEGER)%z, &
1788) ' not found in global domain.', count
1789) call printErrMsg(option)
1790) else if (count > 1) then
1791) write(option%io_buffer,*) 'Region: (coord)', &
1792) region%coordinates(ONE_INTEGER)%x, &
1793) region%coordinates(ONE_INTEGER)%y, &
1794) region%coordinates(ONE_INTEGER)%z, &
1795) ' duplicated across ', count, &
1796) ' procs in global domain.'
1797) call printErrMsg(option)
1798) endif
1799) else ! 2 coordinates
1800) x_min = min(region%coordinates(ONE_INTEGER)%x, &
1801) region%coordinates(TWO_INTEGER)%x)
1802) x_max = max(region%coordinates(ONE_INTEGER)%x, &
1803) region%coordinates(TWO_INTEGER)%x)
1804) y_min = min(region%coordinates(ONE_INTEGER)%y, &
1805) region%coordinates(TWO_INTEGER)%y)
1806) y_max = max(region%coordinates(ONE_INTEGER)%y, &
1807) region%coordinates(TWO_INTEGER)%y)
1808) z_min = min(region%coordinates(ONE_INTEGER)%z, &
1809) region%coordinates(TWO_INTEGER)%z)
1810) z_max = max(region%coordinates(ONE_INTEGER)%z, &
1811) region%coordinates(TWO_INTEGER)%z)
1812)
1813) if (grid%itype == STRUCTURED_GRID) then
1814) ! shift box slightly inward
1815) x_shift = 1.d-8*(grid%x_max_global-grid%x_min_global)
1816) x_min = x_min+x_shift
1817) x_max = x_max-x_shift
1818) y_shift = 1.d-8*(grid%y_max_global-grid%y_min_global)
1819) y_min = y_min+y_shift
1820) y_max = y_max-y_shift
1821) z_shift = 1.d-8*(grid%z_max_global-grid%z_min_global)
1822) z_min = z_min+z_shift
1823) z_max = z_max-z_shift
1824)
1825) ! if plane or line, ensure it is within the grid cells
1826) if (x_max-x_min < 1.d-10) then
1827) x_max = region%coordinates(ONE_INTEGER)%x
1828) x_shift = 1.d-8*(grid%x_max_global-grid%x_min_global)
1829) if (region%iface == WEST_FACE) then
1830) x_max = x_max + x_shift
1831) else if (region%iface == EAST_FACE) then
1832) x_max = x_max - x_shift
1833) ! otherwise, shift upwind, unless at upwind physical boundary
1834) else
1835) if (x_max > grid%x_min_global + x_shift) then
1836) x_max = x_max - x_shift
1837) else
1838) x_max = x_max + x_shift
1839) endif
1840) endif
1841) x_min = x_max
1842) endif
1843) if (y_max-y_min < 1.d-10) then
1844) y_max = region%coordinates(ONE_INTEGER)%y
1845) y_shift = 1.d-8*(grid%y_max_global-grid%y_min_global)
1846) if (region%iface == SOUTH_FACE) then
1847) y_max = y_max + y_shift
1848) else if (region%iface == NORTH_FACE) then
1849) y_max = y_max - y_shift
1850) ! otherwise, shift upwind, unless at upwind physical boundary
1851) else
1852) if (y_max > grid%y_min_global + y_shift) then
1853) y_max = y_max - y_shift
1854) else
1855) y_max = y_max + y_shift
1856) endif
1857) endif
1858) y_min = y_max
1859) endif
1860) if (z_max-z_min < 1.d-10) then
1861) z_max = region%coordinates(ONE_INTEGER)%z
1862) z_shift = 1.d-8*(grid%z_max_global-grid%z_min_global)
1863) if (region%iface == BOTTOM_FACE) then
1864) z_max = z_max + z_shift
1865) else if (region%iface == TOP_FACE) then
1866) z_max = z_max - z_shift
1867) ! otherwise, shift upwind, unless at upwind physical boundary
1868) else
1869) if (z_max > grid%z_min_global + z_shift) then
1870) z_max = z_max - z_shift
1871) else
1872) z_max = z_max + z_shift
1873) endif
1874) endif
1875) z_min = z_max
1876) endif
1877) endif
1878)
1879) ! ensure overlap
1880) if (x_min <= grid%x_max_local .and. &
1881) x_max >= grid%x_min_local .and. &
1882) y_min <= grid%y_max_local .and. &
1883) y_max >= grid%y_min_local .and. &
1884) z_min <= grid%z_max_local .and. &
1885) z_max >= grid%z_min_local) then
1886)
1887) ! get I,J,K bounds
1888) select case(grid%itype)
1889) case(STRUCTURED_GRID)
1890) ! local, non-ghosted i,j,k's are returned
1891) call StructGridGetIJKFromCoordinate(grid%structured_grid, &
1892) max(x_min,grid%x_min_local+x_shift), &
1893) max(y_min,grid%y_min_local+y_shift), &
1894) max(z_min,grid%z_min_local+z_shift), &
1895) i_min,j_min,k_min)
1896) call StructGridGetIJKFromCoordinate(grid%structured_grid, &
1897) min(x_max,grid%x_max_local-x_shift), &
1898) min(y_max,grid%y_max_local-y_shift), &
1899) min(z_max,grid%z_max_local-z_shift), &
1900) i_max,j_max,k_max)
1901) if (i_min > 0 .and. j_min > 0 .and. k_min > 0 .and. &
1902) i_max > 0 .and. j_max > 0 .and. k_max > 0) then
1903) region%num_cells = (i_max-i_min+1)*(j_max-j_min+1)*(k_max-k_min+1)
1904) allocate(region%cell_ids(region%num_cells))
1905) if (region%iface /= 0) then
1906) allocate(region%faces(region%num_cells))
1907) region%faces = region%iface
1908) endif
1909) region%cell_ids = 0
1910) count = 0
1911) do k = k_min, k_max
1912) do j = j_min, j_max
1913) do i = i_min, i_max
1914) count = count+1
1915) region%cell_ids(count) = i + (j-1)*grid%structured_grid%nlx + &
1916) (k-1)*grid%structured_grid%nlxy
1917) enddo
1918) enddo
1919) enddo
1920) else
1921) iflag = 1
1922) endif
1923) case(IMPLICIT_UNSTRUCTURED_GRID,EXPLICIT_UNSTRUCTURED_GRID, &
1924) POLYHEDRA_UNSTRUCTURED_GRID)
1925) del_x = x_max-x_min
1926) del_y = y_max-y_min
1927) del_z = z_max-z_min
1928) ! 3D box
1929) if (del_x > 1.d-10 .and. &
1930) del_y > 1.d-10 .and. &
1931) del_z > 1.d-10) then
1932) ! geh: if the coordinates define a 3D box, add all cell centers
1933) ! that reside within the box
1934) count = 0
1935) do local_id = 1, grid%nlmax
1936) ghosted_id = grid%nL2G(local_id)
1937) if (grid%x(ghosted_id) >= x_min .and. &
1938) grid%x(ghosted_id) <= x_max .and. &
1939) grid%y(ghosted_id) >= y_min .and. &
1940) grid%y(ghosted_id) <= y_max .and. &
1941) grid%z(ghosted_id) >= z_min .and. &
1942) grid%z(ghosted_id) <= z_max) then
1943) count = count + 1
1944) endif
1945) enddo
1946) allocate(region%cell_ids(count))
1947) region%cell_ids = 0
1948) count = 0
1949) do local_id = 1, grid%nlmax
1950) ghosted_id = grid%nL2G(local_id)
1951) if (grid%x(ghosted_id) >= x_min .and. &
1952) grid%x(ghosted_id) <= x_max .and. &
1953) grid%y(ghosted_id) >= y_min .and. &
1954) grid%y(ghosted_id) <= y_max .and. &
1955) grid%z(ghosted_id) >= z_min .and. &
1956) grid%z(ghosted_id) <= z_max) then
1957) count = count + 1
1958) region%cell_ids(count) = local_id
1959) endif
1960) enddo
1961) region%num_cells = count
1962) ! 2D plane
1963) elseif ((del_x < 1.d-10 .and. del_y > 1.d-10 .and. &
1964) del_z > 1.d-10) .or. &
1965) (del_x > 1.d-10 .and. del_y < 1.d-10 .and. &
1966) del_z > 1.d-10) .or. &
1967) (del_x > 1.d-10 .and. del_y > 1.d-10 .and. &
1968) del_z < 1.d-10)) then
1969) if (grid%itype == EXPLICIT_UNSTRUCTURED_GRID) then
1970) option%io_buffer = 'Regions defined with 2D planes are not ' // &
1971) 'supported with explicit unstructured grids.'
1972) call printErrMsg(option)
1973) endif
1974) if (grid%itype == IMPLICIT_UNSTRUCTURED_GRID) then
1975) call UGridGetCellsInRectangle(x_min,x_max,y_min,y_max, &
1976) z_min,z_max, &
1977) grid%unstructured_grid,option, &
1978) region%num_cells,region%cell_ids, &
1979) region%faces)
1980) else
1981) call UGridPolyhedraGetCellsInRectangle(x_min,x_max,y_min,y_max, &
1982) z_min,z_max, &
1983) grid%unstructured_grid,option, &
1984) region%num_cells,region%cell_ids, &
1985) region%faces)
1986) endif
1987)
1988) endif
1989) end select
1990) endif
1991)
1992) call MPI_Allreduce(iflag,i,ONE_INTEGER_MPI,MPIU_INTEGER,MPI_MAX, &
1993) option%mycomm,ierr)
1994) iflag = i
1995) if (iflag > 0) then
1996) option%io_buffer = 'GridLocalizeRegions, between two points'
1997) call printErrMsg(option)
1998) endif
1999) endif
2000)
2001) end subroutine GridLocalizeRegionFromCoordinates
2002)
2003) ! ************************************************************************** !
2004)
2005) subroutine GridMapCellsInPolVol(grid,polygonal_volume, &
2006) region_name,option,cell_ids)
2007) !
2008) ! Maps all global boundary cells within a polygonal volume to a region
2009) !
2010) ! Author: Glenn Hammond
2011) ! Date: 10/16/15
2012) !
2013) use Option_module
2014) use Geometry_module
2015)
2016) implicit none
2017)
2018) type(grid_type) :: grid
2019) type(polygonal_volume_type) :: polygonal_volume
2020) character(len=MAXWORDLENGTH) :: region_name
2021) type(option_type) :: option
2022) PetscInt, pointer :: cell_ids(:)
2023)
2024) PetscInt :: local_id, ghosted_id, icount
2025) PetscBool :: found
2026) PetscInt, allocatable :: temp_int(:)
2027)
2028) allocate(temp_int(grid%nlmax))
2029) temp_int = 0
2030) icount = 0
2031) do local_id = 1, grid%nlmax
2032) ghosted_id = grid%nL2G(local_id)
2033) found = GeometryPointInPolygonalVolume(grid%x(ghosted_id), &
2034) grid%y(ghosted_id), &
2035) grid%z(ghosted_id), &
2036) polygonal_volume,option)
2037) if (found) then
2038) icount = icount + 1
2039) temp_int(icount) = local_id
2040) endif
2041) enddo
2042) allocate(cell_ids(icount))
2043) cell_ids = temp_int(1:icount)
2044) deallocate(temp_int)
2045)
2046) end subroutine GridMapCellsInPolVol
2047)
2048) ! ************************************************************************** !
2049)
2050) subroutine GridGetLocalIDFromCoordinate(grid,coordinate,option,local_id)
2051) !
2052) ! Returns the local id of the grid cell occupied by a coordinate
2053) !
2054) ! Author: Glenn Hammond
2055) ! Date: 10/16/15
2056) !
2057) use Option_module
2058) use Geometry_module
2059)
2060) implicit none
2061)
2062) type(grid_type) :: grid
2063) type(point3d_type) :: coordinate
2064) type(option_type) :: option
2065) PetscInt :: local_id
2066)
2067) PetscReal, parameter :: pert = 1.d-8, tol = 1.d-20
2068) PetscReal :: x_shift, y_shift, z_shift
2069) PetscInt :: i, j, k
2070)
2071) local_id = UNINITIALIZED_INTEGER
2072) if (coordinate%x >= grid%x_min_global .and. &
2073) coordinate%x <= grid%x_max_global .and. &
2074) coordinate%y >= grid%y_min_global .and. &
2075) coordinate%y <= grid%y_max_global .and. &
2076) coordinate%z >= grid%z_min_global .and. &
2077) coordinate%z <= grid%z_max_global) then
2078) ! If a point is on the corner of 4 or 8 patches in AMR, the region
2079) ! will be assigned to all 4/8 patches...a problem. To avoid this,
2080) ! we are going to perturb all point coordinates slightly upwind, as
2081) ! long as they are not on a global boundary (i.e. boundary condition)
2082) ! -- shift the coorindate slightly upwind
2083) x_shift = coordinate%x - &
2084) pert*(grid%x_max_global-grid%x_min_global)
2085) y_shift = coordinate%y - &
2086) pert*(grid%y_max_global-grid%y_min_global)
2087) z_shift = coordinate%z - &
2088) pert*(grid%z_max_global-grid%z_min_global)
2089) ! if the coodinate is shifted out of the global domain or
2090) ! onto an exterior edge, set it back to the original value
2091) if (x_shift - grid%x_min_global < tol) &
2092) x_shift = coordinate%x
2093) if (y_shift - grid%y_min_global < tol) &
2094) y_shift = coordinate%y
2095) if (z_shift - grid%z_min_global < tol) &
2096) z_shift = coordinate%z
2097) select case(grid%itype)
2098) case(STRUCTURED_GRID)
2099) call StructGridGetIJKFromCoordinate(grid%structured_grid, &
2100) x_shift,y_shift,z_shift, &
2101) i,j,k)
2102) if (i > 0 .and. j > 0 .and. k > 0) then
2103) local_id = i + (j-1)*grid%structured_grid%nlx + &
2104) (k-1)*grid%structured_grid%nlxy
2105) endif
2106) case(IMPLICIT_UNSTRUCTURED_GRID)
2107) !geh: must check each cell individually
2108) call UGridGetCellFromPoint(coordinate%x, &
2109) coordinate%y, &
2110) coordinate%z, &
2111) grid%unstructured_grid,option,local_id)
2112) case(EXPLICIT_UNSTRUCTURED_GRID)
2113) if (grid%itype == EXPLICIT_UNSTRUCTURED_GRID) then
2114) option%io_buffer = 'Locating a grid cell through a specified &
2115) &coordinate (GridGetLocalIDFromCoordinate)is not supported for &
2116) &explicit (primal) unstructured grids.'
2117) call printErrMsg(option)
2118) endif
2119) case(POLYHEDRA_UNSTRUCTURED_GRID)
2120) option%io_buffer = &
2121) 'add code POLYHDERA in GridGetLocalIDFromCoordinate'
2122) call printErrMsg(option)
2123) end select
2124) endif
2125)
2126) ! Several of the subroutines above may return local_id = 0, therefore, we
2127) ! need to ensure that local_id is reset back to uninitiazed if 0
2128) if (local_id <= 0) then
2129) local_id = UNINITIALIZED_INTEGER
2130) endif
2131)
2132) end subroutine GridGetLocalIDFromCoordinate
2133)
2134)
2135) end module Grid_module