grid_unstructured_aux.F90 coverage: 100.00 %func 77.57 %block
1) module Grid_Unstructured_Aux_module
2)
3) ! use Connection_module
4) use Grid_Unstructured_Cell_module
5) use Geometry_module
6)
7) use PFLOTRAN_Constants_module
8)
9) implicit none
10)
11) private
12)
13) #include "petsc/finclude/petscsys.h"
14) #include "petsc/finclude/petscvec.h"
15) #include "petsc/finclude/petscvec.h90"
16) #include "petsc/finclude/petscis.h"
17) #include "petsc/finclude/petscis.h90"
18) #if defined(SCORPIO)
19) include "scorpiof.h"
20) #endif
21)
22) type, public :: grid_unstructured_type
23) ! variables for all unstructured grids
24) PetscInt :: num_ghost_cells ! number of ghost cells (only) on processor
25) PetscInt :: global_offset ! offset in petsc ordering for the first cell on a processor???
26) PetscInt :: nmax ! Total number of nodes in global domain
27) PetscInt :: nlmax ! Total number of non-ghosted nodes in local domain.
28) PetscInt :: ngmax ! Number of ghosted & non-ghosted nodes in local domain.
29) PetscInt, pointer :: hash(:,:,:)
30) PetscInt :: num_hash
31) PetscInt, pointer :: cell_ids_natural(:) ! natural 1d right-hand i,j,k ordering
32) PetscInt, pointer :: cell_ids_petsc(:) ! petsc ordering of cell ids
33) PetscInt, pointer :: ghost_cell_ids_petsc(:) ! petsc ordering of ghost cells ids
34) AO :: ao_natural_to_petsc ! mapping of natural to petsc ordering
35) type(unstructured_explicit_type), pointer :: explicit_grid
36) type(unstructured_polyhedra_type), pointer :: polyhedra_grid
37) ! variables for implicit unstructured grids
38) PetscInt :: grid_type ! 3D subsurface (default) or 2D surface grid
39) PetscInt :: num_vertices_global ! number of vertices in entire problem domain
40) PetscInt :: num_vertices_local ! number of vertices in local grid cells
41) PetscInt :: num_vertices_natural ! number of vertices read initially
42) PetscInt :: max_ndual_per_cell
43) PetscInt :: max_nvert_per_cell
44) PetscInt :: max_cells_sharing_a_vertex
45) PetscInt, pointer :: cell_type(:)
46) PetscInt, pointer :: cell_vertices(:,:) ! vertices for each grid cell (NO LONGER zero-based)
47) PetscInt, pointer :: face_to_cell_ghosted(:,:) !
48) PetscInt, pointer :: connection_to_face(:)
49) !geh: Should not need face_to_vertex_nindex() as one could use face_to_vertex()
50) ! and vertex_ids_nindex() to get the same result.
51) !gb: face_to_vertex_natural is required in GridLocalizeRegionsForUGrid() and needs
52) ! to be saved because:
53) ! (i) face_to_vertex() - Removes the duplicate faces, and the algorithm in
54) ! GridLocalizeRegionsForUGrid() assumes presence of ALL faces.
55) ! (ii) More importantly, face_to_vertex() eventually has vertex indices in
56) ! local index (different from natural index, when using multiple processors).
57) ! Note: A region in the inputfile will be described in terms of natural index.
58) PetscInt, pointer :: face_to_vertex_natural(:,:)
59) PetscInt, pointer :: face_to_vertex(:,:)
60) PetscInt, pointer :: cell_to_face_ghosted(:,:)
61) PetscInt, pointer :: vertex_ids_natural(:)
62) PetscInt, pointer :: cell_neighbors_local_ghosted(:,:) ! local neighbors
63) type(point_type), pointer :: vertices(:)
64) type(point_type), pointer :: face_centroid(:)
65) PetscReal, pointer :: face_area(:)
66) end type grid_unstructured_type
67)
68) type, public :: unstructured_explicit_type
69) PetscInt, pointer :: cell_ids(:)
70) PetscReal, pointer :: cell_volumes(:)
71) type(point3d_type), pointer :: cell_centroids(:)
72) PetscInt, pointer :: connections(:,:)
73) PetscReal, pointer :: face_areas(:)
74) type(point3d_type), pointer :: face_centroids(:)
75) PetscInt :: num_cells_global ! Number of cells in the entire domain
76) PetscInt :: num_elems
77) PetscInt :: num_elems_local ! Number of elements locally
78) PetscInt, pointer :: cell_connectivity(:,:)
79) type(point3d_type), pointer :: vertex_coordinates(:)
80) end type unstructured_explicit_type
81)
82) type, public :: unstructured_polyhedra_type
83) PetscInt, pointer :: cell_ids(:)
84) PetscInt, pointer :: cell_nfaces(:)
85) PetscInt, pointer :: cell_nverts(:)
86) PetscInt, pointer :: cell_faceids(:,:)
87) PetscInt, pointer :: cell_vertids(:,:)
88) PetscReal, pointer :: cell_volumes(:)
89) type(point3d_type), pointer :: cell_centroids(:)
90) PetscInt, pointer :: face_ids(:)
91) PetscInt, pointer :: face_cellids(:)
92) PetscInt, pointer :: face_nverts(:)
93) PetscInt, pointer :: face_vertids(:,:)
94) PetscReal, pointer :: face_areas(:)
95) type(point3d_type), pointer :: face_centroids(:)
96) type(point3d_type), pointer :: vertex_coordinates(:)
97) PetscInt :: num_cells_global
98) PetscInt :: num_cells_local
99) PetscInt :: num_faces_global
100) PetscInt :: num_faces_local
101) PetscInt :: num_vertices_global
102) PetscInt :: num_vertices_local
103) PetscInt :: max_nface_per_cell
104) PetscInt :: max_nvert_per_face
105) PetscInt :: max_nvert_per_cell
106) PetscInt :: num_ufaces_local
107) PetscInt :: num_ufaces_global
108) PetscInt :: num_verts_of_ufaces_local
109) PetscInt :: num_verts_of_ufaces_global
110) PetscInt, pointer :: uface_localids(:)
111) PetscInt, pointer :: uface_nverts(:)
112) PetscInt, pointer :: uface_natvertids(:)
113) PetscInt, pointer :: uface_left_natcellids(:)
114) PetscInt, pointer :: uface_right_natcellids(:)
115) PetscInt, pointer :: ugridf2pgridf(:)
116) PetscInt :: ugrid_num_faces_local
117) end type unstructured_polyhedra_type
118)
119) type, public :: ugdm_type
120) ! local: included both local (non-ghosted) and ghosted cells
121) ! global: includes only local (non-ghosted) cells
122) PetscInt :: ndof
123) ! for the below
124) ! ghosted = local (non-ghosted) and ghosted cells
125) ! local = local (non-ghosted) cells
126) IS :: is_ghosted_local ! IS for ghosted cells with local on-processor numbering
127) IS :: is_local_local ! IS for local cells with local on-processor numbering
128) IS :: is_ghosted_petsc ! IS for ghosted cells with petsc numbering
129) IS :: is_local_petsc ! IS for local cells with petsc numbering
130) IS :: is_ghosts_local ! IS for ghost cells with local on-processor numbering
131) IS :: is_ghosts_petsc ! IS for ghost cells with petsc numbering
132) IS :: is_local_natural ! IS for local cells with natural (global) numbering
133) VecScatter :: scatter_ltog ! scatter context for local to global updates
134) VecScatter :: scatter_gtol ! scatter context for global to local updates
135) VecScatter :: scatter_ltol ! scatter context for local to local updates
136) VecScatter :: scatter_gton ! scatter context for global to natural updates
137) ISLocalToGlobalMapping :: mapping_ltog ! petsc vec local to global mapping
138) !geh: deprecated in PETSc in spring 2014
139) ! ISLocalToGlobalMapping :: mapping_ltogb ! block form of mapping_ltog
140) Vec :: global_vec ! global vec (no ghost cells), petsc-ordering
141) Vec :: local_vec ! local vec (includes local and ghosted cells), local ordering
142) VecScatter :: scatter_bet_grids ! scatter context between surface and subsurface
143) ! grids
144) VecScatter :: scatter_bet_grids_1dof ! scatter context between surface and
145) ! subsurface grids for 1-DOF
146) VecScatter :: scatter_bet_grids_ndof ! scatter context between surface and
147) ! subsurface grids for N-DOFs
148) AO :: ao_natural_to_petsc
149) end type ugdm_type
150)
151) ! PetscInt, parameter :: HEX_TYPE = 1
152) ! PetscInt, parameter :: TET_TYPE = 2
153) ! PetscInt, parameter :: WEDGE_TYPE = 3
154) ! PetscInt, parameter :: PYR_TYPE = 4
155) ! PetscInt, parameter :: TRI_FACE_TYPE = 1
156) ! PetscInt, parameter :: QUAD_FACE_TYPE = 2
157) ! PetscInt, parameter :: MAX_VERT_PER_FACE = 4
158)
159) public :: UGridCreate, &
160) UGridExplicitCreate, &
161) UGridPolyhedraCreate, &
162) UGridMapIndices, &
163) UGridDMCreateJacobian, &
164) UGridDMCreateVector, &
165) UGridDestroy, &
166) UGridCreateUGDM, &
167) UGridCreateUGDMShell, &
168) UGridDMDestroy, &
169) UGridPartition, &
170) UGridNaturalToPetsc, &
171) UGridCreateOldVec, &
172) UGridExplicitDestroy
173)
174) contains
175)
176) ! ************************************************************************** !
177)
178) function UGDMCreate()
179) !
180) ! Creates an unstructured grid distributed mesh object
181) !
182) ! Author: Glenn Hammond
183) ! Date: 10/21/09
184) !
185)
186) implicit none
187)
188) type(ugdm_type), pointer :: UGDMCreate
189)
190) type(ugdm_type), pointer :: ugdm
191)
192) allocate(ugdm)
193) ugdm%is_ghosted_local = 0
194) ugdm%is_local_local = 0
195) ugdm%is_ghosted_petsc = 0
196) ugdm%is_local_petsc = 0
197) ugdm%is_ghosts_local = 0
198) ugdm%is_ghosts_petsc = 0
199) ugdm%is_local_natural = 0
200) ugdm%scatter_ltog = 0
201) ugdm%scatter_gtol = 0
202) ugdm%scatter_ltol = 0
203) ugdm%scatter_gton = 0
204) ugdm%mapping_ltog = 0
205) ugdm%global_vec = 0
206) ugdm%local_vec = 0
207) ugdm%scatter_bet_grids = 0
208) ugdm%scatter_bet_grids_1dof = 0
209) ugdm%scatter_bet_grids_ndof = 0
210) ugdm%ao_natural_to_petsc = 0 ! this is solely a pointer, do not destroy
211) UGDMCreate => ugdm
212)
213) end function UGDMCreate
214)
215) ! ************************************************************************** !
216)
217) function UGridCreate()
218) !
219) ! Creates an unstructured grid object
220) !
221) ! Author: Glenn Hammond
222) ! Date: 09/30/09
223) !
224)
225) implicit none
226)
227) type(grid_unstructured_type), pointer :: UGridCreate
228)
229) type(grid_unstructured_type), pointer :: unstructured_grid
230)
231) allocate(unstructured_grid)
232)
233) ! variables for all unstructured grids
234) unstructured_grid%num_ghost_cells = 0
235) unstructured_grid%global_offset = 0
236) unstructured_grid%nmax = 0
237) unstructured_grid%nlmax = 0
238) unstructured_grid%ngmax = 0
239) nullify(unstructured_grid%hash)
240) unstructured_grid%num_hash = 100
241) nullify(unstructured_grid%cell_ids_natural)
242) nullify(unstructured_grid%cell_ids_petsc)
243) nullify(unstructured_grid%ghost_cell_ids_petsc)
244) unstructured_grid%ao_natural_to_petsc = 0
245) nullify(unstructured_grid%explicit_grid)
246) nullify(unstructured_grid%polyhedra_grid)
247)
248) ! variables for implicit unstructured grids
249) unstructured_grid%grid_type = THREE_DIM_GRID
250) unstructured_grid%num_vertices_global = 0
251) unstructured_grid%num_vertices_local = 0
252) unstructured_grid%max_ndual_per_cell = 0
253) unstructured_grid%max_nvert_per_cell = 0
254) unstructured_grid%max_cells_sharing_a_vertex = 24
255) nullify(unstructured_grid%cell_type)
256) nullify(unstructured_grid%cell_vertices)
257) nullify(unstructured_grid%face_to_cell_ghosted)
258) nullify(unstructured_grid%face_to_vertex_natural)
259) nullify(unstructured_grid%face_to_vertex)
260) nullify(unstructured_grid%cell_to_face_ghosted)
261) nullify(unstructured_grid%vertex_ids_natural)
262) nullify(unstructured_grid%vertices)
263) nullify(unstructured_grid%cell_neighbors_local_ghosted)
264) nullify(unstructured_grid%connection_to_face)
265) nullify(unstructured_grid%face_centroid)
266) nullify(unstructured_grid%face_area)
267)
268) UGridCreate => unstructured_grid
269)
270) end function UGridCreate
271)
272) ! ************************************************************************** !
273)
274) function UGridExplicitCreate()
275) !
276) ! Creates an explicit unstructured grid object
277) !
278) ! Author: Glenn Hammond
279) ! Date: 05/14/12
280) !
281)
282) implicit none
283)
284) type(unstructured_explicit_type), pointer :: UGridExplicitCreate
285)
286) type(unstructured_explicit_type), pointer :: explicit_grid
287)
288) allocate(explicit_grid)
289)
290) nullify(explicit_grid%cell_ids)
291) nullify(explicit_grid%cell_volumes)
292) nullify(explicit_grid%cell_centroids)
293) nullify(explicit_grid%connections)
294) nullify(explicit_grid%face_areas)
295) nullify(explicit_grid%face_centroids)
296) nullify(explicit_grid%cell_connectivity)
297) nullify(explicit_grid%vertex_coordinates)
298)
299) explicit_grid%num_cells_global = 0
300) explicit_grid%num_elems = 0
301) explicit_grid%num_elems_local = 0
302)
303) UGridExplicitCreate => explicit_grid
304)
305) end function UGridExplicitCreate
306)
307) ! ************************************************************************** !
308)
309) function UGridPolyhedraCreate()
310) !
311) ! Creates a polyhedra unstructured grid object.
312) !
313) ! Author: Gautam Bisht, LBL
314) ! Date: 09/29/13
315) !
316)
317) implicit none
318)
319) type(unstructured_polyhedra_type), pointer :: UGridPolyhedraCreate
320)
321) type(unstructured_polyhedra_type), pointer :: polyhedra_grid
322)
323) allocate(polyhedra_grid)
324)
325) polyhedra_grid%num_cells_global = 0
326) polyhedra_grid%num_cells_local = 0
327) polyhedra_grid%num_faces_global = 0
328) polyhedra_grid%num_faces_local = 0
329) polyhedra_grid%num_vertices_global = 0
330) polyhedra_grid%num_vertices_local = 0
331) polyhedra_grid%max_nface_per_cell = 0
332) polyhedra_grid%max_nvert_per_face = 0
333) polyhedra_grid%max_nvert_per_cell = 0
334) polyhedra_grid%num_ufaces_local = 0
335) polyhedra_grid%num_ufaces_global = 0
336) polyhedra_grid%num_verts_of_ufaces_local = 0
337) polyhedra_grid%num_verts_of_ufaces_global = 0
338) polyhedra_grid%ugrid_num_faces_local = 0
339)
340) nullify(polyhedra_grid%cell_ids)
341) nullify(polyhedra_grid%cell_nfaces)
342) nullify(polyhedra_grid%cell_vertids)
343) nullify(polyhedra_grid%cell_faceids)
344) nullify(polyhedra_grid%cell_nverts)
345) nullify(polyhedra_grid%cell_volumes)
346) nullify(polyhedra_grid%cell_centroids)
347) nullify(polyhedra_grid%face_ids)
348) nullify(polyhedra_grid%face_cellids)
349) nullify(polyhedra_grid%face_nverts)
350) nullify(polyhedra_grid%face_vertids)
351) nullify(polyhedra_grid%face_areas)
352) nullify(polyhedra_grid%face_centroids)
353) nullify(polyhedra_grid%vertex_coordinates)
354) nullify(polyhedra_grid%uface_localids)
355) nullify(polyhedra_grid%uface_nverts)
356) nullify(polyhedra_grid%uface_natvertids)
357) nullify(polyhedra_grid%uface_left_natcellids)
358) nullify(polyhedra_grid%uface_right_natcellids)
359) nullify(polyhedra_grid%ugridf2pgridf)
360)
361) UGridPolyhedraCreate => polyhedra_grid
362)
363) end function UGridPolyhedraCreate
364)
365) ! ************************************************************************** !
366)
367) subroutine UGridCreateUGDM(unstructured_grid,ugdm,ndof,option)
368) !
369) ! Constructs mappings / scatter contexts for PETSc DM
370) ! object
371) !
372) ! Author: Glenn Hammond
373) ! Date: 09/30/09
374) !
375)
376) use Option_module
377) use Utility_module, only: reallocateIntArray
378)
379) implicit none
380)
381) #include "petsc/finclude/petscvec.h"
382) #include "petsc/finclude/petscvec.h90"
383) #include "petsc/finclude/petscmat.h"
384) #include "petsc/finclude/petscmat.h90"
385) #include "petsc/finclude/petscdm.h"
386) #include "petsc/finclude/petscdm.h90"
387) #include "petsc/finclude/petscis.h"
388) #include "petsc/finclude/petscis.h90"
389) #include "petsc/finclude/petscviewer.h"
390)
391) type(grid_unstructured_type) :: unstructured_grid
392) type(ugdm_type), pointer :: ugdm
393) PetscInt :: ndof
394) type(option_type) :: option
395)
396) PetscInt, pointer :: int_ptr(:)
397) PetscInt :: local_id, ghosted_id
398) PetscInt :: idof
399) IS :: is_tmp
400) Vec :: vec_tmp
401) PetscErrorCode :: ierr
402) character(len=MAXWORDLENGTH) :: ndof_word
403) character(len=MAXSTRINGLENGTH) :: string
404)
405) PetscViewer :: viewer
406)
407) PetscInt, allocatable :: int_array(:)
408)
409) ugdm => UGDMCreate()
410) ugdm%ndof = ndof
411)
412) #if UGRID_DEBUG
413) write(ndof_word,*) ndof
414) ndof_word = adjustl(ndof_word)
415) ndof_word = '_' // trim(ndof_word)
416) string = 'Vectors' // ndof_word
417) call printMsg(option,string)
418) #endif
419)
420) ! create global vec
421) !call VecCreateMPI(option%mycomm,unstructured_grid%nlmax*ndof, &
422) ! PETSC_DETERMINE,ugdm%global_vec,ierr)
423) call VecCreate(option%mycomm,ugdm%global_vec,ierr);CHKERRQ(ierr)
424) call VecSetSizes(ugdm%global_vec,unstructured_grid%nlmax*ndof, &
425) PETSC_DECIDE,ierr);CHKERRQ(ierr)
426) call VecSetBlockSize(ugdm%global_vec,ndof,ierr);CHKERRQ(ierr)
427) call VecSetFromOptions(ugdm%global_vec,ierr);CHKERRQ(ierr)
428)
429) ! create local vec
430) !call VecCreateSeq(PETSC_COMM_SELF,unstructured_grid%ngmax*ndof, &
431) ! ugdm%local_vec,ierr)
432) call VecCreate(PETSC_COMM_SELF,ugdm%local_vec,ierr);CHKERRQ(ierr)
433) call VecSetSizes(ugdm%local_vec,unstructured_grid%ngmax*ndof,PETSC_DECIDE, &
434) ierr);CHKERRQ(ierr)
435) call VecSetBlockSize(ugdm%local_vec,ndof,ierr);CHKERRQ(ierr)
436) call VecSetFromOptions(ugdm%local_vec,ierr);CHKERRQ(ierr)
437)
438) ! IS for global numbering of local, non-ghosted cells
439) !geh call VecGetOwnershipRange(ugdm%global_vec,istart,iend,ierr)
440) ! ISCreateBlock requires block ids, not indices. Therefore, istart should be
441) ! the offset of the block from the beginning of the vector.
442) !geh istart = istart / ndof
443) allocate(int_array(unstructured_grid%nlmax))
444) do local_id = 1, unstructured_grid%nlmax
445) !geh int_array(local_id) = (local_id-1)+istart
446) int_array(local_id) = (local_id-1) + unstructured_grid%global_offset
447) enddo
448)
449) ! arguments for ISCreateBlock():
450) ! option%mycomm - the MPI communicator
451) ! ndof - number of elements in each block
452) ! unstructured_grid%nlmax - the length of the index set
453) ! (the number of blocks
454) ! int_array - the list of integers, one for each block and count
455) ! of block not indices
456) ! PETSC_COPY_VALUES - see PetscCopyMode, only PETSC_COPY_VALUES and
457) ! PETSC_OWN_POINTER are supported in this routine
458) ! ugdm%is_local_petsc - the new index set
459) ! ierr - PETScErrorCode
460) call ISCreateBlock(option%mycomm,ndof,unstructured_grid%nlmax, &
461) int_array,PETSC_COPY_VALUES,ugdm%is_local_petsc, &
462) ierr);CHKERRQ(ierr)
463) deallocate(int_array)
464)
465) #if UGRID_DEBUG
466) string = 'is_local_petsc' // trim(ndof_word) // '.out'
467) call PetscViewerASCIIOpen(option%mycomm,trim(string),viewer, &
468) ierr);CHKERRQ(ierr)
469) call ISView(ugdm%is_local_petsc,viewer,ierr);CHKERRQ(ierr)
470) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
471) #endif
472)
473) ! IS for local numbering of ghosts cells
474) allocate(int_array(unstructured_grid%num_ghost_cells))
475) do ghosted_id = 1, unstructured_grid%num_ghost_cells
476) int_array(ghosted_id) = (ghosted_id+unstructured_grid%nlmax-1)
477) enddo
478) call ISCreateBlock(option%mycomm,ndof,unstructured_grid%num_ghost_cells, &
479) int_array,PETSC_COPY_VALUES,ugdm%is_ghosts_local, &
480) ierr);CHKERRQ(ierr)
481) deallocate(int_array)
482)
483) #if UGRID_DEBUG
484) string = 'is_ghosts_local' // trim(ndof_word) // '.out'
485) call PetscViewerASCIIOpen(option%mycomm,trim(string),viewer, &
486) ierr);CHKERRQ(ierr)
487) call ISView(ugdm%is_ghosts_local,viewer,ierr);CHKERRQ(ierr)
488) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
489) #endif
490)
491) #if UGRID_DEBUG
492) string = 'Index Sets' // ndof_word
493) call printMsg(option,'Index Sets')
494) #endif
495)
496) ! IS for local numbering of ghosts cells
497) allocate(int_array(unstructured_grid%num_ghost_cells))
498) do ghosted_id = 1, unstructured_grid%num_ghost_cells
499) int_array(ghosted_id) = &
500) (unstructured_grid%ghost_cell_ids_petsc(ghosted_id)-1)
501) enddo
502) call ISCreateBlock(option%mycomm,ndof,unstructured_grid%num_ghost_cells, &
503) int_array,PETSC_COPY_VALUES,ugdm%is_ghosts_petsc, &
504) ierr);CHKERRQ(ierr)
505) deallocate(int_array)
506)
507) #if UGRID_DEBUG
508) string = 'is_ghosts_petsc' // trim(ndof_word) // '.out'
509) call PetscViewerASCIIOpen(option%mycomm,trim(string),viewer, &
510) ierr);CHKERRQ(ierr)
511) call ISView(ugdm%is_ghosts_petsc,viewer,ierr);CHKERRQ(ierr)
512) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
513) #endif
514)
515) ! IS for local numbering of local, non-ghosted cells
516) allocate(int_array(unstructured_grid%nlmax))
517) do local_id = 1, unstructured_grid%nlmax
518) int_array(local_id) = (local_id-1)
519) enddo
520) call ISCreateBlock(option%mycomm,ndof,unstructured_grid%nlmax, &
521) int_array,PETSC_COPY_VALUES,ugdm%is_local_local, &
522) ierr);CHKERRQ(ierr)
523) deallocate(int_array)
524)
525) #if UGRID_DEBUG
526) string = 'is_local_local' // trim(ndof_word) // '.out'
527) call PetscViewerASCIIOpen(option%mycomm,trim(string),viewer, &
528) ierr);CHKERRQ(ierr)
529) call ISView(ugdm%is_local_local,viewer,ierr);CHKERRQ(ierr)
530) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
531) #endif
532)
533) ! IS for ghosted numbering of local ghosted cells
534) allocate(int_array(unstructured_grid%ngmax))
535) do ghosted_id = 1, unstructured_grid%ngmax
536) int_array(ghosted_id) = (ghosted_id-1)
537) enddo
538) call ISCreateBlock(option%mycomm,ndof,unstructured_grid%ngmax, &
539) int_array,PETSC_COPY_VALUES,ugdm%is_ghosted_local, &
540) ierr);CHKERRQ(ierr)
541) deallocate(int_array)
542)
543) #if UGRID_DEBUG
544) string = 'is_ghosted_local' // trim(ndof_word) // '.out'
545) call PetscViewerASCIIOpen(option%mycomm,trim(string),viewer, &
546) ierr);CHKERRQ(ierr)
547) call ISView(ugdm%is_ghosted_local,viewer,ierr);CHKERRQ(ierr)
548) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
549) #endif
550)
551) ! IS for petsc numbering of local ghosted cells
552) allocate(int_array(unstructured_grid%ngmax))
553) do local_id = 1, unstructured_grid%nlmax
554) !geh int_array(local_id) = istart+(local_id-1)
555) int_array(local_id) = (local_id-1) + unstructured_grid%global_offset
556) enddo
557) do ghosted_id = 1,unstructured_grid%num_ghost_cells
558) int_array(unstructured_grid%nlmax+ghosted_id) = &
559) (unstructured_grid%ghost_cell_ids_petsc(ghosted_id)-1)
560) enddo
561) call ISCreateBlock(option%mycomm,ndof,unstructured_grid%ngmax, &
562) int_array,PETSC_COPY_VALUES,ugdm%is_ghosted_petsc, &
563) ierr);CHKERRQ(ierr)
564) deallocate(int_array)
565)
566) #if UGRID_DEBUG
567) string = 'is_ghosted_petsc' // trim(ndof_word) // '.out'
568) call PetscViewerASCIIOpen(option%mycomm,trim(string),viewer, &
569) ierr);CHKERRQ(ierr)
570) call ISView(ugdm%is_ghosted_petsc,viewer,ierr);CHKERRQ(ierr)
571) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
572) #endif
573)
574) ! create a local to global mapping
575) #if UGRID_DEBUG
576) string = 'ISLocalToGlobalMapping' // ndof_word
577) call printMsg(option,string)
578) #endif
579)
580) call ISLocalToGlobalMappingCreateIS(ugdm%is_ghosted_petsc, &
581) ugdm%mapping_ltog,ierr);CHKERRQ(ierr)
582)
583) #if UGRID_DEBUG
584) string = 'mapping_ltog' // trim(ndof_word) // '.out'
585) call PetscViewerASCIIOpen(option%mycomm,trim(string),viewer, &
586) ierr);CHKERRQ(ierr)
587) call ISLocalToGlobalMappingView(ugdm%mapping_ltog,viewer,ierr);CHKERRQ(ierr)
588) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
589) #endif
590)
591) #if UGRID_DEBUG
592) string = 'local to global' // ndof_word
593) call printMsg(option,string)
594) #endif
595)
596) ! Create local to global scatter
597) call VecScatterCreate(ugdm%local_vec,ugdm%is_local_local,ugdm%global_vec, &
598) ugdm%is_local_petsc,ugdm%scatter_ltog, &
599) ierr);CHKERRQ(ierr)
600)
601) #if UGRID_DEBUG
602) string = 'scatter_ltog' // trim(ndof_word) // '.out'
603) call PetscViewerASCIIOpen(option%mycomm,trim(string),viewer, &
604) ierr);CHKERRQ(ierr)
605) call VecScatterView(ugdm%scatter_ltog,viewer,ierr);CHKERRQ(ierr)
606) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
607) #endif
608)
609) #if UGRID_DEBUG
610) string = 'global to local' // ndof_word
611) call printMsg(option,string)
612) #endif
613)
614) ! Create global to local scatter
615) call VecScatterCreate(ugdm%global_vec,ugdm%is_ghosted_petsc,ugdm%local_vec, &
616) ugdm%is_ghosted_local,ugdm%scatter_gtol, &
617) ierr);CHKERRQ(ierr)
618)
619) #if UGRID_DEBUG
620) string = 'scatter_gtol' // trim(ndof_word) // '.out'
621) call PetscViewerASCIIOpen(option%mycomm,trim(string),viewer, &
622) ierr);CHKERRQ(ierr)
623) call VecScatterView(ugdm%scatter_gtol,viewer,ierr);CHKERRQ(ierr)
624) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
625) #endif
626)
627) #if UGRID_DEBUG
628) string = 'local to local' // ndof_word
629) call printMsg(option,string)
630) #endif
631)
632) ! Create local to local scatter. Essentially remap the global to local as
633) ! PETSc does in daltol.c
634) call VecScatterCopy(ugdm%scatter_gtol,ugdm%scatter_ltol,ierr);CHKERRQ(ierr)
635) call ISGetIndicesF90(ugdm%is_local_local,int_ptr,ierr);CHKERRQ(ierr)
636) call VecScatterRemap(ugdm%scatter_ltol,int_ptr,PETSC_NULL_INTEGER, &
637) ierr);CHKERRQ(ierr)
638) call ISRestoreIndicesF90(ugdm%is_local_local,int_ptr,ierr);CHKERRQ(ierr)
639)
640) #if UGRID_DEBUG
641) string = 'scatter_ltol' // trim(ndof_word) // '.out'
642) call PetscViewerASCIIOpen(option%mycomm,trim(string),viewer, &
643) ierr);CHKERRQ(ierr)
644) call VecScatterView(ugdm%scatter_ltol,viewer,ierr);CHKERRQ(ierr)
645) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
646) #endif
647)
648) ! Set up global to natural scatter
649) ! Create index set of local non-ghosted Petsc ordering
650) call VecCreateMPI(option%mycomm,unstructured_grid%nlmax, &
651) PETSC_DETERMINE,vec_tmp,ierr);CHKERRQ(ierr)
652) !geh call VecGetOwnershipRange(vec_tmp,istart,iend,ierr)
653) call VecDestroy(vec_tmp,ierr);CHKERRQ(ierr)
654) allocate(int_array(unstructured_grid%nlmax))
655) do local_id = 1, unstructured_grid%nlmax
656) !geh int_array(local_id) = (local_id-1)+istart
657) int_array(local_id) = (local_id-1) + unstructured_grid%global_offset
658) enddo
659) call ISCreateGeneral(option%mycomm,unstructured_grid%nlmax, &
660) int_array,PETSC_COPY_VALUES,is_tmp,ierr);CHKERRQ(ierr)
661) deallocate(int_array)
662) call AOPetscToApplicationIS(unstructured_grid%ao_natural_to_petsc, &
663) is_tmp,ierr);CHKERRQ(ierr)
664) ! remap for ndof > 1 !geh: no longer need to accommodate ndof > 1, but leave
665) ! alone for now.
666) allocate(int_array(unstructured_grid%nlmax))
667) call ISGetIndicesF90(is_tmp,int_ptr,ierr);CHKERRQ(ierr)
668) do local_id = 1, unstructured_grid%nlmax
669) int_array(local_id) = int_ptr(local_id)
670) enddo
671) call ISRestoreIndicesF90(is_tmp,int_ptr,ierr);CHKERRQ(ierr)
672) call ISDestroy(is_tmp,ierr);CHKERRQ(ierr)
673) call ISCreateBlock(option%mycomm,ndof,unstructured_grid%nlmax, &
674) int_array,PETSC_COPY_VALUES,ugdm%is_local_natural, &
675) ierr);CHKERRQ(ierr)
676) deallocate(int_array)
677)
678) #if UGRID_DEBUG
679) string = 'is_local_natural' // trim(ndof_word) // '.out'
680) call PetscViewerASCIIOpen(option%mycomm,trim(string),viewer, &
681) ierr);CHKERRQ(ierr)
682) call ISView(ugdm%is_local_natural,viewer,ierr);CHKERRQ(ierr)
683) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
684) #endif
685)
686) !call VecCreateMPI(option%mycomm,unstructured_grid%nlmax*ndof, &
687) ! PETSC_DETERMINE,vec_tmp,ierr)
688) call VecCreate(option%mycomm,vec_tmp,ierr);CHKERRQ(ierr)
689) call VecSetSizes(vec_tmp,unstructured_grid%nlmax*ndof,PETSC_DECIDE, &
690) ierr);CHKERRQ(ierr)
691) call VecSetBlockSize(vec_tmp,ndof,ierr);CHKERRQ(ierr)
692) call VecSetFromOptions(vec_tmp,ierr);CHKERRQ(ierr)
693) call VecScatterCreate(ugdm%global_vec,ugdm%is_local_petsc,vec_tmp, &
694) ugdm%is_local_natural,ugdm%scatter_gton, &
695) ierr);CHKERRQ(ierr)
696) call VecDestroy(vec_tmp,ierr);CHKERRQ(ierr)
697)
698) #if UGRID_DEBUG
699) string = 'scatter_gton' // trim(ndof_word) // '.out'
700) call PetscViewerASCIIOpen(option%mycomm,trim(string),viewer, &
701) ierr);CHKERRQ(ierr)
702) call VecScatterView(ugdm%scatter_gton,viewer,ierr);CHKERRQ(ierr)
703) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
704)
705) #endif
706)
707) ! set the ao_natural_to_petsc pointer
708) ugdm%ao_natural_to_petsc = unstructured_grid%ao_natural_to_petsc
709)
710) end subroutine UGridCreateUGDM
711)
712) ! ************************************************************************** !
713)
714) subroutine UGridCreateUGDMShell(unstructured_grid,da,ugdm,ndof,option)
715)
716) !
717) ! Sets up PETSc DM Shell for unstructured grid
718) !
719) ! Author: Gautam Bisht, LBNL
720) ! Date: 11/10/15
721) !
722)
723) use Option_module
724) use Utility_module, only: reallocateIntArray
725)
726) implicit none
727)
728) #include "petsc/finclude/petscdm.h90"
729) #include "petsc/finclude/petscdmda.h"
730)
731) type(grid_unstructured_type) :: unstructured_grid
732) DM :: da
733) type(ugdm_type), pointer :: ugdm
734) PetscInt :: ndof
735) type(option_type) :: option
736)
737) Vec :: global_vec, local_vec
738) !Mat :: jac
739) PetscErrorCode :: ierr
740)
741) ! Create UGDM
742) call UGridCreateUGDM(unstructured_grid,ugdm,ndof,option)
743)
744) ! Create the DMShell
745) call DMShellCreate(option%mycomm,da,ierr);CHKERRQ(ierr)
746)
747) ! Set VecScatters
748) call DMShellSetGlobalToLocalVecScatter(da,ugdm%scatter_gtol,ierr);CHKERRQ(ierr)
749) call DMShellSetLocalToGlobalVecScatter(da,ugdm%scatter_ltog,ierr);CHKERRQ(ierr)
750) call DMShellSetLocalToLocalVecScatter(da,ugdm%scatter_ltol,ierr);CHKERRQ(ierr)
751)
752) ! Create vectors
753) call UGridDMCreateVector(unstructured_grid,ugdm,global_vec,GLOBAL,option)
754) call UGridDMCreateVector(unstructured_grid,ugdm,local_vec,LOCAL,option)
755)
756) ! Set vectors
757) call DMShellSetGlobalVector(da,global_vec,ierr);CHKERRQ(ierr)
758) call DMShellSetLocalVector(da,local_vec,ierr);CHKERRQ(ierr)
759)
760) call VecDestroy(global_vec,ierr);CHKERRQ(ierr)
761) call VecDestroy(local_vec,ierr);CHKERRQ(ierr)
762)
763) ! GB: Can mat_type passed as an argument to this subroutine?
764) !call UGridDMCreateJacobian(unstructured_grid,ugdm,mat_type,jac,option)
765) !call DMShellSetMatrix(J,ierr);CHKERRQ(ierr)
766) !call MatDestroy(J,ierr);CHKERRQ(ierr)
767)
768) end subroutine UGridCreateUGDMShell
769)
770) ! ************************************************************************** !
771)
772) subroutine UGridDMCreateJacobian(unstructured_grid,ugdm,mat_type,J,option)
773) !
774) ! Creates a Jacobian matrix based on the unstructured
775) ! grid dual
776) !
777) ! Author: Glenn Hammond
778) ! Date: 11/05/09
779) !
780)
781) use Option_module
782)
783) implicit none
784)
785) type(grid_unstructured_type) :: unstructured_grid
786) type(ugdm_type) :: ugdm
787) MatType :: mat_type
788) Mat :: J
789) type(option_type) :: option
790)
791) PetscInt, allocatable :: d_nnz(:), o_nnz(:)
792) PetscInt :: local_id, ineighbor, neighbor_id
793) PetscInt :: iconn, id_up, id_dn
794) PetscInt :: ndof_local
795) PetscErrorCode :: ierr
796)
797) allocate(d_nnz(unstructured_grid%nlmax))
798) allocate(o_nnz(unstructured_grid%nlmax))
799) d_nnz = 1 ! start 1 since diagonal connection to self
800) o_nnz = 0
801) if (associated(unstructured_grid%explicit_grid)) then
802) do iconn = 1, size(unstructured_grid%explicit_grid%connections,2)
803) id_up = unstructured_grid%explicit_grid%connections(1,iconn)
804) id_dn = unstructured_grid%explicit_grid%connections(2,iconn)
805) if (id_up <= unstructured_grid%nlmax) then ! local
806) if (id_dn <= unstructured_grid%nlmax) then
807) d_nnz(id_up) = d_nnz(id_up) + 1
808) else
809) o_nnz(id_up) = o_nnz(id_up) + 1
810) endif
811) endif
812) if (id_dn <= unstructured_grid%nlmax) then ! local
813) if (id_up <= unstructured_grid%nlmax) then
814) d_nnz(id_dn) = d_nnz(id_dn) + 1
815) else
816) o_nnz(id_dn) = o_nnz(id_dn) + 1
817) endif
818) endif
819) enddo
820) else
821) do local_id = 1, unstructured_grid%nlmax
822) do ineighbor = 1, unstructured_grid% &
823) cell_neighbors_local_ghosted(0,local_id)
824) neighbor_id = unstructured_grid% &
825) cell_neighbors_local_ghosted(ineighbor,local_id)
826) if (neighbor_id > 0) then
827) d_nnz(local_id) = d_nnz(local_id) + 1
828) else
829) o_nnz(local_id) = o_nnz(local_id) + 1
830) endif
831) enddo
832) enddo
833) endif
834)
835) ndof_local = unstructured_grid%nlmax*ugdm%ndof
836)
837) call MatCreate(option%mycomm,J,ierr); CHKERRQ(ierr)
838) call MatSetType(J, mat_type, ierr); CHKERRQ(ierr)
839) call MatSetSizes(J,ndof_local,ndof_local,PETSC_DETERMINE,PETSC_DETERMINE, &
840) ierr) ;CHKERRQ(ierr)
841) call MatXAIJSetPreallocation(J,ugdm%ndof,d_nnz,o_nnz, &
842) PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
843) ierr); CHKERRQ(ierr)
844)
845) call MatSetLocalToGlobalMapping(J,ugdm%mapping_ltog, &
846) ugdm%mapping_ltog,ierr);CHKERRQ(ierr)
847)
848) deallocate(d_nnz)
849) deallocate(o_nnz)
850)
851) end subroutine UGridDMCreateJacobian
852)
853) ! ************************************************************************** !
854)
855) subroutine UGridDMCreateVector(unstructured_grid,ugdm,vec,vec_type,option)
856) !
857) ! Creates a global vector with PETSc ordering
858) !
859) ! Author: Glenn Hammond
860) ! Date: 11/06/09
861) !
862)
863) use Option_module
864)
865) implicit none
866)
867) type(grid_unstructured_type) :: unstructured_grid
868) type(ugdm_type) :: ugdm
869) Vec :: vec
870) PetscInt :: vec_type
871) type(option_type) :: option
872)
873) PetscErrorCode :: ierr
874)
875) select case(vec_type)
876) case(GLOBAL)
877) !call VecCreateMPI(option%mycomm,unstructured_grid%nlmax* &
878) ! ugdm%ndof, &
879) ! PETSC_DETERMINE,vec,ierr)
880) call VecCreate(option%mycomm,vec,ierr);CHKERRQ(ierr)
881) call VecSetSizes(vec,unstructured_grid%nlmax*ugdm%ndof, &
882) PETSC_DECIDE,ierr);CHKERRQ(ierr)
883) call VecSetLocalToGlobalMapping(vec,ugdm%mapping_ltog, &
884) ierr);CHKERRQ(ierr)
885) call VecSetBlockSize(vec,ugdm%ndof,ierr);CHKERRQ(ierr)
886) call VecSetFromOptions(vec,ierr);CHKERRQ(ierr)
887) case(LOCAL)
888) !call VecCreateSeq(PETSC_COMM_SELF,unstructured_grid%ngmax* &
889) ! ugdm%ndof, &
890) ! vec,ierr)
891) call VecCreate(PETSC_COMM_SELF,vec,ierr);CHKERRQ(ierr)
892) call VecSetSizes(vec,unstructured_grid%ngmax*ugdm%ndof, &
893) PETSC_DECIDE,ierr);CHKERRQ(ierr)
894) call VecSetBlockSize(vec,ugdm%ndof,ierr);CHKERRQ(ierr)
895) call VecSetFromOptions(vec,ierr);CHKERRQ(ierr)
896) case(NATURAL)
897) !call VecCreateMPI(option%mycomm,unstructured_grid%nlmax* &
898) ! ugdm%ndof, &
899) ! PETSC_DETERMINE,vec,ierr)
900) call VecCreate(option%mycomm,vec,ierr);CHKERRQ(ierr)
901) call VecSetSizes(vec,unstructured_grid%nlmax*ugdm%ndof, &
902) PETSC_DECIDE,ierr);CHKERRQ(ierr)
903) call VecSetBlockSize(vec,ugdm%ndof,ierr);CHKERRQ(ierr)
904) call VecSetFromOptions(vec,ierr);CHKERRQ(ierr)
905) end select
906)
907) end subroutine UGridDMCreateVector
908)
909) ! ************************************************************************** !
910)
911) subroutine UGridMapIndices(unstructured_grid,ugdm,nG2L,nL2G,nG2A,option)
912) !
913) ! maps global, local and natural indices of cells to each other
914) !
915) ! Author: Glenn Hammond
916) ! Date: 11/06/09
917) !
918)
919) use Option_module
920)
921) implicit none
922)
923) type(grid_unstructured_type) :: unstructured_grid
924) type(ugdm_type) :: ugdm
925) PetscInt, pointer :: nG2L(:)
926) PetscInt, pointer :: nL2G(:)
927) PetscInt, pointer :: nG2A(:)
928) type(option_type) :: option
929)
930) PetscErrorCode :: ierr
931) PetscInt, pointer :: int_ptr(:)
932) PetscInt :: local_id
933) PetscInt :: ghosted_id
934)
935) allocate(nG2L(unstructured_grid%ngmax))
936) allocate(nL2G(unstructured_grid%nlmax))
937) allocate(nG2A(unstructured_grid%ngmax))
938)
939) ! initialize ghosted to 0
940) !geh: any index beyond %nlmax will be 0 indicating that there is no local
941) ! counterpart (i.e., it is a ghost cell)
942) nG2L = 0
943)
944) !geh: Yes, it seems redundant that that we are setting both nL2G and nG2L to
945) ! the same index, but keep in mind that nG2L extends beyond %nlmax and
946) ! we need these arrays to provide seemless integration for structured and
947) ! unstructured
948) do local_id = 1, unstructured_grid%nlmax
949) nL2G(local_id) = local_id
950) nG2L(local_id) = local_id
951) enddo
952)
953) call ISGetIndicesF90(ugdm%is_ghosted_petsc,int_ptr,ierr);CHKERRQ(ierr)
954) do ghosted_id = 1, unstructured_grid%ngmax
955) nG2A(ghosted_id) = int_ptr(ghosted_id)+1
956) enddo
957) call ISRestoreIndicesF90(ugdm%is_ghosted_petsc,int_ptr,ierr);CHKERRQ(ierr)
958) nG2A = nG2A - 1
959) call AOPetscToApplication(unstructured_grid%ao_natural_to_petsc, &
960) unstructured_grid%ngmax, &
961) nG2A,ierr);CHKERRQ(ierr)
962) nG2A = nG2A + 1 ! 1-based
963)
964)
965) end subroutine UGridMapIndices
966)
967) ! ************************************************************************** !
968)
969) subroutine UGridPartition(ugrid,option,Dual_mat,is_new, &
970) num_cells_local_new)
971) !
972) ! UGridGet_Dual_Part_IS: Given an adjacency matrix, calculates the dual
973) ! partitions, and provides a new IS with the ids
974) ! of the local cells on the processor
975) !
976) ! Author: Glenn Hammond
977) ! Date: 10/05/12
978) !
979)
980) use Option_module
981)
982) implicit none
983)
984) #include "petsc/finclude/petscmat.h"
985) #include "petsc/finclude/petscmat.h90"
986) #include "petsc/finclude/petscis.h"
987) #include "petsc/finclude/petscis.h90"
988) #include "petsc/finclude/petscviewer.h"
989)
990) type(grid_unstructured_type) :: ugrid
991) type(option_type) :: option
992) Mat :: Dual_mat
993) IS :: is_new
994) PetscInt :: num_cells_local_new
995)
996) MatPartitioning :: Part
997) PetscInt, allocatable :: cell_counts(:)
998) PetscInt :: iflag
999) PetscViewer :: viewer
1000) PetscInt :: local_vertex_offset
1001) PetscErrorCode :: ierr
1002)
1003) #if UGRID_DEBUG
1004) call printMsg(option,'Partitioning')
1005) #endif
1006)
1007) ! create the partitioning
1008) call MatPartitioningCreate(option%mycomm,Part,ierr);CHKERRQ(ierr)
1009) ! MatPartitioningSetAdjacency sets the adjacency graph (matrix) of the
1010) ! thing to be partitioned. - petsc
1011) call MatPartitioningSetAdjacency(Part,Dual_mat,ierr);CHKERRQ(ierr)
1012) call MatPartitioningSetFromOptions(Part,ierr);CHKERRQ(ierr)
1013) ! MatPartitioningApply gets a partitioning for a matrix. For each local cell
1014) ! this tells the processor number that that cell is assigned to. - petsc
1015) ! is_new holds this information
1016) call MatPartitioningApply(Part,is_new,ierr);CHKERRQ(ierr)
1017) call MatPartitioningDestroy(Part,ierr);CHKERRQ(ierr)
1018)
1019) #if UGRID_DEBUG
1020) if (ugrid%grid_type == THREE_DIM_GRID) then
1021) call PetscViewerASCIIOpen(option%mycomm,'is_subsurf.out',viewer, &
1022) ierr);CHKERRQ(ierr)
1023) else
1024) call PetscViewerASCIIOpen(option%mycomm,'is_surf.out',viewer, &
1025) ierr);CHKERRQ(ierr)
1026) endif
1027) call ISView(is_new,viewer,ierr);CHKERRQ(ierr)
1028) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
1029) #endif
1030)
1031) ! calculate the number of local grid cells on each processor
1032) allocate(cell_counts(option%mycommsize))
1033) ! ISPartitioningCount takes a ISPartitioning and determines the number of
1034) ! resulting elements on each (partition) process - petsc
1035) call ISPartitioningCount(is_new,option%mycommsize,cell_counts, &
1036) ierr);CHKERRQ(ierr)
1037) num_cells_local_new = cell_counts(option%myrank+1)
1038) call MPI_Allreduce(num_cells_local_new,iflag,ONE_INTEGER_MPI,MPIU_INTEGER, &
1039) MPI_MIN,option%mycomm,ierr)
1040) deallocate(cell_counts)
1041) if (iflag < 1) then
1042) option%io_buffer = 'A processor core has been assigned zero cells.'
1043) call printErrMsg(option)
1044) endif
1045)
1046) end subroutine UGridPartition
1047)
1048) ! ************************************************************************** !
1049)
1050) subroutine UGridCreateOldVec(ugrid,option,elements_old, &
1051) num_cells_local_old, &
1052) is_new,is_scatter,stride)
1053) !
1054) ! UGridNaturalToPetsc: Deallocates a unstructured grid
1055) !
1056) ! Author: Glenn Hammond
1057) ! Date: 11/01/09
1058) !
1059) use Option_module
1060)
1061) implicit none
1062)
1063) #include "petsc/finclude/petscis.h"
1064) #include "petsc/finclude/petscis.h90"
1065) #include "petsc/finclude/petscviewer.h"
1066)
1067) type(grid_unstructured_type) :: ugrid
1068) type(option_type) :: option
1069) Vec :: elements_old
1070) PetscInt :: num_cells_local_old
1071) IS :: is_new
1072) IS :: is_scatter
1073) PetscInt :: stride
1074)
1075) PetscViewer :: viewer
1076) IS :: is_num
1077) PetscInt, pointer :: index_ptr(:)
1078) PetscErrorCode :: ierr
1079)
1080) ! calculate the global offsets in the new vector for each grid cell
1081)
1082) ! ISPartitioningToNumbering takes an ISPartitioning and on each processor
1083) ! generates an IS that contains a new global node number for each index
1084) ! based on the partitioning. - petsc
1085) call ISPartitioningToNumbering(is_new,is_num,ierr);CHKERRQ(ierr)
1086) call ISDestroy(is_new,ierr);CHKERRQ(ierr)
1087) call ISGetIndicesF90(is_num,index_ptr,ierr);CHKERRQ(ierr)
1088)
1089) ! Create a mapping of local indices to global strided (use block ids, not
1090) ! indices)
1091) call ISCreateBlock(option%mycomm,stride, &
1092) num_cells_local_old, &
1093) index_ptr,PETSC_COPY_VALUES,is_scatter, &
1094) ierr);CHKERRQ(ierr)
1095) call ISRestoreIndicesF90(is_num,index_ptr,ierr);CHKERRQ(ierr)
1096) call ISDestroy(is_num,ierr);CHKERRQ(ierr)
1097)
1098) #if UGRID_DEBUG
1099) if (ugrid%grid_type == THREE_DIM_GRID) then
1100) call PetscViewerASCIIOpen(option%mycomm,'is_scatter_elem_old_to_new_subsurf.out', &
1101) viewer,ierr);CHKERRQ(ierr)
1102) else
1103) call PetscViewerASCIIOpen(option%mycomm,'is_scatter_elem_old_to_new_surf.out', &
1104) viewer,ierr);CHKERRQ(ierr)
1105) endif
1106) call ISView(is_scatter,viewer,ierr);CHKERRQ(ierr)
1107) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
1108) #endif
1109)
1110) ! create another strided vector with the old cell/element distribution
1111) call VecCreate(option%mycomm,elements_old,ierr);CHKERRQ(ierr)
1112) call VecSetSizes(elements_old,stride*num_cells_local_old,PETSC_DECIDE, &
1113) ierr);CHKERRQ(ierr)
1114) call VecSetFromOptions(elements_old,ierr);CHKERRQ(ierr)
1115)
1116) end subroutine UGridCreateOldVec
1117)
1118) ! ************************************************************************** !
1119)
1120) subroutine UGridNaturalToPetsc(ugrid,option,elements_old,elements_local, &
1121) num_cells_local_new,stride,dual_offset, &
1122) natural_id_offset,is_scatter)
1123) !
1124) ! Deallocates a unstructured grid
1125) !
1126) ! Author: Glenn Hammond
1127) ! Date: 11/01/09
1128) !
1129)
1130) use Option_module
1131) use Utility_module, only: reallocateIntArray, DeallocateArray
1132)
1133) implicit none
1134)
1135) #include "petsc/finclude/petscmat.h"
1136) #include "petsc/finclude/petscmat.h90"
1137) #include "petsc/finclude/petscis.h"
1138) #include "petsc/finclude/petscis.h90"
1139) #include "petsc/finclude/petscviewer.h"
1140)
1141) type(grid_unstructured_type) :: ugrid
1142) type(option_type) :: option
1143) Vec :: elements_old, elements_local
1144) PetscInt :: num_cells_local_new
1145) PetscInt :: stride
1146) PetscInt :: dual_offset
1147) PetscInt :: natural_id_offset
1148) IS :: is_scatter
1149)
1150) Vec :: elements_petsc, elements_natural
1151) PetscViewer :: viewer
1152) VecScatter :: vec_scatter
1153) IS :: is_gather
1154)
1155) character(len=MAXSTRINGLENGTH) :: string
1156) PetscInt :: global_offset_new
1157) PetscInt :: local_id, idual, count
1158) PetscInt :: max_ghost_cell_count
1159) PetscInt :: temp_int
1160) PetscInt :: ghosted_id
1161) PetscInt :: ghost_cell_count
1162) PetscInt :: ghost_offset_new
1163) PetscInt :: dual_id
1164) PetscBool :: found
1165) PetscReal, pointer :: vec_ptr(:)
1166) PetscReal, pointer :: vec_ptr2(:)
1167) PetscInt, allocatable :: int_array(:)
1168) PetscInt, allocatable :: int_array2(:)
1169) PetscInt, allocatable :: int_array3(:)
1170) PetscInt, allocatable :: int_array4(:)
1171) PetscInt, allocatable :: int_array5(:)
1172) PetscInt, pointer :: int_array_pointer(:)
1173) PetscErrorCode :: ierr
1174)
1175) ! create a petsc vec to store all the information for each element
1176) ! based on the stride calculated above.
1177) call VecCreate(option%mycomm,elements_natural,ierr);CHKERRQ(ierr)
1178) call VecSetSizes(elements_natural, &
1179) stride*num_cells_local_new, &
1180) PETSC_DECIDE,ierr);CHKERRQ(ierr)
1181) call VecSetFromOptions(elements_natural,ierr);CHKERRQ(ierr)
1182)
1183) #if UGRID_DEBUG
1184) call printMsg(option,'Before element scatter')
1185) #endif
1186)
1187) ! scatter all the cell data from the old decomposition (as read in in
1188) ! parallel) to the more parmetis-calculated decomposition
1189) call VecScatterCreate(elements_old,PETSC_NULL_OBJECT,elements_natural,is_scatter, &
1190) vec_scatter,ierr);CHKERRQ(ierr)
1191) call ISDestroy(is_scatter,ierr);CHKERRQ(ierr)
1192) call VecScatterBegin(vec_scatter,elements_old,elements_natural, &
1193) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1194) call VecScatterEnd(vec_scatter,elements_old,elements_natural, &
1195) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1196) call VecScatterDestroy(vec_scatter,ierr);CHKERRQ(ierr)
1197)
1198) #if UGRID_DEBUG
1199) call printMsg(option,'After element scatter')
1200) if (ugrid%grid_type == THREE_DIM_GRID) then
1201) call PetscViewerASCIIOpen(option%mycomm,'elements_old_suburf.out',viewer, &
1202) ierr);CHKERRQ(ierr)
1203) else
1204) call PetscViewerASCIIOpen(option%mycomm,'elements_old_surf.out',viewer, &
1205) ierr);CHKERRQ(ierr)
1206) endif
1207) call VecView(elements_old,viewer,ierr);CHKERRQ(ierr)
1208) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
1209) #endif
1210)
1211) call VecDestroy(elements_old,ierr);CHKERRQ(ierr)
1212)
1213) #if UGRID_DEBUG
1214) if (ugrid%grid_type == THREE_DIM_GRID) then
1215) call PetscViewerASCIIOpen(option%mycomm,'elements_natural_subsurf.out',viewer, &
1216) ierr);CHKERRQ(ierr)
1217) else
1218) call PetscViewerASCIIOpen(option%mycomm,'elements_natural_surf.out',viewer, &
1219) ierr);CHKERRQ(ierr)
1220) endif
1221) call VecView(elements_natural,viewer,ierr);CHKERRQ(ierr)
1222) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
1223) #endif
1224)
1225) ! update global offset based on new partitioning
1226) global_offset_new = 0
1227) call MPI_Exscan(num_cells_local_new,global_offset_new, &
1228) ONE_INTEGER_MPI,MPIU_INTEGER,MPI_SUM,option%mycomm,ierr)
1229) ugrid%global_offset = global_offset_new
1230)
1231) allocate(ugrid%cell_ids_natural(num_cells_local_new))
1232) ugrid%cell_ids_natural = 0
1233)
1234) ! look at all connections and determine how many are non-local, and create
1235) ! a listof indices
1236) call VecDuplicate(elements_natural,elements_petsc,ierr);CHKERRQ(ierr)
1237) call VecCopy(elements_natural,elements_petsc,ierr);CHKERRQ(ierr)
1238)
1239) #if UGRID_DEBUG
1240) call printMsg(option,'Lists of ids')
1241) #endif
1242)
1243) ! now we unpack the decomposed cell data
1244)
1245) ! store the natural grid cell id for each local cell as read from the grid
1246) ! file
1247) call VecGetArrayF90(elements_natural,vec_ptr,ierr);CHKERRQ(ierr)
1248) do local_id = 1, num_cells_local_new
1249) ! Cell id read from explicit grid input file is second entry in block
1250) ! It may match the first entry (the calculated natural id based on the
1251) ! order that cells were read), but it need not.
1252) ugrid%cell_ids_natural(local_id) = &
1253) int(abs(vec_ptr((local_id-1)*stride+natural_id_offset)))
1254) enddo
1255) call VecRestoreArrayF90(elements_natural,vec_ptr,ierr);CHKERRQ(ierr)
1256)
1257) #if UGRID_DEBUG
1258) write(string,*) option%myrank
1259) string = 'natural_ids' // trim(adjustl(string)) // '.out'
1260) open(unit=86,file=trim(string))
1261) do local_id = 1, num_cells_local_new
1262) write(86,'(i5)') ugrid%cell_ids_natural(local_id)
1263) enddo
1264) close(86)
1265) #endif
1266)
1267) ! make a list of petsc ids for each local cell (you simply take the global
1268) ! offset and add it to the local contiguous cell ids on each processor
1269) allocate(int_array(num_cells_local_new))
1270) do local_id = 1, num_cells_local_new
1271) int_array(local_id) = local_id+global_offset_new
1272) enddo
1273)
1274) ! make the arrays zero-based
1275) int_array = int_array - 1
1276) ugrid%cell_ids_natural = ugrid%cell_ids_natural - 1
1277) ! create an application ordering (mapping of natural to petsc ordering)
1278) call AOCreateBasic(option%mycomm,num_cells_local_new, &
1279) ugrid%cell_ids_natural,int_array, &
1280) ugrid%ao_natural_to_petsc,ierr);CHKERRQ(ierr)
1281) deallocate(int_array)
1282) ! make cell_ids_natural 1-based again
1283) ugrid%cell_ids_natural = ugrid%cell_ids_natural + 1
1284)
1285) #if UGRID_DEBUG
1286) if (ugrid%grid_type == THREE_DIM_GRID) then
1287) call PetscViewerASCIIOpen(option%mycomm,'ao_subsurf.out',viewer, &
1288) ierr);CHKERRQ(ierr)
1289) else
1290) call PetscViewerASCIIOpen(option%mycomm,'ao_surf.out',viewer, &
1291) ierr);CHKERRQ(ierr)
1292) endif
1293) call AOView(ugrid%ao_natural_to_petsc,viewer,ierr);CHKERRQ(ierr)
1294) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
1295) #endif
1296)
1297) ! The below creates a list of cells ids for the duals and converts them
1298) ! to petsc ordering
1299)
1300) ! count the number of cells and their duals
1301) call VecGetArrayF90(elements_natural,vec_ptr,ierr);CHKERRQ(ierr)
1302) count = 0
1303) do local_id=1, num_cells_local_new
1304) count = count + 1
1305) do idual = 1, ugrid%max_ndual_per_cell
1306) dual_id = int(vec_ptr(idual + dual_offset + (local_id-1)*stride))
1307) if (dual_id < 1) exit ! here we hit the 0 at the end of last dual
1308) count = count + 1
1309) enddo
1310) enddo
1311)
1312) ! allocate and fill an array with the natural cell and dual ids
1313) allocate(int_array(count))
1314) count = 0
1315) do local_id=1, num_cells_local_new
1316) count = count + 1
1317) int_array(count) = ugrid%cell_ids_natural(local_id)
1318) do idual = 1, ugrid%max_ndual_per_cell
1319) dual_id = int(vec_ptr(idual + dual_offset + (local_id-1)*stride))
1320) if (dual_id < 1) exit ! again we hit the 0
1321) count = count + 1
1322) int_array(count) = dual_id
1323) enddo
1324) enddo
1325) call VecRestoreArrayF90(elements_natural,vec_ptr,ierr);CHKERRQ(ierr)
1326)
1327) #if UGRID_DEBUG
1328) call printMsg(option,'Application ordering')
1329) #endif
1330)
1331) ! convert the dual ids in int_array from natural to petsc numbering
1332) int_array = int_array - 1
1333) call AOApplicationToPetsc(ugrid%ao_natural_to_petsc,count, &
1334) int_array,ierr);CHKERRQ(ierr)
1335) int_array = int_array + 1
1336)
1337) #if UGRID_DEBUG
1338) call printMsg(option,'PETSc-ordered duals')
1339) #endif
1340)
1341) ! load mapped petsc-ordered dual ids back into duplicated vector
1342) ! exactly the opposite operation of when we loaded the temporary int_array
1343) ! vector
1344) call VecGetArrayF90(elements_petsc,vec_ptr,ierr);CHKERRQ(ierr)
1345) !geh: do not believe that we need elements_natural here
1346) ! call VecGetArrayF90(elements_natural,vec_ptr2,ierr)
1347) allocate(ugrid%cell_ids_petsc(num_cells_local_new))
1348) count = 0
1349) do local_id=1, num_cells_local_new
1350) count = count + 1
1351) ! extract the petsc id for the cell
1352) ugrid%cell_ids_petsc(local_id) = int_array(count)
1353) ! store it in the elements_petsc vector too
1354) vec_ptr((local_id-1)*stride+1) = int_array(count)
1355) do idual = 1, ugrid%max_ndual_per_cell
1356) !geh dual_id = vec_ptr2(idual + dual_offset + (local_id-1)*stride)
1357) dual_id = int(vec_ptr(idual + dual_offset + (local_id-1)*stride))
1358) if (dual_id < 1) exit
1359) count = count + 1
1360) ! store the petsc numbered duals in the vector also
1361) vec_ptr(idual + dual_offset + (local_id-1)*stride) = int_array(count)
1362) enddo
1363) enddo
1364) call VecRestoreArrayF90(elements_petsc,vec_ptr,ierr);CHKERRQ(ierr)
1365) !geh call VecRestoreArrayF90(elements_natural,vec_ptr2,ierr)
1366) deallocate(int_array)
1367)
1368) #if UGRID_DEBUG
1369) if (ugrid%grid_type == THREE_DIM_GRID) then
1370) call PetscViewerASCIIOpen(option%mycomm,'elements_petsc_subsurf.out',viewer, &
1371) ierr);CHKERRQ(ierr)
1372) else
1373) call PetscViewerASCIIOpen(option%mycomm,'elements_petsc_surf.out',viewer, &
1374) ierr);CHKERRQ(ierr)
1375) endif
1376) call VecView(elements_petsc,viewer,ierr);CHKERRQ(ierr)
1377) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
1378) #endif
1379)
1380) ! make a list of ghosted ids in petsc numbering
1381) #if UGRID_DEBUG
1382) call printMsg(option,'Renumbering ghost ids to petsc numbering')
1383) #endif
1384)
1385) call VecGetArrayF90(elements_petsc,vec_ptr,ierr);CHKERRQ(ierr)
1386) ghost_cell_count = 0
1387) ! allocate a temporarily-sized array
1388) ! geh: this assumes that the number of ghost cells will not exceed the number
1389) ! of local and 100 is used to ensure that if this is not true, the array
1390) ! is still large enough
1391) max_ghost_cell_count = max(num_cells_local_new,100)
1392) allocate(int_array_pointer(max_ghost_cell_count))
1393) int_array_pointer = 0
1394) ! loop over all duals and find the off-processor cells on the other
1395) ! end of a dual
1396) do local_id=1, num_cells_local_new
1397) do idual = 1, ugrid%max_ndual_per_cell
1398) dual_id = int(vec_ptr(idual + dual_offset + (local_id-1)*stride))
1399) found = PETSC_FALSE
1400) if (dual_id < 1) exit
1401) if (dual_id <= global_offset_new .or. &
1402) dual_id > global_offset_new + num_cells_local_new) then
1403) ghost_cell_count = ghost_cell_count + 1
1404) ! reallocate the ghost cell array if necessary
1405) if (ghost_cell_count > max_ghost_cell_count) then
1406) call reallocateIntArray(int_array_pointer,max_ghost_cell_count)
1407) endif
1408) int_array_pointer(ghost_cell_count) = dual_id
1409) vec_ptr(idual + dual_offset + (local_id-1)*stride) = &
1410) -ghost_cell_count
1411) ! temporarily store the index of the int_array_pointer
1412) ! flag negative
1413) else
1414) ! if dual_id is in petsc numbering, then the local ids is:
1415) vec_ptr(idual + dual_offset + (local_id-1)*stride) = &
1416) dual_id - global_offset_new
1417) endif
1418) enddo
1419) enddo
1420) call VecRestoreArrayF90(elements_petsc,vec_ptr,ierr);CHKERRQ(ierr)
1421)
1422) #if UGRID_DEBUG
1423) if (ugrid%grid_type == THREE_DIM_GRID) then
1424) call PetscViewerASCIIOpen(option%mycomm, &
1425) 'elements_local_dual_unsorted_subsurf.out',viewer, &
1426) ierr);CHKERRQ(ierr)
1427) else
1428) call PetscViewerASCIIOpen(option%mycomm, &
1429) 'elements_local_dual_unsorted_surf.out',viewer, &
1430) ierr);CHKERRQ(ierr)
1431) endif
1432) call VecView(elements_petsc,viewer,ierr);CHKERRQ(ierr)
1433) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
1434) #endif
1435)
1436)
1437) #if UGRID_DEBUG
1438) call printMsg(option,' Sorting local ghost ids')
1439) #endif
1440)
1441) if (ghost_cell_count > 0) then
1442) ! sort ghost cell ids
1443) allocate(int_array2(ghost_cell_count))
1444) do ghosted_id = 1, ghost_cell_count
1445) int_array2(ghosted_id) = ghosted_id
1446) enddo
1447)
1448) ! sort with permutation
1449) ! int_array_pointer = ghost ids unsorted before and after
1450) ! int_array2 = permutation
1451) ! convert to 0-based
1452) int_array_pointer = int_array_pointer-1
1453) int_array2 = int_array2 - 1
1454) call PetscSortIntWithPermutation(ghost_cell_count,int_array_pointer, &
1455) int_array2,ierr);CHKERRQ(ierr)
1456) ! convert to 1-based
1457) int_array_pointer = int_array_pointer+1
1458) int_array2 = int_array2 + 1
1459)
1460) ! determine how many duplicates
1461) allocate(int_array3(ghost_cell_count))
1462) allocate(int_array4(ghost_cell_count))
1463) allocate(int_array5(ghost_cell_count))
1464) int_array3 = 0
1465) temp_int = 1
1466) int_array3(temp_int) = int_array_pointer(int_array2(1))
1467) ! do not change ghosted_id = 1 to ghosted_id = 2 as the first value in
1468) ! int_array2() will not be set correctly.
1469) do ghosted_id = 1, ghost_cell_count
1470) if (int_array3(temp_int) < &
1471) int_array_pointer(int_array2(ghosted_id))) then
1472) temp_int = temp_int + 1
1473) int_array3(temp_int) = int_array_pointer(int_array2(ghosted_id))
1474) endif
1475) int_array5(int_array2(ghosted_id)) = ghosted_id
1476) int_array4(ghosted_id) = temp_int
1477) enddo
1478)
1479) ghost_cell_count = temp_int
1480) allocate(ugrid%ghost_cell_ids_petsc(ghost_cell_count))
1481) ugrid%ghost_cell_ids_petsc = 0
1482)
1483) ugrid%ghost_cell_ids_petsc(1:ghost_cell_count) = &
1484) int_array3(1:ghost_cell_count)
1485)
1486) #if UGRID_DEBUG
1487) call printMsg(option,' Remappping ghost ids')
1488) #endif
1489)
1490) ! remap of duals of ghost cells
1491) call VecGetArrayF90(elements_petsc,vec_ptr,ierr);CHKERRQ(ierr)
1492) do local_id=1, num_cells_local_new
1493) do idual = 1, ugrid%max_ndual_per_cell
1494) ! dual_id is now the negative of the local unsorted ghost cell id
1495) dual_id = int(vec_ptr(idual + dual_offset + (local_id-1)*stride))
1496) ! dual_id = 0: not assigned
1497) ! dual_id > 0: assigned to local cell
1498) ! dual_id < 0: assigned to ghost cell
1499) if (dual_id < 0) then
1500) vec_ptr(idual + dual_offset + (local_id-1)*stride) = &
1501) int_array4(int_array5(-dual_id)) + num_cells_local_new
1502) endif
1503) enddo
1504) enddo
1505) call VecRestoreArrayF90(elements_petsc,vec_ptr,ierr);CHKERRQ(ierr)
1506)
1507) deallocate(int_array2)
1508) deallocate(int_array3)
1509) deallocate(int_array4)
1510) deallocate(int_array5)
1511) endif
1512) call DeallocateArray(int_array_pointer)
1513)
1514) ugrid%nlmax = num_cells_local_new
1515) ugrid%num_ghost_cells = ghost_cell_count
1516) ugrid%ngmax = &
1517) num_cells_local_new + ghost_cell_count
1518)
1519)
1520) #if UGRID_DEBUG
1521) if (ugrid%grid_type == THREE_DIM_GRID) then
1522) call PetscViewerASCIIOpen(option%mycomm,'elements_local_dual_subsurf.out', &
1523) viewer,ierr);CHKERRQ(ierr)
1524) else
1525) call PetscViewerASCIIOpen(option%mycomm,'elements_local_dual_surf.out', &
1526) viewer,ierr);CHKERRQ(ierr)
1527) endif
1528) call VecView(elements_petsc,viewer,ierr);CHKERRQ(ierr)
1529) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
1530) #endif
1531)
1532) #if UGRID_DEBUG
1533) call printMsg(option,'Resizing natural cell id array')
1534) #endif
1535)
1536) ! Resize cell_ids_natural to include ghosted cells
1537) allocate(int_array(ugrid%nlmax))
1538) int_array(:) = ugrid%cell_ids_natural(:)
1539) deallocate(ugrid%cell_ids_natural)
1540) allocate(ugrid%cell_ids_natural(ugrid%ngmax))
1541) ugrid%cell_ids_natural(:) = UNINITIALIZED_INTEGER
1542) ugrid%cell_ids_natural(1:ugrid%nlmax) = int_array(:)
1543) deallocate(int_array)
1544) call VecGetArrayF90(elements_petsc,vec_ptr,ierr);CHKERRQ(ierr)
1545) call VecGetArrayF90(elements_natural,vec_ptr2,ierr);CHKERRQ(ierr)
1546) do local_id=1, ugrid%nlmax
1547) do idual = 1, ugrid%max_ndual_per_cell
1548) dual_id = int(vec_ptr(idual + dual_offset + (local_id-1)*stride))
1549) if (dual_id < 1) exit
1550) if (dual_id > ugrid%nlmax) then
1551) ugrid%cell_ids_natural(dual_id) = &
1552) int(vec_ptr2(idual + dual_offset + (local_id-1)*stride))
1553) endif
1554) enddo
1555) enddo
1556) if (minval(ugrid%cell_ids_natural) < 1) then
1557) write(string,*) minval( ugrid%cell_ids_natural)
1558) option%io_buffer = 'Negative natural id: ' // trim(adjustl(string))
1559) call printErrMsgByRank(option)
1560) endif
1561) call VecRestoreArrayF90(elements_petsc,vec_ptr,ierr);CHKERRQ(ierr)
1562) call VecRestoreArrayF90(elements_natural,vec_ptr2,ierr);CHKERRQ(ierr)
1563) call VecDestroy(elements_natural,ierr);CHKERRQ(ierr)
1564)
1565) ! NOW START ON GHOSTING CELLS
1566)
1567) #if UGRID_DEBUG
1568) call printMsg(option,'Ghosting local element vectors')
1569) #endif
1570)
1571) ! load cell neighbors into array
1572) ! start first index at zero to store # duals for a cell
1573) allocate(ugrid%cell_neighbors_local_ghosted( &
1574) 0:ugrid%max_ndual_per_cell,ugrid%nlmax))
1575) ugrid%cell_neighbors_local_ghosted = 0
1576) call VecGetArrayF90(elements_petsc,vec_ptr,ierr);CHKERRQ(ierr)
1577) do local_id=1, ugrid%nlmax
1578) count = 0
1579) do idual = 1, ugrid%max_ndual_per_cell
1580) dual_id = int(vec_ptr(idual + dual_offset + (local_id-1)*stride))
1581) if (dual_id < 1) exit
1582) count = count + 1
1583) ! flag ghosted cells in dual as negative
1584) !geh: these negative dual ids are used later in UGridDMCreateJacobian()
1585) ! to specify off processor connectity in the Jacobian
1586) if (dual_id > ugrid%nlmax) dual_id = -dual_id
1587) ugrid%cell_neighbors_local_ghosted(idual,local_id) = dual_id
1588) enddo
1589) ! set the # of duals in for the cell
1590) ugrid%cell_neighbors_local_ghosted(0,local_id) = count
1591) enddo
1592) call VecRestoreArrayF90(elements_petsc,vec_ptr,ierr);CHKERRQ(ierr)
1593)
1594) ! need to create a local ghosted vector in which we can collect element info
1595) ! including ghost cells
1596) !call VecCreateSeq(PETSC_COMM_SELF,stride*ugrid%ngmax, &
1597) ! elements_local,ierr)
1598) call VecCreate(PETSC_COMM_SELF,elements_local,ierr);CHKERRQ(ierr)
1599) call VecSetSizes(elements_local,stride*ugrid%ngmax,PETSC_DECIDE, &
1600) ierr);CHKERRQ(ierr)
1601) call VecSetBlockSize(elements_local,stride,ierr);CHKERRQ(ierr)
1602) call VecSetFromOptions(elements_local,ierr);CHKERRQ(ierr)
1603) allocate(int_array(ugrid%ngmax))
1604) int_array(1:ugrid%nlmax) = &
1605) ugrid%cell_ids_petsc(:)
1606) if (ugrid%num_ghost_cells > 0) then
1607) int_array(ugrid%nlmax+1:ugrid%ngmax) = &
1608) ugrid%ghost_cell_ids_petsc(:)
1609) endif
1610) int_array = int_array-1
1611) call ISCreateBlock(option%mycomm,stride,ugrid%ngmax, &
1612) int_array,PETSC_COPY_VALUES,is_scatter, &
1613) ierr);CHKERRQ(ierr)
1614) do ghosted_id = 1, ugrid%ngmax
1615) int_array(ghosted_id) = ghosted_id-1
1616) enddo
1617) call ISCreateBlock(option%mycomm,stride,ugrid%ngmax, &
1618) int_array,PETSC_COPY_VALUES,is_gather,ierr);CHKERRQ(ierr)
1619) deallocate(int_array)
1620)
1621) #if UGRID_DEBUG
1622) if (ugrid%grid_type == THREE_DIM_GRID) then
1623) call PetscViewerASCIIOpen(option%mycomm,'is_scatter_elem_local_to_ghost_subsurf.out',viewer, &
1624) ierr);CHKERRQ(ierr)
1625) else
1626) call PetscViewerASCIIOpen(option%mycomm,'is_scatter_elem_local_to_ghost_surf.out',viewer, &
1627) ierr);CHKERRQ(ierr)
1628) endif
1629) call ISView(is_scatter,viewer,ierr);CHKERRQ(ierr)
1630) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
1631) if (ugrid%grid_type == THREE_DIM_GRID) then
1632) call PetscViewerASCIIOpen(option%mycomm,'is_gather_elem_local_to_ghost_subsurf.out',viewer, &
1633) ierr);CHKERRQ(ierr)
1634) else
1635) call PetscViewerASCIIOpen(option%mycomm,'is_gather_elem_local_to_ghost_surf.out',viewer, &
1636) ierr);CHKERRQ(ierr)
1637) endif
1638) call ISView(is_gather,viewer,ierr);CHKERRQ(ierr)
1639) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
1640) #endif
1641)
1642) call VecScatterCreate(elements_petsc,is_scatter,elements_local,is_gather, &
1643) vec_scatter,ierr);CHKERRQ(ierr)
1644) call ISDestroy(is_scatter,ierr);CHKERRQ(ierr)
1645) call ISDestroy(is_gather,ierr);CHKERRQ(ierr)
1646) call VecScatterBegin(vec_scatter,elements_petsc,elements_local, &
1647) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1648) call VecScatterEnd(vec_scatter,elements_petsc,elements_local, &
1649) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1650) call VecScatterDestroy(vec_scatter,ierr);CHKERRQ(ierr)
1651) call VecDestroy(elements_petsc,ierr);CHKERRQ(ierr)
1652)
1653) #if UGRID_DEBUG
1654) write(string,*) option%myrank
1655) if (ugrid%grid_type == THREE_DIM_GRID) then
1656) string = 'elements_local' // trim(adjustl(string)) // '_subsurf.out'
1657) else
1658) string = 'elements_local' // trim(adjustl(string)) // '_surf.out'
1659) endif
1660) call PetscViewerASCIIOpen(PETSC_COMM_SELF,trim(string),viewer, &
1661) ierr);CHKERRQ(ierr)
1662) call VecView(elements_local,viewer,ierr);CHKERRQ(ierr)
1663) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
1664)
1665) call printMsg(option,'Scatter/gathering local ghosted vertices')
1666) #endif
1667)
1668) end subroutine UGridNaturalToPetsc
1669)
1670) ! ************************************************************************** !
1671)
1672) subroutine UGridDestroy(unstructured_grid)
1673) !
1674) ! Deallocates a unstructured grid
1675) !
1676) ! Author: Glenn Hammond
1677) ! Date: 11/01/09
1678) !
1679)
1680) use Utility_module, only : DeallocateArray
1681)
1682) implicit none
1683)
1684) type(grid_unstructured_type), pointer :: unstructured_grid
1685)
1686) PetscErrorCode :: ierr
1687)
1688) if (.not.associated(unstructured_grid)) return
1689)
1690) ! variables for all unstructured grids
1691) call DeallocateArray(unstructured_grid%hash)
1692) call DeallocateArray(unstructured_grid%cell_ids_natural)
1693) call DeallocateArray(unstructured_grid%cell_ids_petsc)
1694) call DeallocateArray(unstructured_grid%ghost_cell_ids_petsc)
1695) call UGridExplicitDestroy(unstructured_grid%explicit_grid)
1696) call UGridPolyhedraDestroy(unstructured_grid%polyhedra_grid)
1697) if (unstructured_grid%ao_natural_to_petsc /= 0) then
1698) call AODestroy(unstructured_grid%ao_natural_to_petsc,ierr);CHKERRQ(ierr)
1699) endif
1700)
1701) ! variables for implicit unstructured grids
1702) call DeallocateArray(unstructured_grid%cell_type)
1703) call DeallocateArray(unstructured_grid%cell_vertices)
1704) call DeallocateArray(unstructured_grid%face_to_cell_ghosted)
1705) call DeallocateArray(unstructured_grid%connection_to_face)
1706) call DeallocateArray(unstructured_grid%face_to_vertex_natural)
1707) call DeallocateArray(unstructured_grid%face_to_vertex)
1708) call DeallocateArray(unstructured_grid%cell_to_face_ghosted)
1709) call DeallocateArray(unstructured_grid%vertex_ids_natural)
1710) call DeallocateArray(unstructured_grid%cell_neighbors_local_ghosted)
1711) if (associated(unstructured_grid%vertices)) &
1712) deallocate(unstructured_grid%vertices)
1713) nullify(unstructured_grid%vertices)
1714) if (associated(unstructured_grid%face_centroid)) &
1715) deallocate(unstructured_grid%face_centroid)
1716) nullify(unstructured_grid%face_centroid)
1717) call DeallocateArray(unstructured_grid%face_area)
1718)
1719) deallocate(unstructured_grid)
1720) nullify(unstructured_grid)
1721)
1722) end subroutine UGridDestroy
1723)
1724) ! ************************************************************************** !
1725)
1726) subroutine UGridDMDestroy(ugdm)
1727) !
1728) ! Deallocates a unstructured grid distributed mesh
1729) !
1730) ! Author: Glenn Hammond
1731) ! Date: 11/01/09
1732) !
1733)
1734) implicit none
1735)
1736) type(ugdm_type), pointer :: ugdm
1737)
1738) PetscErrorCode :: ierr
1739)
1740) if (.not.associated(ugdm)) return
1741)
1742) call ISDestroy(ugdm%is_ghosted_local,ierr);CHKERRQ(ierr)
1743) call ISDestroy(ugdm%is_local_local,ierr);CHKERRQ(ierr)
1744) call ISDestroy(ugdm%is_ghosted_petsc,ierr);CHKERRQ(ierr)
1745) call ISDestroy(ugdm%is_local_petsc,ierr);CHKERRQ(ierr)
1746) call ISDestroy(ugdm%is_ghosts_local,ierr);CHKERRQ(ierr)
1747) call ISDestroy(ugdm%is_ghosts_petsc,ierr);CHKERRQ(ierr)
1748) call ISDestroy(ugdm%is_local_natural,ierr);CHKERRQ(ierr)
1749) call VecScatterDestroy(ugdm%scatter_ltog,ierr);CHKERRQ(ierr)
1750) call VecScatterDestroy(ugdm%scatter_gtol,ierr);CHKERRQ(ierr)
1751) call VecScatterDestroy(ugdm%scatter_ltol,ierr);CHKERRQ(ierr)
1752) call VecScatterDestroy(ugdm%scatter_gton,ierr);CHKERRQ(ierr)
1753) call ISLocalToGlobalMappingDestroy(ugdm%mapping_ltog,ierr);CHKERRQ(ierr)
1754) call VecDestroy(ugdm%global_vec,ierr);CHKERRQ(ierr)
1755) call VecDestroy(ugdm%local_vec,ierr);CHKERRQ(ierr)
1756) call VecScatterDestroy(ugdm%scatter_bet_grids,ierr);CHKERRQ(ierr)
1757) call VecScatterDestroy(ugdm%scatter_bet_grids_1dof,ierr);CHKERRQ(ierr)
1758) call VecScatterDestroy(ugdm%scatter_bet_grids_ndof,ierr);CHKERRQ(ierr)
1759) ! ugdm%ao_natural_to_petsc is a pointer to ugrid%ao_natural_to_petsc. Do
1760) ! not destroy here.
1761) ugdm%ao_natural_to_petsc = 0
1762) deallocate(ugdm)
1763) nullify(ugdm)
1764)
1765) end subroutine UGridDMDestroy
1766)
1767) ! ************************************************************************** !
1768)
1769) subroutine UGridExplicitDestroy(explicit_grid)
1770) !
1771) ! Deallocates an explicit unstructured grid object
1772) !
1773) ! Author: Glenn Hammond
1774) ! Date: 05/14/12
1775) !
1776)
1777) use Utility_module, only : DeallocateArray
1778)
1779) implicit none
1780)
1781) type(unstructured_explicit_type), pointer :: explicit_grid
1782)
1783) PetscErrorCode :: ierr
1784)
1785) if (.not.associated(explicit_grid)) return
1786)
1787) call DeallocateArray(explicit_grid%cell_ids)
1788) call DeallocateArray(explicit_grid%cell_volumes)
1789) if (associated(explicit_grid%cell_centroids)) &
1790) deallocate(explicit_grid%cell_centroids)
1791) nullify(explicit_grid%cell_centroids)
1792) call DeallocateArray(explicit_grid%connections)
1793) call DeallocateArray(explicit_grid%face_areas)
1794) if (associated(explicit_grid%face_centroids)) &
1795) deallocate(explicit_grid%face_centroids)
1796) nullify(explicit_grid%face_centroids)
1797)
1798) deallocate(explicit_grid)
1799) nullify(explicit_grid)
1800)
1801) end subroutine UGridExplicitDestroy
1802)
1803) ! ************************************************************************** !
1804) !> This routine deallocates a polyhedra unstructured grid object
1805) !!
1806) !> @author
1807) !! Gautam Bisht, LBL
1808) !!
1809) !! date: 09/29/13
1810) ! ************************************************************************** !
1811) subroutine UGridPolyhedraDestroy(polyhedra_grid)
1812)
1813) use Utility_module, only : DeallocateArray
1814)
1815) implicit none
1816)
1817) type(unstructured_polyhedra_type), pointer :: polyhedra_grid
1818)
1819) PetscErrorCode :: ierr
1820)
1821) if (.not.associated(polyhedra_grid)) return
1822)
1823) call DeallocateArray(polyhedra_grid%cell_ids)
1824) call DeallocateArray(polyhedra_grid%cell_nfaces)
1825) call DeallocateArray(polyhedra_grid%cell_nverts)
1826) call DeallocateArray(polyhedra_grid%cell_faceids)
1827) call DeallocateArray(polyhedra_grid%cell_vertids)
1828) call DeallocateArray(polyhedra_grid%cell_volumes)
1829) call DeallocateArray(polyhedra_grid%face_ids)
1830) call DeallocateArray(polyhedra_grid%face_cellids)
1831) call DeallocateArray(polyhedra_grid%face_nverts)
1832) call DeallocateArray(polyhedra_grid%face_vertids)
1833) call DeallocateArray(polyhedra_grid%face_areas)
1834)
1835) if (associated(polyhedra_grid%cell_centroids)) &
1836) deallocate(polyhedra_grid%cell_centroids)
1837) nullify(polyhedra_grid%cell_centroids)
1838) if (associated(polyhedra_grid%face_centroids)) &
1839) deallocate(polyhedra_grid%face_centroids)
1840) nullify(polyhedra_grid%face_centroids)
1841) if (associated(polyhedra_grid%vertex_coordinates)) &
1842) deallocate(polyhedra_grid%vertex_coordinates)
1843) nullify(polyhedra_grid%vertex_coordinates)
1844) if (associated(polyhedra_grid%ugridf2pgridf)) deallocate(polyhedra_grid%ugridf2pgridf)
1845) nullify(polyhedra_grid%ugridf2pgridf)
1846)
1847) if (associated(polyhedra_grid%uface_localids)) deallocate(polyhedra_grid%uface_localids)
1848) nullify(polyhedra_grid%uface_localids)
1849) if (associated(polyhedra_grid%uface_nverts)) deallocate(polyhedra_grid%uface_nverts)
1850) nullify(polyhedra_grid%uface_nverts)
1851) if (associated(polyhedra_grid%uface_natvertids)) deallocate(polyhedra_grid%uface_natvertids)
1852) nullify(polyhedra_grid%uface_natvertids)
1853) if (associated(polyhedra_grid%uface_left_natcellids)) deallocate(polyhedra_grid%uface_left_natcellids)
1854) nullify(polyhedra_grid%uface_left_natcellids)
1855) if (associated(polyhedra_grid%uface_right_natcellids)) deallocate(polyhedra_grid%uface_right_natcellids)
1856) nullify(polyhedra_grid%uface_right_natcellids)
1857)
1858) deallocate(polyhedra_grid)
1859) nullify(polyhedra_grid)
1860)
1861) end subroutine UGridPolyhedraDestroy
1862)
1863) end module Grid_Unstructured_Aux_module