grid_unstructured_polyhedra.F90 coverage: 87.50 %func 73.65 %block
1) module Grid_Unstructured_Polyhedra_module
2)
3) use Geometry_module
4) use Grid_Unstructured_Aux_module
5) use Grid_Unstructured_Cell_module
6) use PFLOTRAN_Constants_module
7)
8) implicit none
9)
10) private
11)
12) #include "petsc/finclude/petscsys.h"
13) #include "petsc/finclude/petscvec.h"
14) #include "petsc/finclude/petscvec.h90"
15) #include "petsc/finclude/petscis.h"
16) #include "petsc/finclude/petscis.h90"
17) #if defined(SCORPIO)
18) include "scorpiof.h"
19) #endif
20)
21) public :: UGridPolyhedraRead, &
22) UGridPolyhedraDecompose, &
23) UGridPolyhedraSetCellCentroids, &
24) UGridPolyhedraComputeInternConnect, &
25) UGridPolyhedraComputeVolumes, &
26) UGridPolyhedraPopulateConnection, &
27) UGridPolyhedraGetCellsInRectangle, &
28) UGridPolyhedraComputeOutputInfo
29)
30) contains
31)
32) ! ************************************************************************** !
33)
34) subroutine UGridPolyhedraRead(ugrid, filename, option)
35) !
36) ! This routine reads unstructured polyhedra grid in ASCII.
37) !
38) ! Author: Gautam Bisht, LBL
39) ! Date: 09/29/13
40) !
41)
42) use Input_Aux_module
43) use Option_module
44) use String_module
45)
46) implicit none
47)
48) type(grid_unstructured_type) :: ugrid
49) character(len=MAXSTRINGLENGTH) :: filename
50) type(option_type) :: option
51)
52) type(unstructured_polyhedra_type), pointer :: pgrid
53) type(input_type), pointer :: input
54) character(len=MAXSTRINGLENGTH) :: string, hint
55) character(len=MAXWORDLENGTH) :: word, card
56) PetscInt :: fileid, icell, iface, irank, ivert
57) PetscInt :: remainder, temp_int, num_to_read
58)
59) PetscInt :: num_cells
60) PetscInt :: num_cells_local
61) PetscInt :: num_cells_local_save
62) PetscInt :: num_faces
63) PetscInt :: num_faces_local
64) PetscInt :: num_faces_local_save
65) PetscInt :: num_vertices
66) PetscInt :: num_vertices_local
67) PetscInt :: num_vertices_local_save
68) PetscInt :: max_nvert_per_cell
69) PetscInt :: max_nface_per_cell
70) PetscInt :: max_nvert_per_face
71) PetscMPIInt :: status_mpi(MPI_STATUS_SIZE)
72) PetscMPIInt :: int_mpi
73) PetscErrorCode :: ierr
74) PetscReal, allocatable :: temp_real_array(:,:)
75) PetscInt, allocatable :: nfaces_per_proc(:)
76)
77) pgrid => ugrid%polyhedra_grid
78) allocate(nfaces_per_proc(option%mycommsize))
79) nfaces_per_proc = 0
80) num_faces_local_save = 0
81)
82) max_nvert_per_cell = -1
83) if (option%myrank == option%io_rank) then
84) fileid = 86
85) input => InputCreate(fileid,filename,option)
86)
87) call InputReadPflotranString(input,option)
88) ! read CELL card, though we already know the
89) call InputReadWord(input,option,card,PETSC_TRUE)
90) word = 'CELLS'
91) call InputErrorMsg(input,option,word,card)
92) if (.not.StringCompare(word,card)) then
93) option%io_buffer = 'Unrecognized keyword "' // trim(card) // &
94) '" in explicit grid file.'
95) call printErrMsgByRank(option)
96) endif
97)
98) hint = 'Polyhedra Unstructured Grid CELLS'
99) call InputReadInt(input,option,temp_int)
100) call InputErrorMsg(input,option,'number of cells',hint)
101) endif
102)
103) call MPI_Bcast(temp_int,ONE_INTEGER_MPI,MPI_INTEGER,option%io_rank, &
104) option%mycomm,ierr)
105) num_cells = temp_int
106) pgrid%num_cells_global = num_cells
107) ugrid%nmax = num_cells
108)
109) ! divide cells across ranks
110) num_cells_local = num_cells/option%mycommsize
111) num_cells_local_save = num_cells_local
112) remainder = num_cells - &
113) num_cells_local*option%mycommsize
114) if (option%myrank < remainder) num_cells_local = &
115) num_cells_local + 1
116)
117) ! allocate memory
118) allocate(pgrid%cell_ids(num_cells_local))
119) allocate(pgrid%cell_nfaces(num_cells_local))
120) allocate(pgrid%cell_nverts(num_cells_local))
121) allocate(pgrid%cell_volumes(num_cells_local))
122) allocate(pgrid%cell_centroids(num_cells_local))
123) pgrid%cell_ids = 0
124) pgrid%cell_nfaces = 0
125) pgrid%cell_nverts = 0
126) pgrid%cell_volumes = 0.d0
127) do icell = 1, num_cells_local
128) pgrid%cell_centroids(icell)%x = 0.d0
129) pgrid%cell_centroids(icell)%y = 0.d0
130) pgrid%cell_centroids(icell)%z = 0.d0
131) enddo
132)
133) ! Read all cells from ASCII file through io_rank and communicate
134) ! to other ranks
135) max_nface_per_cell = 0
136) if (option%myrank == option%io_rank) then
137)
138) allocate(temp_real_array(7,num_cells_local_save+1))
139) ! read for other processors
140) do irank = 0, option%mycommsize-1
141) temp_real_array = UNINITIALIZED_DOUBLE
142) num_to_read = num_cells_local_save
143) if (irank < remainder) num_to_read = num_to_read + 1
144) do icell = 1, num_to_read
145) call InputReadPflotranString(input,option)
146) call InputReadStringErrorMsg(input,option,hint)
147)
148) call InputReadInt(input,option,temp_int)
149) call InputErrorMsg(input,option,'cell id',hint)
150) temp_real_array(1,icell) = dble(temp_int)
151)
152) call InputReadInt(input,option,temp_int)
153) call InputErrorMsg(input,option,'num faces',hint)
154) temp_real_array(2,icell) = dble(temp_int)
155) nfaces_per_proc(irank+1) = nfaces_per_proc(irank+1) + &
156) temp_int
157) if (temp_int > max_nface_per_cell) &
158) max_nface_per_cell = temp_int
159)
160) call InputReadInt(input,option,temp_int)
161) call InputErrorMsg(input,option,'num vertices',hint)
162) temp_real_array(3,icell) = dble(temp_int)
163) if (temp_int>max_nvert_per_cell) max_nvert_per_cell = temp_int
164)
165) call InputReadDouble(input,option,temp_real_array(4,icell))
166) call InputErrorMsg(input,option,'cell x coordinate',hint)
167)
168) call InputReadDouble(input,option,temp_real_array(5,icell))
169) call InputErrorMsg(input,option,'cell y coordinate',hint)
170)
171) call InputReadDouble(input,option,temp_real_array(6,icell))
172) call InputErrorMsg(input,option,'cell z coordinate',hint)
173)
174) call InputReadDouble(input,option,temp_real_array(7,icell))
175) call InputErrorMsg(input,option,'cell volume',hint)
176) enddo
177) if (nfaces_per_proc(irank+1)>num_faces_local_save) &
178) num_faces_local_save = nfaces_per_proc(irank+1)
179)
180) if (irank == option%io_rank) then
181) ! cells reside on io_rank
182) num_faces_local = 0
183) pgrid%num_cells_local = num_cells_local
184) do icell = 1, num_cells_local
185) pgrid%cell_ids(icell) = int(temp_real_array(1,icell))
186) pgrid%cell_nfaces(icell) = int(temp_real_array(2,icell))
187) pgrid%cell_nverts(icell) = int(temp_real_array(3,icell))
188) pgrid%cell_centroids(icell)%x = temp_real_array(4,icell)
189) pgrid%cell_centroids(icell)%y = temp_real_array(5,icell)
190) pgrid%cell_centroids(icell)%z = temp_real_array(6,icell)
191) pgrid%cell_volumes(icell) = temp_real_array(7,icell)
192)
193) num_faces_local = num_faces_local + &
194) pgrid%cell_nfaces(icell)
195) enddo
196) pgrid%num_faces_local = num_faces_local
197) else
198) ! otherwise communicate to other ranks
199) int_mpi = num_to_read*7
200) call MPI_Send(temp_real_array,int_mpi,MPI_DOUBLE_PRECISION,irank, &
201) num_to_read,option%mycomm,ierr)
202) endif
203) enddo
204)
205) else
206)
207) ! other ranks pos the recv
208) allocate(temp_real_array(7,num_cells_local))
209) int_mpi = num_cells_local*7
210) call MPI_Recv(temp_real_array,int_mpi, &
211) MPI_DOUBLE_PRECISION,option%io_rank, &
212) MPI_ANY_TAG,option%mycomm,status_mpi,ierr)
213) num_faces_local = 0
214) pgrid%num_cells_local = num_cells_local
215) do icell = 1, num_cells_local
216) pgrid%cell_ids(icell) = int(temp_real_array(1,icell))
217) pgrid%cell_nfaces(icell) = int(temp_real_array(2,icell))
218) pgrid%cell_nverts(icell) = int(temp_real_array(3,icell))
219) pgrid%cell_centroids(icell)%x = temp_real_array(4,icell)
220) pgrid%cell_centroids(icell)%y = temp_real_array(5,icell)
221) pgrid%cell_centroids(icell)%z = temp_real_array(6,icell)
222) pgrid%cell_volumes(icell) = temp_real_array(7,icell)
223)
224) num_faces_local = num_faces_local + &
225) pgrid%cell_nfaces(icell)
226) enddo
227) num_faces_local_save = num_faces_local
228) pgrid%num_faces_local = num_faces_local
229) endif
230) deallocate(temp_real_array)
231)
232) call MPI_Bcast(max_nvert_per_cell,ONE_INTEGER_MPI,MPI_INTEGER,option%io_rank, &
233) option%mycomm,ierr)
234) ugrid%max_nvert_per_cell = max_nvert_per_cell
235) pgrid%max_nvert_per_cell = max_nvert_per_cell
236) allocate(pgrid%cell_vertids(max_nvert_per_cell,num_cells_local))
237)
238) if (option%myrank == option%io_rank) then
239)
240) call InputReadPflotranString(input,option)
241) call InputReadWord(input,option,card,PETSC_TRUE)
242) word = 'FACES'
243) call InputErrorMsg(input,option,word,card)
244) if (.not.StringCompare(word,card)) then
245) option%io_buffer = 'Unrecongnized keyword "' // trim(card) // &
246) '" in polyhedra grid file.'
247) call printErrMsgByRank(option)
248) endif
249)
250) hint = 'Polyhedra Unstructured Grid FACES'
251) call InputReadInt(input,option,temp_int)
252) call InputErrorMsg(input,option,'number of faces',hint)
253) endif
254)
255) int_mpi = 1
256) call MPI_Bcast(temp_int,ONE_INTEGER_MPI,MPI_INTEGER,option%io_rank, &
257) option%mycomm,ierr)
258) num_faces = temp_int
259) pgrid%num_faces_global = num_faces
260)
261) ! divide faces across ranks
262) allocate(pgrid%face_ids(num_faces_local))
263) allocate(pgrid%face_cellids(num_faces_local))
264) allocate(pgrid%face_nverts(num_faces_local))
265) allocate(pgrid%face_vertids(max_nvert_per_cell,num_faces_local))
266) allocate(pgrid%face_areas(num_faces_local))
267) allocate(pgrid%face_centroids(num_faces_local))
268) do iface = 1, num_faces_local
269) pgrid%face_centroids(iface)%x = 0.0d0
270) pgrid%face_centroids(iface)%y = 0.0d0
271) pgrid%face_centroids(iface)%z = 0.0d0
272) enddo
273)
274) ! read all faces from ASCII file through io_rank and communicate
275) ! to other ranks
276) max_nvert_per_face = 0
277) if (option%myrank == option%io_rank) then
278) allocate(temp_real_array(7+max_nvert_per_cell,num_faces_local_save))
279) do irank = 0, option%mycommsize-1
280) temp_real_array = UNINITIALIZED_DOUBLE
281) num_to_read = nfaces_per_proc(irank+1)
282) do iface = 1, num_to_read
283) call InputReadPflotranString(input,option)
284) call InputReadStringErrorMsg(input,option,hint)
285)
286) call InputReadInt(input,option,temp_int)
287) call InputErrorMsg(input,option,'face id',hint)
288) temp_real_array(1,iface) = dble(temp_int)
289)
290) call InputReadInt(input,option,temp_int)
291) call InputErrorMsg(input,option,'cell id',hint)
292) temp_real_array(2,iface) = dble(temp_int)
293)
294) call InputReadInt(input,option,temp_int)
295) call InputErrorMsg(input,option,'number of vertices',hint)
296) temp_real_array(3,iface) = dble(temp_int)
297) if (temp_int > max_nvert_per_face) &
298) max_nvert_per_face = temp_int
299)
300) do ivert = 1, int(temp_real_array(3,iface))
301) call InputReadInt(input,option,temp_int)
302) call InputErrorMsg(input,option,'face - vertex id',hint)
303) temp_real_array(3+ivert,iface) = dble(temp_int)
304) enddo
305)
306) call InputReadDouble(input,option,temp_real_array(4+max_nvert_per_cell,iface))
307) call InputErrorMsg(input,option,'face x coordinate',hint)
308) call InputReadDouble(input,option,temp_real_array(5+max_nvert_per_cell,iface))
309) call InputErrorMsg(input,option,'face y coordinate',hint)
310) call InputReadDouble(input,option,temp_real_array(6+max_nvert_per_cell,iface))
311) call InputErrorMsg(input,option,'face z coordinate',hint)
312) call InputReadDouble(input,option,temp_real_array(7+max_nvert_per_cell,iface))
313) call InputErrorMsg(input,option,'face area',hint)
314)
315) enddo
316)
317) if (irank == option%io_rank) then
318) do iface = 1, num_faces_local
319) pgrid%face_ids(iface) = int(temp_real_array(1,iface))
320) pgrid%face_cellids(iface) = int(temp_real_array(2,iface))
321) pgrid%face_nverts(iface) = int(temp_real_array(3,iface))
322) pgrid%face_vertids(:,iface) = UNINITIALIZED_INTEGER
323) do ivert = 1, pgrid%face_nverts(iface)
324) pgrid%face_vertids(ivert,iface) = &
325) int(temp_real_array(3+ivert,iface))
326) enddo
327) pgrid%face_centroids(iface)%x = &
328) temp_real_array(4+max_nvert_per_cell,iface)
329) pgrid%face_centroids(iface)%y = &
330) temp_real_array(5+max_nvert_per_cell,iface)
331) pgrid%face_centroids(iface)%z = &
332) temp_real_array(6+max_nvert_per_cell,iface)
333) pgrid%face_areas(iface) = temp_real_array(7+max_nvert_per_cell,iface)
334) enddo
335) else
336)
337) int_mpi = num_to_read*(7+max_nvert_per_cell)
338) call MPI_Send(temp_real_array,int_mpi,MPI_DOUBLE_PRECISION,irank, &
339) num_to_read,option%mycomm,ierr)
340) endif
341)
342) enddo
343) else
344) ! other ranks post the recv
345) allocate(temp_real_array(7+max_nvert_per_cell,num_faces_local+1))
346) int_mpi = num_faces_local*(7+max_nvert_per_cell)
347) call MPI_Recv(temp_real_array,int_mpi, &
348) MPI_DOUBLE_PRECISION,option%io_rank, &
349) MPI_ANY_TAG,option%mycomm,status_mpi,ierr)
350)
351) do iface = 1, num_faces_local
352) pgrid%face_ids(iface) = int(temp_real_array(1,iface))
353) pgrid%face_cellids(iface) = int(temp_real_array(2,iface))
354) pgrid%face_nverts(iface) = int(temp_real_array(3,iface))
355) pgrid%face_vertids(:,iface) = UNINITIALIZED_INTEGER
356) do ivert = 1, pgrid%face_nverts(iface)
357) pgrid%face_vertids(ivert,iface) = &
358) int(temp_real_array(3+ivert,iface))
359) enddo
360) pgrid%face_centroids(iface)%x = &
361) temp_real_array(4+max_nvert_per_cell,iface)
362) pgrid%face_centroids(iface)%y = &
363) temp_real_array(5+max_nvert_per_cell,iface)
364) pgrid%face_centroids(iface)%z = &
365) temp_real_array(6+max_nvert_per_cell,iface)
366) pgrid%face_areas(iface) = temp_real_array(7+max_nvert_per_cell,iface)
367)
368) enddo
369) endif
370) deallocate(temp_real_array)
371)
372) if (option%myrank == option%io_rank) then
373) call InputReadPflotranString(input,option)
374) call InputReadWord(input,option,card,PETSC_TRUE)
375) word = 'VERTICES'
376) call InputErrorMsg(input,option,word,card)
377) if (.not.StringCompare(word,card)) then
378) option%io_buffer = 'Unrecongnized keyword "' // trim(card) // &
379) '" in polyhedra grid file.'
380) call printErrMsgByRank(option)
381) endif
382) hint = 'Polyhedra Unstructured Grid VERTICES'
383) call InputReadInt(input,option,temp_int)
384) call InputErrorMsg(input,option,'number of vertices',hint)
385) endif
386)
387) int_mpi = 1
388) call MPI_Bcast(temp_int,ONE_INTEGER_MPI,MPI_INTEGER,option%io_rank, &
389) option%mycomm,ierr)
390) num_vertices = temp_int
391) pgrid%num_vertices_global = num_vertices
392) ugrid%num_vertices_global = num_vertices
393)
394) ! divide cells across ranks
395) num_vertices_local = num_vertices/option%mycommsize
396) num_vertices_local_save = num_vertices_local
397) remainder = num_vertices - &
398) num_vertices_local*option%mycommsize
399) if (option%myrank < remainder) num_vertices_local = &
400) num_vertices_local + 1
401)
402) allocate(pgrid%vertex_coordinates(num_vertices_local))
403) do ivert = 1, num_vertices_local
404) pgrid%vertex_coordinates(ivert)%x = 0.d0
405) pgrid%vertex_coordinates(ivert)%y = 0.d0
406) pgrid%vertex_coordinates(ivert)%z = 0.d0
407) enddo
408)
409) ! read all vertices from ASCII file through io_rank and communicate
410) ! to other ranks
411) if (option%myrank == option%io_rank) then
412) allocate(temp_real_array(3,num_vertices_local_save+1))
413) ! read for all processors
414) do irank = 0, option%mycommsize-1
415) temp_real_array = UNINITIALIZED_DOUBLE
416) num_to_read = num_vertices_local_save
417) if (irank < remainder) num_to_read = num_to_read + 1
418) do ivert = 1, num_to_read
419) call InputReadPflotranString(input,option)
420) call InputReadStringErrorMsg(input,option,hint)
421) call InputReadDouble(input,option,temp_real_array(1,ivert))
422) call InputErrorMsg(input,option,'vertex x coordinate',hint)
423) call InputReadDouble(input,option,temp_real_array(2,ivert))
424) call InputErrorMsg(input,option,'vertex y coordinate',hint)
425) call InputReadDouble(input,option,temp_real_array(3,ivert))
426) call InputErrorMsg(input,option,'vertex z coordinate',hint)
427) enddo
428)
429) if (irank == option%io_rank) then
430) pgrid%num_vertices_local = num_vertices_local
431) do ivert = 1, num_vertices_local
432) pgrid%vertex_coordinates(ivert)%x = temp_real_array(1,ivert)
433) pgrid%vertex_coordinates(ivert)%y = temp_real_array(2,ivert)
434) pgrid%vertex_coordinates(ivert)%z = temp_real_array(3,ivert)
435) enddo
436) else
437) int_mpi = num_to_read*3
438) call MPI_Send(temp_real_array,int_mpi,MPI_DOUBLE_PRECISION,irank, &
439) num_to_read,option%mycomm,ierr)
440) endif
441)
442) enddo
443) else
444)
445) pgrid%num_vertices_local = num_vertices_local
446) allocate(temp_real_array(3,num_vertices_local))
447) int_mpi = num_vertices_local*3
448) call MPI_Recv(temp_real_array,int_mpi, &
449) MPI_DOUBLE_PRECISION,option%io_rank, &
450) MPI_ANY_TAG,option%mycomm,status_mpi,ierr)
451) do ivert = 1, num_vertices_local
452) pgrid%vertex_coordinates(ivert)%x = temp_real_array(1,ivert)
453) pgrid%vertex_coordinates(ivert)%y = temp_real_array(2,ivert)
454) pgrid%vertex_coordinates(ivert)%z = temp_real_array(3,ivert)
455) enddo
456) endif
457)
458) deallocate(temp_real_array)
459) deallocate(nfaces_per_proc)
460)
461) temp_int = max_nface_per_cell
462) call MPI_Bcast(temp_int,ONE_INTEGER_MPI,MPI_INTEGER,option%io_rank, &
463) option%mycomm,ierr)
464) pgrid%max_nface_per_cell = temp_int
465)
466) temp_int = max_nvert_per_face
467) call MPI_Bcast(temp_int,ONE_INTEGER_MPI,MPI_INTEGER,option%io_rank, &
468) option%mycomm,ierr)
469) pgrid%max_nvert_per_face = temp_int
470)
471) if (option%myrank == option%io_rank) then
472) call InputDestroy(input)
473) endif
474)
475) end subroutine UGridPolyhedraRead
476)
477) ! ************************************************************************** !
478)
479) subroutine UGridPolyhedraDecompose(ugrid, option)
480) !
481) ! This routine decomposes a polyhedra grid across ranks.
482) !
483) ! Author: Gautam Bisht, LBL
484) ! Date: 09/29/13
485) !
486)
487) use Input_Aux_module
488) use Option_module
489) use String_module
490) use Grid_Unstructured_Cell_module
491) use Utility_module, only: reallocateIntArray, SearchOrderedArray
492)
493) implicit none
494)
495) #include "petsc/finclude/petscvec.h"
496) #include "petsc/finclude/petscvec.h90"
497) #include "petsc/finclude/petscmat.h"
498) #include "petsc/finclude/petscmat.h90"
499) #include "petsc/finclude/petscdm.h"
500) #include "petsc/finclude/petscdm.h90"
501) #include "petsc/finclude/petscis.h"
502) #include "petsc/finclude/petscis.h90"
503) #include "petsc/finclude/petscviewer.h"
504)
505) type(grid_unstructured_type) :: ugrid
506) type(option_type) :: option
507)
508) type(unstructured_polyhedra_type), pointer :: pgrid
509) PetscInt :: icell
510) PetscInt :: iface
511) PetscInt :: ivertex
512) PetscInt :: ivertex2
513) PetscInt :: iface_beg
514) PetscInt :: iface_end
515) PetscInt :: local_id
516) PetscInt :: min_vertex_id
517) PetscInt :: index_format_flag
518) PetscInt :: num_cells_local_old
519) PetscInt :: global_offset_old
520) PetscInt :: num_common_vertices
521) PetscInt :: num_cells_local_new
522) PetscInt :: count
523) PetscInt :: max_nvert_per_cell
524) PetscInt :: max_ndual_per_cell
525) PetscInt :: max_nvert_per_face
526) PetscInt :: max_nface_per_cell
527) PetscInt :: vertex_ids_offset
528) PetscInt :: temp_int
529) PetscInt :: dual_offset
530) PetscInt :: face_offset
531) PetscInt :: natural_id_offset
532) PetscInt :: cell_stride
533) PetscInt :: face_count
534) PetscInt :: max_int_count
535) PetscInt :: vertex_count
536) PetscInt :: vertex_id
537) PetscInt :: ghosted_id
538) PetscInt :: nface
539) PetscInt :: iface_vert
540) PetscInt :: num_faces_local
541) PetscInt :: idx
542) PetscInt, allocatable :: local_vertices(:)
543) PetscInt, allocatable :: local_vertex_offset(:)
544) PetscInt, allocatable :: int_array(:)
545) PetscInt, allocatable :: int_array1(:)
546) PetscInt, allocatable :: int_array2(:)
547) PetscInt, allocatable :: int_array3(:)
548) PetscInt, allocatable :: int_array4(:)
549) PetscInt, allocatable :: needed_vertices_petsc(:)
550) PetscInt, allocatable :: vert_n2g(:,:)
551) PetscInt, pointer :: int_array_pointer(:)
552) PetscErrorCode :: ierr
553) PetscInt, pointer :: ia_ptr(:), ja_ptr(:)
554) PetscReal, pointer :: vec_ptr(:)
555) PetscBool :: success
556) PetscInt :: num_rows, num_cols, istart, iend, icol
557) character(len=MAXSTRINGLENGTH) :: string
558)
559) PetscViewer :: viewer
560) Mat :: Adj_mat
561) Mat :: Dual_mat
562) MatPartitioning :: Part
563) Vec :: elements_natural
564) Vec :: elements_local
565) Vec :: elements_old
566) Vec :: vertices_old
567) Vec :: vertices_new
568) IS :: is_new
569) IS :: is_scatter
570) IS :: is_gather
571)
572) VecScatter :: vec_scatter
573)
574) pgrid => ugrid%polyhedra_grid
575) max_nvert_per_cell = ugrid%max_nvert_per_cell
576)
577) pgrid%cell_vertids = UNINITIALIZED_INTEGER
578)
579) min_vertex_id = 2 ! min value should be either 0 or 1 after global
580) ! reduction to idenitfy if vertex ids is 0-based
581) ! or 1-based.
582) iface_beg = 1
583) iface_end = 0
584) ! find vertex ids forming a cell
585) do icell = 1, pgrid%num_cells_local
586) iface_end = iface_beg + pgrid%cell_nfaces(icell) - 1
587) do iface = iface_beg, iface_end
588) if (pgrid%face_cellids(iface) /= &
589) pgrid%cell_ids(icell)) then
590) option%io_buffer = 'Face ID does not correspond to cell'
591) call printErrMsgByRank(option)
592) endif
593) do ivertex = 1, pgrid%face_nverts(iface)
594) do ivertex2 = 1, max_nvert_per_cell
595) if (pgrid%cell_vertids(ivertex2,icell) == UNINITIALIZED_INTEGER) then
596) pgrid%cell_vertids(ivertex2,icell) = &
597) pgrid%face_vertids(ivertex,iface)
598)
599) if (pgrid%cell_vertids(ivertex2,icell) < min_vertex_id) &
600) min_vertex_id = pgrid%cell_vertids(ivertex2,icell)
601) exit
602) endif
603) if (pgrid%cell_vertids(ivertex2,icell) == &
604) pgrid%face_vertids(ivertex,iface)) exit
605) enddo
606) enddo
607)
608) enddo
609) iface_beg = iface_end + 1
610) enddo
611)
612) allocate(int_array1(max_nvert_per_cell))
613) allocate(int_array2(max_nvert_per_cell))
614)
615) ! for a given cell, sort vertices in ascending order
616) do icell = 1, pgrid%num_cells_local
617) do ivertex = 1, pgrid%cell_nverts(icell)
618) int_array1(ivertex) = pgrid%cell_vertids(ivertex,icell)
619) int_array2(ivertex) = ivertex - 1
620) enddo
621) call PetscSortIntWithPermutation(pgrid%cell_nverts(icell), &
622) int_array1,int_array2,ierr);CHKERRQ(ierr)
623) int_array2 = int_array2 + 1
624) do ivertex = 1, pgrid%cell_nverts(icell)
625) pgrid%cell_vertids(ivertex,icell) = &
626) int_array1(int_array2(ivertex))
627) enddo
628) enddo
629) deallocate(int_array1)
630) deallocate(int_array2)
631)
632) call MPI_Allreduce(min_vertex_id,index_format_flag, &
633) ONE_INTEGER_MPI,MPIU_INTEGER,MPI_MIN,option%mycomm,ierr)
634)
635) if (index_format_flag /= 0 .and. index_format_flag /= 1) then
636) call printErrMsg(option,'Min. vertex id is neither 0 nor 1. Check input mesh.')
637) endif
638)
639) num_cells_local_old = pgrid%num_cells_local
640) ! let's make it Fortran indexing (i.e. 1-based)
641) do local_id = 1, num_cells_local_old
642) do ivertex = 1, max_nvert_per_cell
643) ! at this point we may be zero-based
644) if (pgrid%cell_vertids(ivertex,local_id) < 0) then
645) ! change no_value (UNINITIALIZED_INTEGER) to '0'
646) pgrid%cell_vertids(ivertex,local_id) = 0
647) else
648) if (index_format_flag == 0) then
649) ! let's make it Fortran indexing
650) pgrid%cell_vertids(ivertex,local_id) = &
651) pgrid%cell_vertids(ivertex,local_id) + 1
652) endif
653) endif
654) enddo
655) enddo
656)
657) #if UGRID_DEBUG
658) write(string,*) ugrid%max_nvert_per_cell
659) option%io_buffer = 'Maximum number of vertices per cell: ' // adjustl(string)
660) call printMsg(option)
661) write(string,*) index_format_flag
662) option%io_buffer = 'Vertex indexing starts at: ' // adjustl(string)
663) call printMsg(option)
664) if (index_format_flag == 0) then
665) option%io_buffer = 'Changing vertex indexing to 1-based.'
666) call printMsg(option)
667) endif
668) #endif
669)
670) num_cells_local_old = pgrid%num_cells_local
671) allocate(local_vertices(max_nvert_per_cell*num_cells_local_old))
672) allocate(local_vertex_offset(num_cells_local_old+1))
673) local_vertices = 0
674) local_vertex_offset = 0
675) count = 0
676) local_vertex_offset(1) = 0
677) do local_id = 1, num_cells_local_old
678) do ivertex = 1, ugrid%max_nvert_per_cell
679) if (pgrid%cell_vertids(ivertex,local_id) == 0) exit
680) count = count + 1
681) ! local vertices must be zero-based for MatCreateMPIAdj; thus subtract 1
682) local_vertices(count) = &
683) pgrid%cell_vertids(ivertex,local_id) - 1
684) enddo
685) local_vertex_offset(local_id+1) = count
686) enddo
687)
688) select case (ugrid%grid_type)
689) case(TWO_DIM_GRID)
690) num_common_vertices = 2 ! cells must share at least this number of vertices
691) case(THREE_DIM_GRID)
692) num_common_vertices = 3 ! cells must share at least this number of vertices
693) case default
694) option%io_buffer = 'Grid type not recognized '
695) call printErrMsg(option)
696) end select
697)
698) ! determine the global offset from 0 for cells on this rank
699) global_offset_old = 0
700) call MPI_Exscan(num_cells_local_old,global_offset_old, &
701) ONE_INTEGER_MPI,MPIU_INTEGER,MPI_SUM,option%mycomm,ierr)
702)
703) ! create an adjacency matrix for calculating the duals (connnections)
704) #if UGRID_DEBUG
705) call printMsg(option,'Adjacency matrix')
706) #endif
707)
708) call MatCreateMPIAdj(option%mycomm,num_cells_local_old, &
709) pgrid%num_vertices_global, &
710) local_vertex_offset, &
711) local_vertices,PETSC_NULL_INTEGER,Adj_mat, &
712) ierr);CHKERRQ(ierr)
713)
714) ! do not free local_vertices; MatAdjDestroy will do it
715) ! do not free local_vertex_offset; MatAdjDestroy will do it
716)
717) #if UGRID_DEBUG
718) if (ugrid%grid_type == THREE_DIM_GRID) then
719) call PetscViewerASCIIOpen(option%mycomm,'Adj_subsurf.out',viewer, &
720) ierr);CHKERRQ(ierr)
721) else
722) call PetscViewerASCIIOpen(option%mycomm,'Adj_surf.out',viewer, &
723) ierr);CHKERRQ(ierr)
724) endif
725) call MatView(Adj_mat,viewer,ierr);CHKERRQ(ierr)
726) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
727) #endif
728)
729) #if UGRID_DEBUG
730) call printMsg(option,'Dual matrix')
731) #endif
732)
733) #if defined(PETSC_HAVE_PARMETIS)
734) call MatMeshToCellGraph(Adj_mat,num_common_vertices,Dual_mat, &
735) ierr);CHKERRQ(ierr)
736) #endif
737) call MatDestroy(Adj_mat,ierr);CHKERRQ(ierr)
738)
739) #if UGRID_DEBUG
740) if (ugrid%grid_type == THREE_DIM_GRID) then
741) call PetscViewerASCIIOpen(option%mycomm,'Dual_subsurf.out',viewer, &
742) ierr);CHKERRQ(ierr)
743) else
744) call PetscViewerASCIIOpen(option%mycomm,'Dual_surf.out',viewer, &
745) ierr);CHKERRQ(ierr)
746) endif
747) call MatView(Dual_mat,viewer,ierr);CHKERRQ(ierr)
748) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
749) #endif
750)
751) call UGridPartition(ugrid,option,Dual_mat,is_new, &
752) num_cells_local_new)
753)
754) if (allocated(local_vertices)) deallocate(local_vertices)
755) if (allocated(local_vertex_offset)) deallocate(local_vertex_offset)
756)
757) ! second argument of ZERO_INTEGER means to use 0-based indexing
758) ! MagGetRowIJF90 returns row and column pointers for compressed matrix data
759) call MatGetRowIJF90(Dual_mat,ZERO_INTEGER,PETSC_FALSE,PETSC_FALSE,num_rows, &
760) ia_ptr,ja_ptr,success,ierr);CHKERRQ(ierr)
761)
762) if (.not.success .or. num_rows /= num_cells_local_old) then
763) print *, option%myrank, num_rows, success, num_cells_local_old
764) option%io_buffer = 'Error getting IJ row indices from dual matrix'
765) call printErrMsg(option)
766) endif
767)
768)
769) if (.not.success .or. num_rows /= num_cells_local_old) then
770) print *, option%myrank, num_rows, success, num_cells_local_old
771) option%io_buffer = 'Error getting IJ row indices from dual matrix'
772) call printErrMsg(option)
773) endif
774)
775) ! calculate maximum number of connections for any given cell
776) max_ndual_per_cell = 0
777) do local_id = 1, num_cells_local_old
778) istart = ia_ptr(local_id)
779) iend = ia_ptr(local_id+1)-1
780) num_cols = iend-istart+1
781) if (num_cols > max_ndual_per_cell) &
782) max_ndual_per_cell = num_cols
783) enddo
784) temp_int = max_ndual_per_cell
785) call MPI_Allreduce(temp_int,max_ndual_per_cell, &
786) ONE_INTEGER_MPI,MPIU_INTEGER,MPI_MAX,option%mycomm,ierr)
787) ugrid%max_ndual_per_cell = max_ndual_per_cell
788)
789) #if UGRID_DEBUG
790) write(string,*) max_ndual_per_cell
791) option%io_buffer = 'Maximum number of duals per cell: ' // adjustl(string)
792) call printMsg(option)
793) #endif
794)
795) call MatRestoreRowIJF90(Dual_mat,ZERO_INTEGER,PETSC_FALSE,PETSC_FALSE, &
796) num_rows,ia_ptr,ja_ptr,success,ierr);CHKERRQ(ierr)
797)
798) ! in order to redistributed vertex/cell data among ranks, I package it
799) ! in a crude way within a strided petsc vec and pass it. The stride
800) ! determines the size of each cells "packaged" data
801) max_nface_per_cell = pgrid%max_nface_per_cell
802) max_nvert_per_face = pgrid%max_nvert_per_face
803)
804) face_offset = 8 + 1 ! +1 for -666
805) vertex_ids_offset = face_offset + &
806) ( &
807) 1 + & ! natural face id
808) 1 + & ! natural cell id
809) 1 + & ! number of vertices
810) max_nvert_per_face + & ! maximum vertices forming a face
811) 3 + & ! face centroid (xc,yc,zc)
812) 1 & ! face area
813) ) * max_nface_per_cell + &
814) 1 ! for -777
815) dual_offset = vertex_ids_offset + max_nvert_per_cell + 1 ! +1 for -888
816) cell_stride = dual_offset+ max_ndual_per_cell + 1 ! +1 for -999999
817) natural_id_offset = 2
818)
819) ! Information for each cell is packed in a strided petsc vec
820) ! The information is ordered within each stride as follows:
821) ! -cell_N ! global cell id (negative indicates 1-based)
822) ! natural cell id
823) ! number of faces
824) ! number of vertices
825) ! cell x coordinate
826) ! cell y coordinate
827) ! cell z coordinate
828) ! cell volume
829) ! -666 ! separator between cell info and face info
830) ! face1_id
831) ! face1_cellid
832) ! face1_nverts
833) ! face1_vert1
834) ! face1_vert2
835) ! ...
836) ! face1_vertN
837) ! face1_xc
838) ! face1_yc
839) ! face1_zc
840) ! face1_area
841) ! face2
842) ! face2_...
843) ! ...
844) ! face2_area
845) ! ...
846) ! ...
847) ! faceM_area
848) ! -777 ! separator between face info and vertices
849) ! vertex1 ! in cell_N
850) ! vertex2
851) ! ...
852) ! vertexN
853) ! -888 ! separator between vertex and dual ids
854) ! dual1 ! dual ids between cell_N and others
855) ! dual2
856) ! ...
857) ! dualN
858) ! -999999 ! separator indicating end of information for cell_N
859)
860) ! the purpose of -777, -888, and -999999 is to allow one to use cells of
861) ! various geometry. Currently, the max # vertices = 8 and max # duals = 6.
862) ! But this will be generalized in the future.
863) call UGridCreateOldVec(ugrid,option, elements_old, &
864) num_cells_local_old, is_new, is_scatter, cell_stride)
865)
866) ! 0 = 0-based indexing
867) ! MagGetRowIJF90 returns row and column pointers for compressed matrix data
868) call MatGetRowIJF90(Dual_mat,ZERO_INTEGER,PETSC_FALSE,PETSC_FALSE,num_rows, &
869) ia_ptr,ja_ptr,success,ierr);CHKERRQ(ierr)
870)
871) call VecGetArrayF90(elements_old,vec_ptr,ierr);CHKERRQ(ierr)
872) count = 0
873) face_count = 0
874) do local_id = 1, num_cells_local_old
875) count = count + 1
876) vec_ptr(count) = -(global_offset_old+local_id)
877) count = count + 1
878) vec_ptr(count) = pgrid%cell_ids(local_id)
879) count = count + 1
880) vec_ptr(count) = pgrid%cell_nfaces(local_id)
881) count = count + 1
882) vec_ptr(count) = pgrid%cell_nverts(local_id)
883) count = count + 1
884) vec_ptr(count) = pgrid%cell_centroids(local_id)%x
885) count = count + 1
886) vec_ptr(count) = pgrid%cell_centroids(local_id)%y
887) count = count + 1
888) vec_ptr(count) = pgrid%cell_centroids(local_id)%z
889) count = count + 1
890) vec_ptr(count) = pgrid%cell_volumes(local_id)
891) count = count + 1
892) vec_ptr(count) = -666
893)
894) do iface = 1, pgrid%cell_nfaces(local_id)
895) face_count = face_count + 1
896) count = count + 1
897) vec_ptr(count) = pgrid%face_ids(face_count)
898) count = count + 1
899) vec_ptr(count) = pgrid%face_cellids(face_count)
900) count = count + 1
901) vec_ptr(count) = pgrid%face_nverts(face_count)
902)
903) do ivertex = 1, pgrid%face_nverts(face_count)
904) count = count + 1
905) vec_ptr(count) = pgrid%face_vertids(ivertex,face_count)
906) enddo
907)
908) do ivertex = pgrid%face_nverts(face_count)+1, max_nvert_per_face
909) count = count + 1
910) vec_ptr(count) = 0
911) enddo
912)
913) count = count + 1
914) vec_ptr(count) = pgrid%face_centroids(face_count)%x
915) count = count + 1
916) vec_ptr(count) = pgrid%face_centroids(face_count)%y
917) count = count + 1
918) vec_ptr(count) = pgrid%face_centroids(face_count)%z
919) count = count + 1
920) vec_ptr(count) = pgrid%face_areas(face_count)
921)
922) enddo
923)
924) do iface = pgrid%cell_nfaces(local_id)+1, max_nface_per_cell
925)
926) count = count + 1
927) vec_ptr(count) = 0 ! face_id
928) count = count + 1
929) vec_ptr(count) = 0 ! face_cellid
930) count = count + 1
931) vec_ptr(count) = 0 ! face_nverts
932)
933) do ivertex = 1, max_nvert_per_face
934) count = count + 1
935) vec_ptr(count) = 0 ! face_vertids
936) enddo
937)
938) count = count + 1
939) vec_ptr(count) = 0 ! face_centoid_x
940) count = count + 1
941) vec_ptr(count) = 0 ! face_centoid_y
942) count = count + 1
943) vec_ptr(count) = 0 ! face_centoid_z
944) count = count + 1
945) vec_ptr(count) = 0 ! face_area
946)
947) enddo
948)
949) count = count + 1
950) vec_ptr(count) = -777
951)
952) do ivertex = 1, max_nvert_per_cell
953) count = count + 1
954) if (ivertex <= pgrid%cell_nverts(local_id)) then
955) vec_ptr(count) = pgrid%cell_vertids(ivertex,local_id)
956) else
957) vec_ptr(count) = 0
958) endif
959) enddo
960)
961) count = count + 1
962) vec_ptr(count) = -888
963)
964) ! add the dual ids
965) istart = ia_ptr(local_id)
966) iend = ia_ptr(local_id+1)-1
967) num_cols = iend-istart+1
968) if (num_cols > max_ndual_per_cell) then
969) option%io_buffer = &
970) 'Number of columns in Dual matrix is larger then max_ndual_per_cell.'
971) call printErrMsgByRank(option)
972) endif
973) do icol = 1, max_ndual_per_cell
974) count = count + 1
975) if (icol <= num_cols) then
976) ! increment for 1-based ordering
977) vec_ptr(count) = ja_ptr(icol+istart) + 1
978) else
979) vec_ptr(count) = 0
980) endif
981) enddo
982) count = count + 1
983) ! final separator
984) vec_ptr(count) = -999999 ! help differentiate
985)
986) enddo
987) call VecRestoreArrayF90(elements_old,vec_ptr,ierr);CHKERRQ(ierr)
988)
989) call MatRestoreRowIJF90(Dual_mat, ZERO_INTEGER, PETSC_FALSE, PETSC_FALSE, &
990) num_rows, ia_ptr, ja_ptr, success, &
991) ierr);CHKERRQ(ierr)
992) call MatDestroy(Dual_mat,ierr);CHKERRQ(ierr)
993)
994) call UGridNaturalToPetsc(ugrid, option, &
995) elements_old, elements_local, &
996) num_cells_local_new, cell_stride, dual_offset, &
997) natural_id_offset, is_scatter)
998)
999) ! make a list of local vertices
1000) max_int_count = 2*ugrid%ngmax
1001) allocate(int_array_pointer(max_int_count))
1002) int_array_pointer = 0
1003) vertex_count = 0
1004) ! yep - load them all into a petsc vector
1005) ! note that the vertices are still in natural numbering
1006) call VecGetArrayF90(elements_local,vec_ptr,ierr);CHKERRQ(ierr)
1007) do ghosted_id=1, ugrid%ngmax
1008) do ivertex = 1, ugrid%max_nvert_per_cell
1009) vertex_id = int(vec_ptr(ivertex + vertex_ids_offset + (ghosted_id-1)*cell_stride))
1010) if (vertex_id < 1) exit
1011) vertex_count = vertex_count + 1
1012) if (vertex_count > max_int_count) then
1013) call reallocateIntArray(int_array_pointer,max_int_count)
1014) endif
1015) vec_ptr(ivertex + vertex_ids_offset + (ghosted_id-1)*cell_stride) = vertex_count
1016) int_array_pointer(vertex_count) = vertex_id
1017) enddo
1018) enddo
1019) call VecRestoreArrayF90(elements_local,vec_ptr,ierr);CHKERRQ(ierr)
1020)
1021) ! sort the vertex ids
1022) allocate(int_array(vertex_count))
1023) int_array(1:vertex_count) = int_array_pointer(1:vertex_count)
1024) allocate(int_array2(vertex_count))
1025) do ivertex = 1, vertex_count
1026) int_array2(ivertex) = ivertex
1027) enddo
1028) deallocate(int_array_pointer)
1029) nullify(int_array_pointer)
1030) int_array2 = int_array2-1
1031) call PetscSortIntWithPermutation(vertex_count,int_array,int_array2, &
1032) ierr);CHKERRQ(ierr)
1033) int_array2 = int_array2+1
1034)
1035) ! remove duplicates
1036) allocate(int_array3(vertex_count))
1037) allocate(int_array4(vertex_count))
1038) int_array3 = 0
1039) int_array4 = 0
1040) int_array3(1) = int_array(int_array2(1))
1041) count = 1
1042) int_array4(int_array2(1)) = count
1043) do ivertex = 2, vertex_count
1044) vertex_id = int_array(int_array2(ivertex))
1045) if (vertex_id > int_array3(count)) then
1046) count = count + 1
1047) int_array3(count) = vertex_id
1048) endif
1049) int_array4(int_array2(ivertex)) = count
1050) enddo
1051) vertex_count = count
1052)
1053) allocate(ugrid%vertex_ids_natural(vertex_count))
1054) ugrid%vertex_ids_natural = int_array3(1:vertex_count)
1055)
1056) ! now load all the vertices needed to define all the local cells
1057) ! on the processor
1058) allocate(needed_vertices_petsc(vertex_count))
1059) needed_vertices_petsc(1:vertex_count) = int_array3(1:vertex_count)
1060)
1061) ! allocate the array that will store the vertex ids for each cell.
1062) ! remember that max_nvert_per_cell is the max # of vertices in a cell
1063) ! currently hardwired to 8.
1064) !deallocate(ugrid%cell_vertices)
1065) allocate(ugrid%cell_vertices( &
1066) 0:ugrid%max_nvert_per_cell,ugrid%ngmax))
1067) ugrid%cell_vertices = 0
1068)
1069) ! permute the local ids calculated earlier in the int_array4
1070) allocate(vert_n2g(ugrid%max_nvert_per_cell,2))
1071)
1072) call VecGetArrayF90(elements_local,vec_ptr,ierr);CHKERRQ(ierr)
1073) do ghosted_id = 1, ugrid%ngmax
1074) vert_n2g = 0
1075) do ivertex = 1, ugrid%max_nvert_per_cell
1076) ! extract the original vertex id
1077) vertex_id = int(vec_ptr(ivertex + vertex_ids_offset + (ghosted_id-1)*cell_stride))
1078) if (vertex_id < 1) exit
1079) count = ugrid%cell_vertices(0,ghosted_id)+1
1080) ugrid%cell_vertices(count,ghosted_id) = &
1081) int_array4(vertex_id)
1082) ugrid%cell_vertices(0,ghosted_id) = count
1083) ! load the permuted value back into the petsc vector
1084) vec_ptr(ivertex + vertex_ids_offset + (ghosted_id-1)*cell_stride) = &
1085) int_array4(vertex_id)
1086)
1087) vert_n2g(ivertex,1) = int_array(vertex_id)
1088) vert_n2g(ivertex,2) = int_array4(vertex_id)
1089) enddo
1090)
1091) nface = int(vec_ptr((ghosted_id-1)*cell_stride+3))
1092) do iface = 1,nface
1093) do iface_vert = 1,max_nvert_per_face
1094) vertex_id = (ghosted_id-1)*cell_stride + &
1095) face_offset + &
1096) (iface-1)*(7 + max_nvert_per_face) + &
1097) 3 + iface_vert
1098) if (vec_ptr(vertex_id) < 1) exit
1099) do ivertex = 1, max_nvert_per_cell
1100) if (vec_ptr(vertex_id) == vert_n2g(ivertex,1)) then
1101) vec_ptr(vertex_id) = vert_n2g(ivertex,2)
1102) exit
1103) endif
1104) enddo
1105) enddo
1106) enddo
1107)
1108) enddo
1109) call VecRestoreArrayF90(elements_local,vec_ptr,ierr);CHKERRQ(ierr)
1110) deallocate(int_array)
1111) deallocate(int_array2)
1112) deallocate(int_array3)
1113) deallocate(int_array4)
1114)
1115) #if UGRID_DEBUG
1116) write(string,*) option%myrank
1117) if (ugrid%grid_type == THREE_DIM_GRID) then
1118) string = 'elements_vert_local' // trim(adjustl(string)) // '_subsurf.out'
1119) else
1120) string = 'elements_vert_local' // trim(adjustl(string)) // '_surf.out'
1121) endif
1122) call PetscViewerASCIIOpen(PETSC_COMM_SELF,trim(string),viewer, &
1123) ierr);CHKERRQ(ierr)
1124) call VecView(elements_local,viewer,ierr);CHKERRQ(ierr)
1125) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
1126) #endif
1127)
1128) deallocate(pgrid%cell_ids)
1129) deallocate(pgrid%cell_nfaces)
1130) deallocate(pgrid%cell_vertids)
1131) deallocate(pgrid%cell_nverts)
1132) deallocate(pgrid%cell_volumes)
1133) deallocate(pgrid%cell_centroids)
1134)
1135) deallocate(pgrid%face_ids)
1136) deallocate(pgrid%face_cellids)
1137) deallocate(pgrid%face_nverts)
1138) deallocate(pgrid%face_vertids)
1139) deallocate(pgrid%face_areas)
1140) deallocate(pgrid%face_centroids)
1141)
1142) allocate(pgrid%cell_ids(ugrid%ngmax))
1143) allocate(pgrid%cell_nfaces(ugrid%ngmax))
1144) allocate(pgrid%cell_vertids(max_nvert_per_cell,ugrid%ngmax))
1145) allocate(pgrid%cell_faceids(max_nface_per_cell,ugrid%ngmax))
1146) allocate(pgrid%cell_nverts(ugrid%ngmax))
1147) allocate(pgrid%cell_volumes(ugrid%ngmax))
1148) allocate(pgrid%cell_centroids(ugrid%ngmax))
1149)
1150) pgrid%cell_faceids = 0
1151) pgrid%cell_vertids = 0
1152)
1153) call VecGetArrayF90(elements_local,vec_ptr,ierr);CHKERRQ(ierr)
1154)
1155) num_faces_local = 0
1156) do ghosted_id = 1, ugrid%ngmax
1157) idx = (ghosted_id-1)*cell_stride
1158) pgrid%cell_ids(ghosted_id) = int(vec_ptr(idx + 2))
1159) pgrid%cell_nfaces(ghosted_id) = int(vec_ptr(idx + 3))
1160) pgrid%cell_nverts(ghosted_id) = int(vec_ptr(idx + 4))
1161) pgrid%cell_centroids(ghosted_id)%x = vec_ptr(idx + 5)
1162) pgrid%cell_centroids(ghosted_id)%y = vec_ptr(idx + 6)
1163) pgrid%cell_centroids(ghosted_id)%z = vec_ptr(idx + 7)
1164) pgrid%cell_volumes(ghosted_id) = vec_ptr(idx + 8)
1165)
1166) pgrid%cell_vertids(1:max_nvert_per_cell,ghosted_id) = &
1167) ugrid%cell_vertices(1:max_nvert_per_cell,ghosted_id)
1168)
1169) do iface = 1,pgrid%cell_nfaces(ghosted_id)
1170) num_faces_local = num_faces_local + 1
1171) pgrid%cell_faceids(iface,ghosted_id) = num_faces_local
1172) enddo
1173) enddo
1174)
1175) allocate(pgrid%face_ids(num_faces_local))
1176) allocate(pgrid%face_cellids(num_faces_local))
1177) allocate(pgrid%face_nverts(num_faces_local))
1178) allocate(pgrid%face_vertids(max_nvert_per_face,num_faces_local))
1179) allocate(pgrid%face_areas(num_faces_local))
1180) allocate(pgrid%face_centroids(num_faces_local))
1181) pgrid%num_faces_local = num_faces_local
1182) pgrid%face_vertids = 0
1183)
1184) count = 0
1185) do ghosted_id = 1,ugrid%ngmax
1186) idx = (ghosted_id-1)*cell_stride + face_offset
1187) do iface = 1,pgrid%cell_nfaces(ghosted_id)
1188) count = count + 1
1189) pgrid%face_ids(count) = int(vec_ptr(idx + 1))
1190) !pgrid%face_cellids(count) = int(vec_ptr(idx + 2))
1191) pgrid%face_cellids(count) = ghosted_id
1192) pgrid%face_nverts(count) = int(vec_ptr(idx + 3))
1193)
1194) do ivertex = 1,pgrid%face_nverts(count)
1195) pgrid%face_vertids(ivertex,count) = int(vec_ptr(idx + 3 + ivertex))
1196) enddo
1197)
1198) idx = idx + 3 + max_nvert_per_face
1199)
1200) pgrid%face_centroids(count)%x = vec_ptr(idx + 1)
1201) pgrid%face_centroids(count)%y = vec_ptr(idx + 2)
1202) pgrid%face_centroids(count)%z = vec_ptr(idx + 3)
1203) pgrid%face_areas(count) = vec_ptr(idx + 4)
1204)
1205) idx = idx + 4
1206) enddo
1207) enddo
1208)
1209) call VecRestoreArrayF90(elements_local,vec_ptr,ierr);CHKERRQ(ierr)
1210)
1211) call VecDestroy(elements_local,ierr);CHKERRQ(ierr)
1212)
1213) ! now we need to work on aligning the original vertex coordinates with
1214) ! the current ordering or permuted/rearranged ordering.
1215)
1216) ugrid%num_vertices_local = pgrid%num_vertices_local
1217)
1218) ! IS for gather operation - need local numbering
1219) allocate(int_array(vertex_count))
1220) ! vertex_count = # of local vertices (I believe ghosted+non-ghosted)
1221) do ivertex = 1, vertex_count
1222) int_array(ivertex) = ivertex-1
1223) enddo
1224)
1225) ! include cell ids (use block ids, not indices)
1226) call ISCreateBlock(option%mycomm,3,vertex_count, &
1227) int_array,PETSC_COPY_VALUES,is_gather,ierr);CHKERRQ(ierr)
1228) deallocate(int_array)
1229)
1230) ! create a parallel petsc vector with a stride of 3.
1231) call VecCreate(option%mycomm,vertices_old,ierr);CHKERRQ(ierr)
1232) call VecSetSizes(vertices_old,ugrid%num_vertices_local*3, &
1233) PETSC_DECIDE,ierr);CHKERRQ(ierr)
1234) call VecSetBlockSize(vertices_old,3,ierr);CHKERRQ(ierr)
1235) call VecSetFromOptions(vertices_old,ierr);CHKERRQ(ierr)
1236)
1237) ! create serial petsc vector with a stride of 3
1238) call VecCreate(PETSC_COMM_SELF,vertices_new,ierr);CHKERRQ(ierr)
1239) call VecSetSizes(vertices_new,vertex_count*3,PETSC_DECIDE, &
1240) ierr);CHKERRQ(ierr)
1241) call VecSetBlockSize(vertices_new,3,ierr);CHKERRQ(ierr)
1242) call VecSetFromOptions(vertices_new,ierr);CHKERRQ(ierr)
1243)
1244) call VecCreate(option%mycomm,vertices_old,ierr);CHKERRQ(ierr)
1245) call VecSetSizes(vertices_old,ugrid%num_vertices_local*3, &
1246) PETSC_DECIDE,ierr);CHKERRQ(ierr)
1247) call VecSetBlockSize(vertices_old,3,ierr);CHKERRQ(ierr)
1248) call VecSetFromOptions(vertices_old,ierr);CHKERRQ(ierr)
1249)
1250) ! create serial petsc vector with a stride of 3
1251) call VecCreate(PETSC_COMM_SELF,vertices_new,ierr);CHKERRQ(ierr)
1252) call VecSetSizes(vertices_new,vertex_count*3,PETSC_DECIDE, &
1253) ierr);CHKERRQ(ierr)
1254) call VecSetBlockSize(vertices_new,3,ierr);CHKERRQ(ierr)
1255) call VecSetFromOptions(vertices_new,ierr);CHKERRQ(ierr)
1256)
1257) ! load up the coordinates
1258) call VecGetArrayF90(vertices_old,vec_ptr,ierr);CHKERRQ(ierr)
1259) do ivertex = 1, pgrid%num_vertices_local
1260) vec_ptr((ivertex-1)*3+1) = pgrid%vertex_coordinates(ivertex)%x
1261) vec_ptr((ivertex-1)*3+2) = pgrid%vertex_coordinates(ivertex)%y
1262) vec_ptr((ivertex-1)*3+3) = pgrid%vertex_coordinates(ivertex)%z
1263) enddo
1264) call VecRestoreArrayF90(vertices_old,vec_ptr,ierr);CHKERRQ(ierr)
1265) deallocate(pgrid%vertex_coordinates)
1266) nullify(pgrid%vertex_coordinates)
1267)
1268) ! IS for scatter - provide petsc global numbering
1269) allocate(int_array(vertex_count))
1270) do ivertex = 1, vertex_count
1271) int_array(ivertex) = (needed_vertices_petsc(ivertex)-1)
1272) enddo
1273) ! include cell ids
1274) call ISCreateBlock(option%mycomm,3,vertex_count, &
1275) int_array,PETSC_COPY_VALUES,is_scatter, &
1276) ierr);CHKERRQ(ierr)
1277) deallocate(int_array)
1278)
1279) ! resize vertex array to new size
1280) ugrid%num_vertices_natural = ugrid%num_vertices_local
1281) ugrid%num_vertices_local = vertex_count
1282) allocate(ugrid%vertices(vertex_count))
1283) do ivertex = 1, vertex_count
1284) ugrid%vertices(ivertex)%x = 0.d0
1285) ugrid%vertices(ivertex)%y = 0.d0
1286) ugrid%vertices(ivertex)%z = 0.d0
1287) enddo
1288) pgrid%num_vertices_local = vertex_count
1289) allocate(pgrid%vertex_coordinates(vertex_count))
1290) do ivertex = 1, vertex_count
1291) pgrid%vertex_coordinates(ivertex)%x = 0.d0
1292) pgrid%vertex_coordinates(ivertex)%y = 0.d0
1293) pgrid%vertex_coordinates(ivertex)%z = 0.d0
1294) enddo
1295)
1296) #if UGRID_DEBUG
1297) if (ugrid%grid_type == THREE_DIM_GRID) then
1298) call PetscViewerASCIIOpen(option%mycomm,'is_scatter_vert_old_to_new_subsurf.out',viewer, &
1299) ierr);CHKERRQ(ierr)
1300) else
1301) call PetscViewerASCIIOpen(option%mycomm,'is_scatter_vert_old_to_new_surf.out',viewer, &
1302) ierr);CHKERRQ(ierr)
1303) endif
1304) call ISView(is_scatter,viewer,ierr);CHKERRQ(ierr)
1305) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
1306) if (ugrid%grid_type == THREE_DIM_GRID) then
1307) call PetscViewerASCIIOpen(option%mycomm,'is_gather_vert_old_to_new_subsurf.out',viewer, &
1308) ierr);CHKERRQ(ierr)
1309) else
1310) call PetscViewerASCIIOpen(option%mycomm,'is_gather_vert_old_to_new_surf.out',viewer, &
1311) ierr);CHKERRQ(ierr)
1312) endif
1313) call ISView(is_gather,viewer,ierr);CHKERRQ(ierr)
1314) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
1315) #endif
1316)
1317) call VecScatterCreate(vertices_old,is_scatter,vertices_new,is_gather, &
1318) vec_scatter,ierr);CHKERRQ(ierr)
1319) call ISDestroy(is_scatter,ierr);CHKERRQ(ierr)
1320) call ISDestroy(is_gather,ierr);CHKERRQ(ierr)
1321) call VecScatterBegin(vec_scatter,vertices_old,vertices_new, &
1322) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1323) call VecScatterEnd(vec_scatter,vertices_old,vertices_new, &
1324) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1325) call VecScatterDestroy(vec_scatter,ierr);CHKERRQ(ierr)
1326)
1327) #if UGRID_DEBUG
1328) if (ugrid%grid_type == THREE_DIM_GRID) then
1329) call PetscViewerASCIIOpen(option%mycomm,'vertex_coord_old_subsurf.out',viewer, &
1330) ierr);CHKERRQ(ierr)
1331) else
1332) call PetscViewerASCIIOpen(option%mycomm,'vertex_coord_old_surf.out',viewer, &
1333) ierr);CHKERRQ(ierr)
1334) endif
1335) call VecView(vertices_old,viewer,ierr);CHKERRQ(ierr)
1336) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
1337) #endif
1338) call VecDestroy(vertices_old,ierr);CHKERRQ(ierr)
1339)
1340) call VecGetArrayF90(vertices_new,vec_ptr,ierr);CHKERRQ(ierr)
1341) do ivertex = 1, ugrid%num_vertices_local
1342) ugrid%vertices(ivertex)%id = needed_vertices_petsc(ivertex)
1343) ugrid%vertices(ivertex)%x = vec_ptr((ivertex-1)*3+1)
1344) ugrid%vertices(ivertex)%y = vec_ptr((ivertex-1)*3+2)
1345) ugrid%vertices(ivertex)%z = vec_ptr((ivertex-1)*3+3)
1346) pgrid%vertex_coordinates(ivertex)%x = vec_ptr((ivertex-1)*3+1)
1347) pgrid%vertex_coordinates(ivertex)%y = vec_ptr((ivertex-1)*3+2)
1348) pgrid%vertex_coordinates(ivertex)%z = vec_ptr((ivertex-1)*3+3)
1349) enddo
1350) call VecRestoreArrayF90(vertices_new,vec_ptr,ierr);CHKERRQ(ierr)
1351)
1352) #if UGRID_DEBUG
1353) write(string,*) option%myrank
1354) if (ugrid%grid_type == THREE_DIM_GRID) then
1355) string = 'vertex_coord_new' // trim(adjustl(string)) // '_subsurf.out'
1356) else
1357) string = 'vertex_coord_new' // trim(adjustl(string)) // '_surf.out'
1358) endif
1359) call PetscViewerASCIIOpen(PETSC_COMM_SELF,trim(string),viewer, &
1360) ierr);CHKERRQ(ierr)
1361) call VecView(vertices_new,viewer,ierr);CHKERRQ(ierr)
1362) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
1363) #endif
1364)
1365) call VecDestroy(vertices_new,ierr);CHKERRQ(ierr)
1366)
1367) #if UGRID_DEBUG
1368) call printMsg(option,'Setting cell types')
1369) #endif
1370)
1371) allocate(ugrid%cell_type(ugrid%ngmax))
1372)
1373) select case(ugrid%grid_type)
1374) case(THREE_DIM_GRID)
1375) do ghosted_id = 1, ugrid%ngmax
1376) ugrid%cell_type(ghosted_id) = POLY_TYPE
1377) enddo
1378) case default
1379) option%io_buffer = 'Grid type not recognized in UGridPolyhedraDecompose.'
1380) call printErrMsg(option)
1381) end select
1382)
1383) end subroutine UGridPolyhedraDecompose
1384)
1385) ! ************************************************************************** !
1386)
1387) subroutine UGridPolyhedraSetCellCentroids(pgrid,x,y,z, &
1388) x_min,x_max,y_min,y_max,z_min,z_max,option)
1389) !
1390) ! This routine set cell centroids for local+ghosted control volumes.
1391) !
1392) ! Author: Gautam Bisht, LBL
1393) ! Date: 12/28/13
1394) !
1395)
1396) use Option_module
1397) type(unstructured_polyhedra_type), pointer :: pgrid
1398) PetscReal :: x(:), y(:), z(:)
1399) PetscReal :: x_min, x_max, y_min, y_max, z_min, z_max
1400) type(option_type) :: option
1401)
1402) PetscInt :: icell
1403) PetscInt :: ivertex
1404)
1405) do icell = 1,size(pgrid%cell_centroids)
1406) x(icell) = pgrid%cell_centroids(icell)%x
1407) y(icell) = pgrid%cell_centroids(icell)%y
1408) z(icell) = pgrid%cell_centroids(icell)%z
1409) enddo
1410)
1411) do ivertex = 1, pgrid%num_vertices_local
1412) if (x_max < pgrid%vertex_coordinates(ivertex)%x) &
1413) x_max = pgrid%vertex_coordinates(ivertex)%x
1414) if (x_min > pgrid%vertex_coordinates(ivertex)%x) &
1415) x_min = pgrid%vertex_coordinates(ivertex)%x
1416) if (y_max < pgrid%vertex_coordinates(ivertex)%y) &
1417) y_max = pgrid%vertex_coordinates(ivertex)%y
1418) if (y_min > pgrid%vertex_coordinates(ivertex)%y) &
1419) y_min = pgrid%vertex_coordinates(ivertex)%y
1420) if (z_max < pgrid%vertex_coordinates(ivertex)%z) &
1421) z_max = pgrid%vertex_coordinates(ivertex)%z
1422) if (z_min > pgrid%vertex_coordinates(ivertex)%z) &
1423) z_min = pgrid%vertex_coordinates(ivertex)%z
1424) enddo
1425)
1426) end subroutine UGridPolyhedraSetCellCentroids
1427)
1428) ! ************************************************************************** !
1429)
1430) function UGridPolyhedraComputeInternConnect(ugrid, grid_x, &
1431) grid_y, grid_z, option)
1432) !
1433) ! This routine compute internal connectivity of an unstrucutred polyhedra
1434) ! grid.
1435) !
1436) ! Author: Gautam Bisht, LBL
1437) ! Date: 12/28/13
1438) !
1439)
1440) use Connection_module
1441) use Option_module
1442) use Grid_Unstructured_module
1443) use Utility_module, only : DotProduct, CrossProduct
1444)
1445) implicit none
1446)
1447) type(connection_set_type), pointer :: UGridPolyhedraComputeInternConnect
1448) type(grid_unstructured_type) :: ugrid
1449) PetscReal :: grid_x(*), grid_y(*), grid_z(*)
1450) type(option_type) :: option
1451)
1452) type(connection_set_type), pointer :: connections
1453) type(unstructured_polyhedra_type), pointer :: pgrid
1454)
1455) PetscInt :: nconn, iconn
1456) PetscInt :: idual, dual_id
1457)
1458) PetscInt, allocatable :: face_to_vertex(:,:)
1459) PetscInt, allocatable :: cell_to_face(:,:)
1460) PetscInt, allocatable :: face_to_cell(:,:)
1461) PetscInt, allocatable :: vertex_to_cell(:,:)
1462) PetscInt, allocatable :: temp_int(:)
1463) PetscInt, allocatable :: temp_int_2d(:,:)
1464) PetscInt, allocatable :: face_nverts(:)
1465) PetscInt, allocatable :: vertex_ids(:)
1466) PetscInt, allocatable :: dup_face_id(:)
1467) PetscBool, allocatable :: local_boundary_face(:)
1468)
1469) PetscInt :: max_face_per_cell
1470) PetscInt :: max_vert_per_cell
1471) PetscInt :: max_vert_per_face
1472)
1473) PetscInt :: num_match
1474) PetscInt :: found_count
1475) PetscInt :: face_count
1476) PetscInt :: count
1477) PetscInt :: iside
1478) PetscInt :: icell
1479) PetscInt :: dual_local_id
1480) PetscInt :: cell_nvert
1481)
1482) PetscInt :: iface, iface2
1483) PetscInt :: ivertex, ivertex2
1484) PetscInt :: face_id, face_id2
1485) PetscInt :: ghosted_id, ghosted_id2
1486) PetscInt :: local_id, local_id2
1487) PetscInt :: cell_id, cell_id2
1488) PetscInt :: vertex_id, vertex_id2
1489) PetscInt :: cell_type, cell_type2
1490) PetscInt :: face_type, face_type2
1491) PetscInt :: nfaces, nfaces2
1492) PetscInt :: nvertices, nvertices2
1493)
1494) PetscInt :: nintercp
1495) PetscInt :: iintercp
1496) PetscInt :: idx
1497) PetscBool :: face_found
1498) PetscBool :: vertex_found
1499) PetscBool :: cell_found
1500)
1501) PetscReal :: v1(3), v2(3), v3(3)
1502) PetscReal :: n1(3), n2(3), n_up_dn(3)
1503) PetscReal :: dist_up, dist_dn
1504)
1505) type(plane_type) :: plane1
1506) type(point_type) :: point1, point2, point3
1507) type(point_type) :: point_up, point_dn
1508) type(point_type) :: intercept1, intercept2, intercept
1509)
1510) PetscErrorCode :: ierr
1511)
1512) character(len=MAXSTRINGLENGTH) :: string
1513)
1514) pgrid => ugrid%polyhedra_grid
1515)
1516) max_face_per_cell = pgrid%max_nface_per_cell
1517) max_vert_per_face = pgrid%max_nvert_per_face
1518) max_vert_per_cell = ugrid%max_nvert_per_cell
1519)
1520) allocate(vertex_ids(max_vert_per_face))
1521)
1522) ! create mappings of [cells,faces,vertices] to [cells,faces,vertices]
1523) allocate(face_to_vertex(max_vert_per_face, &
1524) max_face_per_cell*ugrid%ngmax))
1525) face_to_vertex = 0
1526)
1527) allocate(cell_to_face(max_face_per_cell, &
1528) ugrid%ngmax))
1529) cell_to_face = 0
1530) allocate(face_to_cell(2,max_face_per_cell* &
1531) ugrid%ngmax))
1532) face_to_cell = 0
1533) allocate(vertex_to_cell(0:ugrid%max_cells_sharing_a_vertex, &
1534) ugrid%num_vertices_local))
1535) vertex_to_cell = 0
1536)
1537) allocate(ugrid%face_to_vertex_natural(max_vert_per_face, &
1538) max_face_per_cell*ugrid%ngmax))
1539) ugrid%face_to_vertex_natural = 0
1540)
1541) allocate(dup_face_id(max_face_per_cell*ugrid%ngmax))
1542) dup_face_id = 0
1543)
1544) face_count = 0
1545) do ghosted_id = 1, ugrid%ngmax
1546) nfaces = pgrid%cell_nfaces(ghosted_id)
1547) do iface = 1, nfaces
1548) face_count = face_count + 1
1549) cell_to_face(iface,ghosted_id) = face_count
1550) face_to_cell(1,face_count) = ghosted_id
1551) nvertices = pgrid%face_nverts(face_count)
1552) vertex_ids = pgrid%face_vertids(:,face_count)
1553) do ivertex = 1, nvertices
1554) face_to_vertex(ivertex,face_count) = vertex_ids(ivertex)
1555) if (face_to_vertex(ivertex,face_count) > 0) then
1556) ugrid%face_to_vertex_natural(ivertex,face_count) = &
1557) ugrid%vertex_ids_natural(face_to_vertex(ivertex,face_count))
1558) endif
1559) enddo
1560) enddo
1561) enddo
1562)
1563) !
1564) ! Remove duplicate faces:
1565) !
1566) ! A cell (cell_id) and Neighboring-Cell (cell_id2) will only share ONE face.
1567) ! Find the face that cell_id ane cell_id2 share and remove it.
1568) !
1569) ! Method:
1570) ! - Pick i-th face (iface) of cell_id and check if ALL the vertices of
1571) ! the iface are present in cell_id2. If all the vertices of iface are
1572) ! not present in cell_id2, move to the next face.
1573) ! - After finding the iface, now find iface2 in cell_id2 that
1574) ! corresponds to iface.
1575) ! - Check to ensure that atleast on face of cell_id is shared
1576) ! with cell_id2.
1577) iconn = 0
1578) do local_id = 1, ugrid%nlmax
1579) ! Select a cell and find number of vertices
1580) cell_id = local_id
1581) nfaces = pgrid%cell_nfaces(local_id)
1582) do idual = 1, ugrid%cell_neighbors_local_ghosted(0,local_id)
1583) cell_id2 = &
1584) abs(ugrid%cell_neighbors_local_ghosted(idual,local_id))
1585) ! If cell-id is neighbor is lower, skip it
1586) if (cell_id2 <= cell_id) cycle
1587) nfaces2 = pgrid%cell_nfaces(cell_id2)
1588) face_found = PETSC_FALSE
1589)
1590) do iface = 1, nfaces
1591) face_id = cell_to_face(iface,cell_id)
1592) nvertices = pgrid%face_nverts(face_id)
1593)
1594) do ivertex = 1, nvertices
1595) vertex_id = face_to_vertex(ivertex,face_id)
1596) vertex_found = PETSC_FALSE
1597) do ivertex2 = 1, ugrid%cell_vertices(0,cell_id2)
1598) vertex_id2 = ugrid%cell_vertices(ivertex2,cell_id2)
1599) if (vertex_id == vertex_id2) then
1600) vertex_found = PETSC_TRUE
1601) exit
1602) endif
1603) enddo
1604) !
1605) ! If ivertex of iface of the Cell is not present as vertices of the
1606) ! Neighboring-Cell, then iface is not the shared face. Skip iterating
1607) ! over the remaing vertices of iface
1608) if (.not.vertex_found) exit
1609) enddo ! do-loop 'ivertex'
1610)
1611) if (vertex_found) then
1612) ! All the vertices of iface are present in the Neighboring cells.
1613) ! Thus, iface is the shared face.
1614) face_found = PETSC_TRUE
1615)
1616) ! Now, we have to find iface2 that corresponds to iface
1617) do iface2 = 1, nfaces2
1618) face_id2 = cell_to_face(iface2,cell_id2)
1619) nvertices2 = pgrid%face_nverts(face_id2)
1620) ! iface and iface2 need to have same number of vertices
1621) if (nvertices == nvertices2) then
1622) ! Count the number of vertices of iface which match vertices of iface2
1623) num_match = 0
1624) do ivertex = 1, nvertices
1625) vertex_id = face_to_vertex(ivertex,face_id)
1626) vertex_found = PETSC_FALSE
1627) do ivertex2 = 1, nvertices2
1628) vertex_id2 = face_to_vertex(ivertex2,face_id2)
1629) if (vertex_id == vertex_id2) then
1630) vertex_found = PETSC_TRUE
1631) num_match = num_match + 1
1632) exit
1633) endif
1634) enddo
1635) !
1636) ! If vertex_id of face_id not found as one of the vertices of face_id2,
1637) ! face_id2 is not shared between cells
1638) if (.not.vertex_found) exit
1639) enddo
1640)
1641) if (num_match == nvertices) then
1642) ! remove duplicate face
1643) if (face_id2 > face_id) then
1644) #ifdef UGRID_DEBUG
1645) write(string,*) option%myrank, face_id2, ' -> ', face_id
1646) option%io_buffer = 'Duplicated face removed:' // trim(string)
1647) call printMsg(option)
1648) #endif
1649) cell_to_face(iface2,cell_id2) = face_id
1650) ! flag face_id2 as removed
1651) face_to_cell(1,face_id2) = -face_to_cell(1,face_id2)
1652) face_to_cell(2,face_id2) = cell_id
1653) ! add cell_id2 to face_ids list
1654) face_to_cell(2,face_id) = cell_id2
1655) dup_face_id(face_id) = face_id2
1656) else
1657) #ifdef UGRID_DEBUG
1658) write(string,*) option%myrank, face_id, ' -> ', face_id2
1659) option%io_buffer = 'Duplicated face removed:' // trim(string)
1660) call printMsg(option)
1661) #endif
1662) cell_to_face(iface,cell_id) = face_id2
1663) ! flag face_id as removed
1664) face_to_cell(1,face_id) = -face_to_cell(1,face_id)
1665) face_to_cell(2,face_id) = cell_id2
1666) ! add cell_id to face_ids2 list
1667) face_to_cell(2,face_id2) = cell_id
1668) dup_face_id(face_id2)= face_id
1669) endif ! if (face_id2 > face_id)
1670) exit
1671) endif ! if (num_match == nvertices)
1672) endif ! if (nvertices == nvertices2)
1673)
1674) ! Check that one shared face was found between the Cell and Neighboring-Cell
1675) if (.not.face_found) then
1676) write(string,*) option%myrank
1677) string = '(' // trim(adjustl(string)) // ')'
1678) write(*,'(a,'' local_id = '',i3,'' natural_id = '',i3,'' vertices: '',8i3)') &
1679) trim(string), &
1680) cell_id,ugrid%cell_ids_natural(cell_id), &
1681) (ugrid%vertex_ids_natural( &
1682) ugrid%cell_vertices(ivertex,cell_id)), &
1683) ivertex=1,ugrid%cell_vertices(0,cell_id))
1684) write(*,'(a,'' local_id2 = '',i3,'' natural_id2 = '',i3,'' vertices2: '',8i3)') &
1685) trim(string), &
1686) cell_id2,ugrid%cell_ids_natural(cell_id2), &
1687) (ugrid%vertex_ids_natural( &
1688) ugrid%cell_vertices(ivertex2,cell_id2)), &
1689) ivertex2=1,ugrid%cell_vertices(0,cell_id2))
1690) option%io_buffer='No shared face found.'
1691) call printErrMsgByRank(option)
1692) endif
1693)
1694) enddo ! do-loop iface2
1695)
1696) exit
1697) endif ! if (vertex_found)
1698) enddo ! do-loop 'iface'
1699)
1700) ! Check that one shared face was found between the Cell and Neighboring-Cell
1701) if (.not.face_found) then
1702) write(string,*) option%myrank
1703) string = '(' // trim(adjustl(string)) // ')'
1704) write(*,'(a,'' local_id = '',i3,'' natural_id = '',i3,'' vertices: '',8i3)') &
1705) trim(string), &
1706) cell_id,ugrid%cell_ids_natural(cell_id), &
1707) (ugrid%vertex_ids_natural( &
1708) ugrid%cell_vertices(ivertex,cell_id)), &
1709) ivertex=1,ugrid%cell_vertices(0,cell_id))
1710) write(*,'(a,'' local_id2 = '',i3,'' natural_id2 = '',i3,'' vertices2: '',8i3)') &
1711) trim(string), &
1712) cell_id2,ugrid%cell_ids_natural(cell_id2), &
1713) (ugrid%vertex_ids_natural( &
1714) ugrid%cell_vertices(ivertex2,cell_id2)), &
1715) ivertex2=1,ugrid%cell_vertices(0,cell_id2))
1716) option%io_buffer='No shared face found.'
1717) call printErrMsgByRank(option)
1718) endif
1719)
1720) enddo ! do-loop 'idual'
1721) enddo ! do-loop 'local_id'
1722)
1723) ! GB: Add a check for dup_face_id
1724)
1725) ! count up the # of faces
1726) face_count = 0
1727) do iface = 1, size(face_to_cell,2)
1728) if (face_to_cell(1,iface) > 0) &
1729) face_count = face_count + 1
1730) enddo
1731) allocate(ugrid%face_to_vertex(max_vert_per_face,face_count))
1732) face_count = 0
1733) do iface = 1, size(face_to_cell,2)
1734) if (face_to_cell(1,iface) > 0) then
1735) face_count = face_count + 1
1736) ugrid%face_to_vertex(:,face_count) = face_to_vertex(:,iface)
1737) endif
1738) enddo
1739) deallocate(face_to_vertex)
1740) ! reallocate face_to_cell to proper size
1741) allocate(temp_int_2d(2,face_count))
1742) allocate(face_nverts(face_count))
1743) allocate(temp_int(size(face_to_cell,2)))
1744) temp_int = 0
1745) face_count = 0
1746) do iface = 1, size(face_to_cell,2)
1747) if (face_to_cell(1,iface) > 0) then
1748) face_count = face_count + 1
1749) temp_int_2d(:,face_count) = face_to_cell(:,iface)
1750) face_nverts(face_count) = pgrid%face_nverts(iface)
1751) temp_int(iface) = face_count
1752) endif
1753) enddo
1754) deallocate(face_to_cell)
1755) allocate(face_to_cell(2,face_count))
1756) face_to_cell = temp_int_2d
1757) deallocate(temp_int_2d)
1758)
1759) ! remap faces in cells using temp_int from above
1760) do iface = 1, size(face_to_cell,2)
1761) face_id = iface
1762) do icell = 1,2
1763) cell_id = face_to_cell(icell,face_id)
1764) ! check for exterior face
1765) if (cell_id < 1) cycle
1766) cell_found = PETSC_FALSE
1767) nfaces = pgrid%cell_nfaces(cell_id)
1768) do iface2 = 1, nfaces
1769) face_id2 = cell_to_face(iface2,cell_id)
1770) if (face_id < 0) cycle
1771) if (face_id == temp_int(face_id2)) then
1772) cell_found = PETSC_TRUE
1773) cell_to_face(iface2,cell_id) = face_id
1774) exit
1775) endif
1776) enddo
1777)
1778) if (.not.cell_found) then
1779) option%io_buffer = 'POLYHEDRA_UGRID: Remapping of cell face id unsuccessful'
1780) call printErrMsg(option)
1781) endif
1782) enddo
1783) enddo
1784) deallocate(temp_int)
1785)
1786) do ghosted_id = 1, ugrid%ngmax
1787) do ivertex = 1, ugrid%cell_vertices(0,ghosted_id)
1788) vertex_id = ugrid%cell_vertices(ivertex,ghosted_id)
1789) if ( vertex_id <= 0) cycle
1790) count = vertex_to_cell(0,vertex_id) + 1
1791) if (count > ugrid%max_cells_sharing_a_vertex) then
1792) write(string,*) 'Vertex can be shared by at most by ', &
1793) ugrid%max_cells_sharing_a_vertex, &
1794) ' cells. Rank = ', option%myrank, ' vertex_id = ', vertex_id, ' exceeds it.'
1795) option%io_buffer = string
1796) call printErrMsg(option)
1797) endif
1798) vertex_to_cell(count,vertex_id) = ghosted_id
1799) vertex_to_cell(0,vertex_id) = count
1800) enddo
1801) enddo
1802)
1803) nconn = 0
1804) do local_id = 1, ugrid%nlmax
1805) do idual = 1, ugrid%cell_neighbors_local_ghosted(0,local_id)
1806) dual_id = ugrid%cell_neighbors_local_ghosted(idual,local_id)
1807) ! count all ghosted connections (dual_id < 0)
1808) ! only count connection with cells of larger ids to avoid double counts
1809) if (dual_id < 0 .or. local_id < dual_id) then
1810) nconn = nconn + 1
1811) endif
1812) enddo
1813) enddo
1814)
1815) connections => ConnectionCreate(nconn,INTERNAL_CONNECTION_TYPE)
1816)
1817) allocate(ugrid%connection_to_face(nconn))
1818) ugrid%connection_to_face = 0
1819)
1820) ! loop over connection again
1821) iconn = 0
1822) do local_id = 1, ugrid%nlmax
1823) do idual = 1, ugrid%cell_neighbors_local_ghosted(0,local_id)
1824) dual_local_id = &
1825) ugrid%cell_neighbors_local_ghosted(idual,local_id)
1826) if (local_id < abs(dual_local_id)) then
1827) iconn = iconn + 1
1828) ! find face
1829) face_found = PETSC_FALSE
1830) do iface = 1, ugrid%cell_vertices(0,local_id)
1831) face_id = cell_to_face(iface,local_id)
1832) do iside = 1,2
1833) cell_id2 = face_to_cell(iside,face_id)
1834) if (cell_id2 == abs(dual_local_id)) then
1835) face_found = PETSC_TRUE
1836) exit
1837) endif
1838) enddo
1839) if (face_found) exit
1840) enddo
1841) if (face_found) then
1842) ugrid%connection_to_face(iconn) = face_id
1843) else
1844) write(string,*) option%myrank,local_id,dual_local_id
1845) option%io_buffer = 'face not found in connection loop' // trim(string)
1846) call printErrMsg(option)
1847) endif
1848)
1849) face_found = PETSC_FALSE
1850) do iface2 = 1, ugrid%cell_vertices(0,cell_id2)
1851) if (cell_to_face(iface,local_id) == &
1852) cell_to_face(iface2,cell_id2)) then
1853) face_found = PETSC_TRUE
1854) exit
1855) endif
1856) enddo
1857)
1858) if (.not.face_found) then
1859) write(string,*) option%myrank, iface, cell_id2
1860) option%io_buffer = 'global face not found' // trim(string)
1861) call printErrMsg(option)
1862) endif
1863)
1864) connections%id_up(iconn) = local_id
1865) connections%id_dn(iconn) = abs(dual_local_id)
1866) connections%face_id(iconn) = cell_to_face(iface,local_id)
1867)
1868) point_up%x = grid_x(local_id)
1869) point_up%y = grid_y(local_id)
1870) point_up%z = grid_z(local_id)
1871) point_dn%x = grid_x(abs(dual_local_id))
1872) point_dn%y = grid_y(abs(dual_local_id))
1873) point_dn%z = grid_z(abs(dual_local_id))
1874)
1875) ! Find intercept
1876) nintercp = face_nverts(cell_to_face(iface,local_id)) - 2
1877)
1878) intercept%x = 0.d0
1879) intercept%y = 0.d0
1880) intercept%z = 0.d0
1881)
1882) do iintercp = 0, nintercp - 1
1883) idx = ugrid%face_to_vertex(1 + iintercp,face_id)
1884) point1 = ugrid%vertices(idx)
1885) idx = ugrid%face_to_vertex(2 + iintercp,face_id)
1886) point2 = ugrid%vertices(idx)
1887) idx = ugrid%face_to_vertex(3 + iintercp,face_id)
1888) point3 = ugrid%vertices(idx)
1889)
1890) call UCellComputePlane(plane1,point1,point2,point3)
1891) call UCellGetPlaneIntercept(plane1,point_up,point_dn,intercept1)
1892)
1893) intercept%x = intercept%x + intercept1%x
1894) intercept%y = intercept%y + intercept1%y
1895) intercept%z = intercept%z + intercept1%z
1896)
1897) enddo
1898)
1899) intercept%x = intercept%x/nintercp
1900) intercept%y = intercept%y/nintercp
1901) intercept%z = intercept%z/nintercp
1902)
1903) ! This is very crude, but for now use average location of intercept
1904) v1(1) = intercept%x-point_up%x
1905) v1(2) = intercept%y-point_up%y
1906) v1(3) = intercept%z-point_up%z
1907) v2(1) = point_dn%x-intercept%x
1908) v2(2) = point_dn%y-intercept%y
1909) v2(3) = point_dn%z-intercept%z
1910) dist_up = sqrt(DotProduct(v1,v1))
1911) dist_dn = sqrt(DotProduct(v2,v2))
1912)
1913) connections%dist(-1:3,iconn) = 0.d0
1914) connections%dist(-1,iconn) = dist_up/(dist_up + dist_dn)
1915) connections%dist(0,iconn) = dist_up + dist_dn
1916) v3 = v1 + v2
1917) connections%dist(1:3,iconn) = v3/sqrt(DotProduct(v3,v3))
1918) connections%area(iconn) = pgrid%face_areas(connections%face_id(iconn))
1919) connections%intercp(1,iconn) = intercept%x
1920) connections%intercp(2,iconn) = intercept%y
1921) connections%intercp(3,iconn) = intercept%z
1922)
1923) endif ! (local_id < abs(dual_local_id))
1924)
1925) enddo
1926) enddo ! 'local_id'
1927)
1928) allocate(ugrid%face_area(face_count))
1929) allocate(ugrid%face_centroid(face_count))
1930)
1931) do local_id = 1, ugrid%nlmax
1932) do iface = 1, pgrid%cell_nfaces(local_id)
1933) face_id = cell_to_face(iface,local_id)
1934)
1935) ugrid%face_centroid(face_id)%x = &
1936) pgrid%face_centroids(face_id)%x
1937) ugrid%face_centroid(face_id)%y = &
1938) pgrid%face_centroids(face_id)%y
1939) ugrid%face_centroid(face_id)%z = &
1940) pgrid%face_centroids(face_id)%z
1941)
1942) ugrid%face_area(face_id) = &
1943) pgrid%face_areas(face_id)
1944)
1945) enddo
1946) enddo
1947)
1948) allocate(ugrid%face_to_cell_ghosted(size(face_to_cell,1), &
1949) size(face_to_cell,2)))
1950) ugrid%face_to_cell_ghosted = face_to_cell
1951) allocate(ugrid%cell_to_face_ghosted(size(cell_to_face,1), &
1952) size(cell_to_face,2)))
1953) ugrid%cell_to_face_ghosted(:,:) = cell_to_face(:,:)
1954)
1955)
1956) #if UGRID_DEBUG
1957) write(string,*) option%myrank
1958) string = 'face_to_cell' // trim(adjustl(string)) // '.out'
1959) open(unit=86,file=trim(string))
1960) do iface = 1, size(face_to_cell,2)
1961) write(86,'(i5)',advance='no') iface
1962) write(86,'(i5)',advance='no') face_to_cell(1,iface)
1963) write(86,'(i5)',advance='no') face_to_cell(2,iface)
1964) write(86,'(a)') ""
1965) enddo
1966) close(86)
1967)
1968) write(string,*) option%myrank
1969) string = 'cell_to_face' // trim(adjustl(string)) // '.out'
1970) open(unit=86,file=trim(string))
1971) do ghosted_id = 1, ugrid%ngmax
1972) write(86,'(i5)',advance='no') ghosted_id
1973) do iface = 1,max_face_per_cell
1974) write(86,'(i5)',advance='no') cell_to_face(iface,ghosted_id)
1975) enddo
1976) write(86,'(a)') ""
1977) enddo
1978) close(86)
1979)
1980) write(string,*) option%myrank
1981) string = 'poly_cell_vertids' // trim(adjustl(string)) // '.out'
1982) open(unit=86,file=trim(string))
1983) do ghosted_id = 1, ugrid%ngmax
1984) write(86,'(i5)',advance='no') ghosted_id
1985) do ivertex = 1,ugrid%max_nvert_per_cell
1986) write(86,'(i5)',advance='no') pgrid%cell_vertids(ivertex,ghosted_id)
1987) enddo
1988) write(86,'(a)') ""
1989) enddo
1990) close(86)
1991)
1992) write(string,*) option%myrank
1993) string = 'poly_cell_faceids' // trim(adjustl(string)) // '.out'
1994) open(unit=86,file=trim(string))
1995) do ghosted_id = 1, ugrid%ngmax
1996) write(86,'(i5)',advance='no') ghosted_id
1997) do ivertex = 1,pgrid%max_nface_per_cell
1998) write(86,'(i5)',advance='no') pgrid%cell_faceids(ivertex,ghosted_id)
1999) enddo
2000) write(86,'(a)') ""
2001) enddo
2002) close(86)
2003)
2004) write(string,*) option%myrank
2005) string = 'poly_face_vertids' // trim(adjustl(string)) // '.out'
2006) open(unit=86,file=trim(string))
2007) do iface = 1, pgrid%num_faces_local
2008) write(86,'(i5)',advance='no') iface
2009) do ivertex = 1,pgrid%max_nvert_per_face
2010) write(86,'(i5)',advance='no') pgrid%face_vertids(ivertex,iface)
2011) enddo
2012) write(86,'(a)') ""
2013) enddo
2014) close(86)
2015)
2016) #endif
2017)
2018) deallocate(cell_to_face)
2019) deallocate(face_to_cell)
2020) deallocate(vertex_to_cell)
2021) deallocate(vertex_ids)
2022) deallocate(dup_face_id)
2023)
2024) UGridPolyhedraComputeInternConnect => connections
2025)
2026) end function UGridPolyhedraComputeInternConnect
2027)
2028) ! ************************************************************************** !
2029)
2030) subroutine UGridPolyhedraComputeVolumes(ugrid, option, volume)
2031) !
2032) ! This routine sets volumes of local control volumes.
2033) !
2034) ! Author: Gautam Bisht, LBL
2035) ! Date: 12/28/13
2036) !
2037)
2038) use Option_module
2039)
2040) implicit none
2041)
2042) type(grid_unstructured_type) :: ugrid
2043) type(option_type) :: option
2044) Vec :: volume
2045)
2046) type(unstructured_polyhedra_type), pointer :: pgrid
2047)
2048) PetscInt :: icell
2049) PetscReal, pointer :: vec_ptr(:)
2050) PetscErrorCode :: ierr
2051)
2052) pgrid => ugrid%polyhedra_grid
2053)
2054) call VecGetArrayF90(volume,vec_ptr,ierr);CHKERRQ(ierr)
2055) do icell = 1, ugrid%nlmax
2056) vec_ptr(icell) = pgrid%cell_volumes(icell)
2057) enddo
2058) call VecRestoreArrayF90(volume,vec_ptr,ierr);CHKERRQ(ierr)
2059)
2060) end subroutine UGridPolyhedraComputeVolumes
2061)
2062) ! ************************************************************************** !
2063)
2064) subroutine UGridPolyhedraPopulateConnection(ugrid, connection, iface_cell, &
2065) iconn, ghosted_id, option)
2066) !
2067) ! This routine computes details about boundary connections (area, dist, etc)
2068) !
2069) ! Author: Gautam Bisht, LBL
2070) ! Date: 12/28/13
2071) !
2072)
2073) use Connection_module
2074) use Utility_module, only : DotProduct
2075) use Option_module
2076) use Grid_Unstructured_Cell_module
2077)
2078) implicit none
2079)
2080) type(grid_unstructured_type) :: ugrid
2081) type(connection_set_type) :: connection
2082) PetscInt :: iface_cell
2083) PetscInt :: iconn
2084) PetscInt :: ghosted_id
2085) type(option_type) :: option
2086)
2087) PetscInt :: face_id
2088) PetscInt :: ivert,vert_id
2089) PetscInt :: face_type
2090) PetscReal :: v1(3),v2(3),n_dist(3), dist
2091) type(point_type) :: vertex_8(8)
2092) type(plane_type) :: plane
2093) type(point_type) :: point, vertex1, vertex2, vertex3, intercept
2094) type(unstructured_polyhedra_type), pointer :: pgrid
2095) character(len=MAXWORDLENGTH) :: word
2096) PetscErrorCode :: ierr
2097)
2098) pgrid => ugrid%polyhedra_grid
2099)
2100) select case(connection%itype)
2101) case(BOUNDARY_CONNECTION_TYPE)
2102) if (iface_cell == 0) then
2103) write(word,*) ghosted_id
2104) option%io_buffer = 'Face id undefined for cell ' // &
2105) trim(adjustl(word)) // &
2106) ' in boundary condition. Should this be a source/sink?'
2107) call printErrMsgByRank(option)
2108) endif
2109) ! Compute cell centeroid
2110) v2(1) = pgrid%cell_centroids(ghosted_id)%x
2111) v2(2) = pgrid%cell_centroids(ghosted_id)%y
2112) v2(3) = pgrid%cell_centroids(ghosted_id)%z
2113)
2114) ! Instead of connecting centroid with face center, calculate the shortest
2115) ! distance between the centroid and face and use that distance - geh
2116)
2117) point%x = v2(1)
2118) point%y = v2(2)
2119) point%z = v2(3)
2120)
2121) !face_id = ugrid%cell_to_face_ghosted(iface_cell, ghosted_id)
2122) face_id = pgrid%cell_faceids(iface_cell, ghosted_id)
2123) intercept%x = pgrid%face_centroids(face_id)%x
2124) intercept%y = pgrid%face_centroids(face_id)%y
2125) intercept%z = pgrid%face_centroids(face_id)%z
2126)
2127) ! Compute distance vector: cell_center - face_centroid
2128) v1(1) = v2(1) - intercept%x
2129) v1(2) = v2(2) - intercept%y
2130) v1(3) = v2(3) - intercept%z
2131)
2132) dist = sqrt(DotProduct(v1, v1))
2133) n_dist = v1/dist
2134) connection%dist(0, iconn) = dist
2135) connection%dist(1, iconn) = n_dist(1)
2136) connection%dist(2, iconn) = n_dist(2)
2137) connection%dist(3, iconn) = n_dist(3)
2138) connection%area(iconn) = pgrid%face_areas(face_id)
2139) connection%intercp(1,iconn)= intercept%x
2140) connection%intercp(2,iconn)= intercept%y
2141) connection%intercp(3,iconn)= intercept%z
2142) connection%face_id(iconn) = face_id
2143)
2144) end select
2145)
2146) end subroutine UGridPolyhedraPopulateConnection
2147)
2148) ! ************************************************************************** !
2149)
2150) subroutine UGridPolyhedraGetCellsInRectangle(x_min, x_max, y_min, y_max, z_min, z_max, &
2151) ugrid, option, num_cells, &
2152) cell_ids, cell_face_ids)
2153) !
2154) ! This routine returns cells that are within a cube.
2155) !
2156) ! Author: Gautam Bisht, LBL
2157) ! Date: 12/28/13
2158) !
2159) use Option_module
2160) use Utility_module, only : reallocateIntArray
2161)
2162) implicit none
2163)
2164) PetscReal :: x_min, x_max, y_min, y_max, z_min, z_max
2165) type(grid_unstructured_type) :: ugrid
2166) type(option_type) :: option
2167) PetscInt :: num_cells
2168) PetscInt, pointer :: cell_ids(:)
2169) PetscInt, pointer :: cell_face_ids(:)
2170)
2171) type(unstructured_polyhedra_type), pointer :: pgrid
2172) PetscInt :: cell_type, num_faces, iface, face_type
2173) PetscInt :: vertex_id
2174) PetscInt :: num_vertices, ivertex
2175) PetscInt :: local_id, ghosted_id
2176) type(point_type) :: point
2177)
2178) PetscReal :: x_min_adj, x_max_adj, y_min_adj, y_max_adj, z_min_adj, z_max_adj
2179) PetscReal :: pert
2180) PetscBool :: in_rectangle
2181)
2182) PetscInt, pointer :: temp_cell_array(:), temp_face_array(:)
2183) PetscInt :: temp_array_size
2184) PetscInt :: face_id
2185)
2186) pgrid => ugrid%polyhedra_grid
2187)
2188) temp_array_size = 100
2189) allocate(temp_cell_array(temp_array_size))
2190) allocate(temp_face_array(temp_array_size))
2191) temp_cell_array = 0
2192) temp_face_array = 0
2193)
2194) ! enlarge box slightly
2195) pert = max(1.d-8*(x_max-x_min),1.d-8)
2196) x_min_adj = x_min - pert
2197) x_max_adj = x_max + pert
2198) pert = max(1.d-8*(y_max-y_min),1.d-8)
2199) y_min_adj = y_min - pert
2200) y_max_adj = y_max + pert
2201) pert = max(1.d-8*(z_max-z_min),1.d-8)
2202) z_min_adj = z_min - pert
2203) z_max_adj = z_max + pert
2204)
2205) do local_id = 1, ugrid%nlmax
2206) ghosted_id = local_id ! ghosted ids are same for first nlocal cells
2207) num_faces = pgrid%cell_nfaces(ghosted_id)
2208) do iface = 1, num_faces
2209) face_id = ugrid%cell_to_face_ghosted(iface, ghosted_id)
2210) num_vertices = pgrid%face_nverts(ghosted_id)
2211) in_rectangle = PETSC_TRUE
2212) do ivertex = 1, num_vertices
2213) !vertex_id = pgrid%face_vertids(ivertex,face_id)
2214) vertex_id = ugrid%face_to_vertex(ivertex,face_id)
2215) point = ugrid%vertices(vertex_id)
2216) if (point%x < x_min_adj .or. &
2217) point%x > x_max_adj .or. &
2218) point%y < y_min_adj .or. &
2219) point%y > y_max_adj .or. &
2220) point%z < z_min_adj .or. &
2221) point%z > z_max_adj) then
2222) in_rectangle = PETSC_FALSE
2223) exit
2224) endif
2225) enddo
2226)
2227) if (in_rectangle) then
2228) num_cells = num_cells + 1
2229) if (num_cells > temp_array_size) then
2230) call reallocateIntArray(temp_cell_array,temp_array_size)
2231) temp_array_size = temp_array_size / 2 ! convert back for next call
2232) call reallocateIntArray(temp_face_array,temp_array_size)
2233) endif
2234) temp_cell_array(num_cells) = local_id
2235) temp_face_array(num_cells) = iface
2236) endif
2237)
2238) enddo
2239) enddo
2240)
2241) allocate(cell_ids(num_cells))
2242) allocate(cell_face_ids(num_cells))
2243) cell_ids = temp_cell_array(1:num_cells)
2244) cell_face_ids = temp_face_array(1:num_cells)
2245) deallocate(temp_cell_array)
2246) nullify(temp_cell_array)
2247) deallocate(temp_face_array)
2248) nullify(temp_face_array)
2249)
2250) end subroutine UGridPolyhedraGetCellsInRectangle
2251)
2252) ! ************************************************************************** !
2253)
2254) subroutine UGridPolyhedraComputeOutputInfo(ugrid, nL2G, nG2L, nG2A, option)
2255) !
2256) ! This routine computes informations later required to write tecplot output
2257) !
2258) ! Author: Gautam Bisht, LBL
2259) ! Date: 12/29/13
2260) !
2261)
2262) use Option_module
2263)
2264) implicit none
2265)
2266) #include "petsc/finclude/petscvec.h"
2267) #include "petsc/finclude/petscvec.h90"
2268) #include "petsc/finclude/petscis.h"
2269) #include "petsc/finclude/petscis.h90"
2270) #include "petsc/finclude/petscviewer.h"
2271)
2272) type(grid_unstructured_type) :: ugrid
2273) type(option_type) :: option
2274) PetscInt, pointer :: nL2G(:)
2275) PetscInt, pointer :: nG2L(:)
2276) PetscInt, pointer :: nG2A(:)
2277)
2278) type(unstructured_polyhedra_type), pointer :: pgrid
2279)
2280) Vec :: nat_cv_proc_rank
2281) Vec :: ghosted_cv_proc_rank
2282) VecScatter :: vec_scat
2283) IS :: is_scatter
2284) IS :: is_gather
2285)
2286) Vec :: ghosted_vec
2287) Vec :: natural_vec
2288)
2289) PetscInt :: istart
2290) PetscInt :: iend
2291) PetscInt :: iface
2292) PetscInt :: ivertex
2293) PetscInt :: count
2294)
2295) PetscInt :: local_id
2296) PetscInt :: ghosted_id
2297) PetscInt :: dual_id
2298) PetscInt :: face_id
2299) PetscInt :: vertex_id
2300)
2301) PetscInt :: pgridface_id
2302)
2303) PetscInt, allocatable :: int_array(:)
2304) PetscReal, allocatable :: real_array(:)
2305) PetscScalar,pointer :: v_loc_p(:)
2306)
2307) PetscErrorCode :: ierr
2308)
2309) pgrid => ugrid%polyhedra_grid
2310)
2311) pgrid%ugrid_num_faces_local = size(ugrid%face_to_cell_ghosted,2)
2312) allocate(pgrid%ugridf2pgridf(pgrid%ugrid_num_faces_local))
2313) pgrid%ugridf2pgridf = 0
2314)
2315) ! Initialize mapping of faces between unstructured_grid and polyhedra_grid
2316) do ghosted_id = 1, ugrid%ngmax
2317) do iface = 1, pgrid%cell_nfaces(ghosted_id)
2318) face_id = ugrid%cell_to_face_ghosted(iface, ghosted_id)
2319)
2320) if (pgrid%ugridf2pgridf(face_id) == 0) then ! mapping not initialized
2321) pgrid%ugridf2pgridf(face_id) = pgrid%cell_faceids(iface, ghosted_id)
2322) else
2323) ! iface-th of ghosted_id-cell is an internal face shared by:
2324) ! - two local control volumes, or
2325) ! - a local and ghosted control volume
2326)
2327) ! Determine if ghosted_id is the upwind control volume.
2328) if (ugrid%face_to_cell_ghosted(1, face_id) == ghosted_id) then
2329) ! Map ugrid face to corresponding face of upwind control volume
2330) ! in pgrid. Upwind control volume is choosen because the order of
2331) ! vertices forming a face is required for output. Also, unit normal
2332) ! vector point from upwind-to-downwind control volume.
2333) pgrid%ugridf2pgridf(face_id) = pgrid%cell_faceids(iface, ghosted_id)
2334) endif
2335) endif
2336) enddo
2337) enddo
2338)
2339) ! Find number of global unique faces. This is required for output
2340) call VecCreateMPI(option%mycomm, ugrid%nlmax, PETSC_DETERMINE, nat_cv_proc_rank, &
2341) ierr);CHKERRQ(ierr)
2342) call VecCreateMPI(option%mycomm, ugrid%ngmax, PETSC_DETERMINE, ghosted_cv_proc_rank, &
2343) ierr);CHKERRQ(ierr)
2344)
2345) ! Populate a vector that contains rank of procoessor on which a given
2346) ! control volume is active.
2347) allocate(int_array(ugrid%nlmax))
2348) allocate(real_array(ugrid%nlmax))
2349) do local_id = 1,ugrid%nlmax
2350) int_array(local_id) = nG2A(nL2G(local_id))
2351) real_array(local_id) = option%myrank
2352) enddo
2353) int_array = int_array - 1
2354)
2355) call VecSetValues(nat_cv_proc_rank, ugrid%nlmax, int_array, real_array, INSERT_VALUES, &
2356) ierr);CHKERRQ(ierr)
2357) call VecAssemblyBegin(nat_cv_proc_rank, ierr);CHKERRQ(ierr)
2358) call VecAssemblyEnd(nat_cv_proc_rank, ierr);CHKERRQ(ierr)
2359) deallocate(int_array)
2360) deallocate(real_array)
2361)
2362) ! Find processor rank for ghost control volumes by scattering data stored in
2363) ! vector nat_cv_proc_rank.
2364) allocate(int_array(ugrid%ngmax))
2365)
2366) do ghosted_id = 1, ugrid%ngmax
2367) int_array(ghosted_id) = nG2A(ghosted_id)
2368) enddo
2369) int_array = int_array - 1
2370) call ISCreateBlock(option%mycomm, 1, ugrid%ngmax, int_array, PETSC_COPY_VALUES, &
2371) is_scatter, ierr);CHKERRQ(ierr)
2372)
2373) call VecGetOwnershipRange(ghosted_cv_proc_rank, istart, iend, &
2374) ierr);CHKERRQ(ierr)
2375) do ghosted_id = 1, ugrid%ngmax
2376) int_array(ghosted_id) = ghosted_id + istart
2377) enddo
2378) int_array = int_array - 1
2379) call ISCreateBlock(option%mycomm, 1, ugrid%ngmax, int_array, PETSC_COPY_VALUES, &
2380) is_gather, ierr);CHKERRQ(ierr)
2381)
2382) call VecScatterCreate(nat_cv_proc_rank, is_scatter, ghosted_cv_proc_rank, is_gather, &
2383) vec_scat, ierr);CHKERRQ(ierr)
2384) call ISDestroy(is_scatter, ierr);CHKERRQ(ierr)
2385) call ISDestroy(is_gather, ierr);CHKERRQ(ierr)
2386) deallocate(int_array)
2387)
2388) call VecScatterBegin(vec_scat, nat_cv_proc_rank, ghosted_cv_proc_rank, &
2389) INSERT_VALUES, SCATTER_FORWARD, ierr);CHKERRQ(ierr)
2390) call VecScatterEnd(vec_scat, nat_cv_proc_rank, ghosted_cv_proc_rank, &
2391) INSERT_VALUES, SCATTER_FORWARD, ierr);CHKERRQ(ierr)
2392) call VecScatterDestroy(vec_scat, ierr);CHKERRQ(ierr)
2393)
2394) ! Find the number of unique faces
2395) allocate(pgrid%uface_localids(pgrid%ugrid_num_faces_local))
2396) allocate(pgrid%uface_nverts(pgrid%ugrid_num_faces_local))
2397) allocate(pgrid%uface_left_natcellids(pgrid%ugrid_num_faces_local))
2398) allocate(pgrid%uface_right_natcellids(pgrid%ugrid_num_faces_local))
2399)
2400) pgrid%uface_localids = -1
2401) pgrid%uface_nverts = -1
2402) pgrid%uface_left_natcellids = -1
2403) pgrid%uface_right_natcellids = -1
2404)
2405) call VecGetArrayF90(ghosted_cv_proc_rank, v_loc_p, ierr);CHKERRQ(ierr)
2406) pgrid%num_ufaces_local = 0
2407) pgrid%uface_nverts = 0
2408) do iface = 1, pgrid%ugrid_num_faces_local
2409) ghosted_id = ugrid%face_to_cell_ghosted(1, iface)
2410) dual_id = ugrid%face_to_cell_ghosted(2, iface)
2411) local_id = nG2L(ghosted_id)
2412)
2413) if (ghosted_id > ugrid%nlmax) cycle
2414)
2415) if (dual_id == 0) then
2416) pgrid%num_ufaces_local = pgrid%num_ufaces_local + 1
2417) pgridface_id = pgrid%ugridf2pgridf(iface)
2418) pgrid%uface_localids(pgrid%num_ufaces_local) = pgridface_id
2419) pgrid%uface_nverts(pgrid%num_ufaces_local) = pgrid%face_nverts(pgridface_id)
2420) pgrid%uface_left_natcellids(pgrid%num_ufaces_local) = nG2A(ghosted_id)
2421) pgrid%uface_right_natcellids(pgrid%num_ufaces_local) = 0
2422) pgrid%num_verts_of_ufaces_local = pgrid%num_verts_of_ufaces_local + pgrid%face_nverts(pgridface_id)
2423) else
2424) if (nG2L(dual_id) == 0) then
2425) if (v_loc_p(dual_id) >= option%myrank) then
2426) pgrid%num_ufaces_local = pgrid%num_ufaces_local + 1
2427) pgridface_id = pgrid%ugridf2pgridf(iface)
2428) pgrid%uface_localids(pgrid%num_ufaces_local) = pgridface_id
2429) pgrid%uface_nverts(pgrid%num_ufaces_local) = pgrid%face_nverts(pgridface_id)
2430) pgrid%uface_left_natcellids(pgrid%num_ufaces_local) = nG2A(ghosted_id)
2431) pgrid%uface_right_natcellids(pgrid%num_ufaces_local) = nG2A(dual_id)
2432) pgrid%num_verts_of_ufaces_local = pgrid%num_verts_of_ufaces_local + pgrid%face_nverts(pgridface_id)
2433) endif
2434) else
2435) pgrid%num_ufaces_local = pgrid%num_ufaces_local + 1
2436) pgridface_id = pgrid%ugridf2pgridf(iface)
2437) pgrid%uface_localids(pgrid%num_ufaces_local) = pgridface_id
2438) pgrid%uface_nverts(pgrid%num_ufaces_local) = pgrid%face_nverts(pgridface_id)
2439) pgrid%num_verts_of_ufaces_local = pgrid%num_verts_of_ufaces_local + pgrid%face_nverts(pgridface_id)
2440) pgrid%uface_left_natcellids(pgrid%num_ufaces_local) = nG2A(ghosted_id)
2441) pgrid%uface_right_natcellids(pgrid%num_ufaces_local) = nG2A(dual_id)
2442) endif
2443) endif
2444) enddo
2445) call VecRestoreArrayF90(ghosted_cv_proc_rank, v_loc_p, ierr);CHKERRQ(ierr)
2446)
2447) call VecDestroy(ghosted_cv_proc_rank, ierr);CHKERRQ(ierr)
2448)
2449)
2450) call MPI_Allreduce(pgrid%num_ufaces_local, pgrid%num_ufaces_global, &
2451) ONE_INTEGER_MPI, MPI_INTEGER, MPI_SUM, option%mycomm, ierr)
2452) call MPI_Allreduce(pgrid%num_verts_of_ufaces_local, pgrid%num_verts_of_ufaces_global, &
2453) ONE_INTEGER_MPI, MPI_INTEGER, MPI_SUM, option%mycomm, ierr)
2454)
2455) allocate(pgrid%uface_natvertids(pgrid%num_verts_of_ufaces_local))
2456) count = 0
2457) do iface = 1, pgrid%num_ufaces_local
2458)
2459) face_id = pgrid%uface_localids(iface)
2460)
2461) do ivertex = 1, pgrid%uface_nverts(iface)
2462) vertex_id = pgrid%face_vertids(ivertex, face_id)
2463) count = count + 1
2464) pgrid%uface_natvertids(count) = ugrid%vertex_ids_natural(vertex_id)
2465) enddo
2466) enddo
2467)
2468) end subroutine UGridPolyhedraComputeOutputInfo
2469)
2470) end module Grid_Unstructured_Polyhedra_module