grid_unstructured_explicit.F90 coverage: 100.00 %func 77.24 %block
1) module Grid_Unstructured_Explicit_module
2)
3) use Geometry_module
4) use Grid_Unstructured_Aux_module
5)
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 :: UGridExplicitRead, &
22) UGridExplicitDecompose, &
23) UGridExplicitSetInternConnect, &
24) UGridExplicitSetCellCentroids, &
25) UGridExplicitComputeVolumes, &
26) UGridExplicitSetBoundaryConnect, &
27) UGridExplicitSetConnections
28)
29) contains
30)
31) ! ************************************************************************** !
32)
33) subroutine UGridExplicitRead(unstructured_grid,filename,option)
34) !
35) ! Reads an explicit unstructured grid in parallel
36) !
37) ! Author: Glenn Hammond
38) ! Date: 10/03/12
39) !
40)
41) use Input_Aux_module
42) use Option_module
43) use String_module
44)
45) implicit none
46)
47) type(grid_unstructured_type) :: unstructured_grid
48) type(unstructured_explicit_type), pointer :: explicit_grid
49) character(len=MAXSTRINGLENGTH) :: filename
50) type(option_type) :: option
51)
52) type(input_type), pointer :: input
53) character(len=MAXSTRINGLENGTH) :: string, hint
54) character(len=MAXWORDLENGTH) :: word, card
55) PetscInt :: fileid, icell, iconn, irank, remainder, temp_int, num_to_read
56)
57) PetscInt :: num_cells, num_connections, num_elems
58) PetscInt :: num_cells_local, num_cells_local_save
59) PetscInt :: num_connections_local, num_connections_local_save
60) PetscInt :: num_elems_local, num_elems_local_save
61) PetscMPIInt :: status_mpi(MPI_STATUS_SIZE)
62) PetscMPIInt :: int_mpi
63) PetscErrorCode :: ierr
64) PetscReal, allocatable :: temp_real_array(:,:)
65) PetscInt, allocatable :: temp_int_array(:,:)
66) PetscInt :: ivertex, num_vertices
67)
68) explicit_grid => unstructured_grid%explicit_grid
69) ! Format of explicit unstructured grid file
70) ! id_, id_up_, id_dn_ = integer
71) ! x_, y_, z_, area_, volume_ = real
72) ! definitions
73) ! id_ = id of grid cell
74) ! id_up_ = id of upwind grid cell in connection
75) ! id_dn_ = id of downwind grid cell in connection
76) ! x_ = x coordinate of cell center
77) ! y_ = y coordinate of cell center
78) ! z_ = z coordinate of cell center
79) ! volume_ = volume of grid cell
80) ! -----------------------------------------------------------------
81) ! CELLS <integer> integer = # cells (N)
82) ! id_1 x_1 y_1 z_1 volume_1
83) ! id_2 x_2 y_2 z_2 volume_2
84) ! ...
85) ! ...
86) ! id_N x_N y_N z_N volume_N
87) ! CONNECTIONS <integer> integer = # connections (M)
88) ! id_up_1 id_dn_1 x_1 y_1 z_1 area_1
89) ! id_up_2 id_dn_2 x_2 y_2 z_2 area_2
90) ! ...
91) ! ...
92) ! id_up_M id_dn_M x_M y_M z_M area_M
93) ! -----------------------------------------------------------------
94)
95) if (option%myrank == option%io_rank) then
96)
97) fileid = 86
98) input => InputCreate(fileid,filename,option)
99)
100) call InputReadPflotranString(input,option)
101) ! read CELL card, though we already know the
102) call InputReadWord(input,option,card,PETSC_TRUE)
103) word = 'CELLS'
104) call InputErrorMsg(input,option,word,card)
105) if (.not.StringCompare(word,card)) then
106) option%io_buffer = 'Unrecognized keyword "' // trim(card) // &
107) '" in explicit grid file.'
108) call printErrMsgByRank(option)
109) endif
110)
111) hint = 'Explicit Unstructured Grid CELLS'
112) call InputReadInt(input,option,temp_int)
113) call InputErrorMsg(input,option,'number of cells',hint)
114) endif
115)
116) call MPI_Bcast(temp_int,ONE_INTEGER_MPI,MPI_INTEGER,option%io_rank, &
117) option%mycomm,ierr)
118) num_cells = temp_int
119) explicit_grid%num_cells_global = num_cells
120)
121) ! divide cells across ranks
122) num_cells_local = num_cells/option%mycommsize
123) num_cells_local_save = num_cells_local
124) remainder = num_cells - &
125) num_cells_local*option%mycommsize
126) if (option%myrank < remainder) num_cells_local = &
127) num_cells_local + 1
128)
129) allocate(explicit_grid%cell_ids(num_cells_local))
130) explicit_grid%cell_ids = 0
131) allocate(explicit_grid%cell_volumes(num_cells_local))
132) explicit_grid%cell_volumes = 0
133) allocate(explicit_grid%cell_centroids(num_cells_local))
134) do icell = 1, num_cells_local
135) explicit_grid%cell_centroids(icell)%x = 0.d0
136) explicit_grid%cell_centroids(icell)%y = 0.d0
137) explicit_grid%cell_centroids(icell)%z = 0.d0
138) enddo
139)
140) ! for now, read all cells from ASCII file through io_rank and communicate
141) ! to other ranks
142) if (option%myrank == option%io_rank) then
143) allocate(temp_real_array(5,num_cells_local_save+1))
144) ! read for other processors
145) do irank = 0, option%mycommsize-1
146) temp_real_array = UNINITIALIZED_DOUBLE
147) num_to_read = num_cells_local_save
148) if (irank < remainder) num_to_read = num_to_read + 1
149) do icell = 1, num_to_read
150) call InputReadPflotranString(input,option)
151) call InputReadStringErrorMsg(input,option,hint)
152) call InputReadInt(input,option,temp_int)
153) call InputErrorMsg(input,option,'cell id',hint)
154) temp_real_array(1,icell) = dble(temp_int)
155) call InputReadDouble(input,option,temp_real_array(2,icell))
156) call InputErrorMsg(input,option,'cell x coordinate',hint)
157) call InputReadDouble(input,option,temp_real_array(3,icell))
158) call InputErrorMsg(input,option,'cell y coordinate',hint)
159) call InputReadDouble(input,option,temp_real_array(4,icell))
160) call InputErrorMsg(input,option,'cell z coordinate',hint)
161) call InputReadDouble(input,option,temp_real_array(5,icell))
162) call InputErrorMsg(input,option,'cell volume',hint)
163) enddo
164)
165) ! if the cells reside on io_rank
166) if (irank == option%io_rank) then
167) #if UGRID_DEBUG
168) write(string,*) num_cells_local
169) string = trim(adjustl(string)) // ' cells stored on p0'
170) print *, trim(string)
171) #endif
172) do icell = 1, num_cells_local
173) explicit_grid%cell_ids(icell) = int(temp_real_array(1,icell))
174) explicit_grid%cell_centroids(icell)%x = temp_real_array(2,icell)
175) explicit_grid%cell_centroids(icell)%y = temp_real_array(3,icell)
176) explicit_grid%cell_centroids(icell)%z = temp_real_array(4,icell)
177) explicit_grid%cell_volumes(icell) = temp_real_array(5,icell)
178) enddo
179) else
180) ! otherwise communicate to other ranks
181) #if UGRID_DEBUG
182) write(string,*) num_to_read
183) write(word,*) irank
184) string = trim(adjustl(string)) // ' cells sent from p0 to p' // &
185) trim(adjustl(word))
186) print *, trim(string)
187) #endif
188) int_mpi = num_to_read*5
189) call MPI_Send(temp_real_array,int_mpi,MPI_DOUBLE_PRECISION,irank, &
190) num_to_read,option%mycomm,ierr)
191) endif
192) enddo
193) else
194) ! other ranks post the recv
195) #if UGRID_DEBUG
196) write(string,*) num_cells_local
197) write(word,*) option%myrank
198) string = trim(adjustl(string)) // ' cells received from p0 at p' // &
199) trim(adjustl(word))
200) print *, trim(string)
201) #endif
202) allocate(temp_real_array(5,num_cells_local))
203) int_mpi = num_cells_local*5
204) call MPI_Recv(temp_real_array,int_mpi, &
205) MPI_DOUBLE_PRECISION,option%io_rank, &
206) MPI_ANY_TAG,option%mycomm,status_mpi,ierr)
207) do icell = 1, num_cells_local
208) explicit_grid%cell_ids(icell) = int(temp_real_array(1,icell))
209) explicit_grid%cell_centroids(icell)%x = temp_real_array(2,icell)
210) explicit_grid%cell_centroids(icell)%y = temp_real_array(3,icell)
211) explicit_grid%cell_centroids(icell)%z = temp_real_array(4,icell)
212) explicit_grid%cell_volumes(icell) = temp_real_array(5,icell)
213) enddo
214)
215) endif
216) deallocate(temp_real_array)
217)
218) if (option%myrank == option%io_rank) then
219)
220)
221) call InputReadPflotranString(input,option)
222) ! read CONNECTIONS card, though we already know the
223) call InputReadWord(input,option,card,PETSC_TRUE)
224) word = 'CONNECTIONS'
225) call InputErrorMsg(input,option,word,card)
226) if (.not.StringCompare(word,card)) then
227) option%io_buffer = 'Unrecognized keyword "' // trim(card) // &
228) '" in explicit grid file.'
229) call printErrMsgByRank(option)
230) endif
231)
232) hint = 'Explicit Unstructured Grid CONNECTIONS'
233) call InputReadInt(input,option,temp_int)
234) call InputErrorMsg(input,option,'number of connections',hint)
235) endif
236)
237) int_mpi = 1
238) call MPI_Bcast(temp_int,ONE_INTEGER_MPI,MPI_INTEGER,option%io_rank, &
239) option%mycomm,ierr)
240) num_connections = temp_int
241)
242) ! divide cells across ranks
243) num_connections_local = num_connections/option%mycommsize
244) num_connections_local_save = num_connections_local
245) remainder = num_connections - &
246) num_connections_local*option%mycommsize
247) if (option%myrank < remainder) num_connections_local = &
248) num_connections_local + 1
249)
250) allocate(explicit_grid%connections(2,num_connections_local))
251) explicit_grid%connections = 0
252) allocate(explicit_grid%face_areas(num_connections_local))
253) explicit_grid%face_areas = 0
254) allocate(explicit_grid%face_centroids(num_connections_local))
255) do iconn = 1, num_connections_local
256) explicit_grid%face_centroids(iconn)%x = 0.d0
257) explicit_grid%face_centroids(iconn)%y = 0.d0
258) explicit_grid%face_centroids(iconn)%z = 0.d0
259) enddo
260)
261) ! for now, read all cells from ASCII file through io_rank and communicate
262) ! to other ranks
263) if (option%myrank == option%io_rank) then
264) allocate(temp_real_array(6,num_connections_local_save+1))
265) ! read for other processors
266) do irank = 0, option%mycommsize-1
267) temp_real_array = UNINITIALIZED_DOUBLE
268) num_to_read = num_connections_local_save
269) if (irank < remainder) num_to_read = num_to_read + 1
270) do iconn = 1, num_to_read
271) call InputReadPflotranString(input,option)
272) call InputReadStringErrorMsg(input,option,hint)
273) call InputReadInt(input,option,temp_int)
274) call InputErrorMsg(input,option,'cell id upwind',hint)
275) temp_real_array(1,iconn) = dble(temp_int)
276) call InputReadInt(input,option,temp_int)
277) call InputErrorMsg(input,option,'cell id downwind',hint)
278) temp_real_array(2,iconn) = dble(temp_int)
279) call InputReadDouble(input,option,temp_real_array(3,iconn))
280) call InputErrorMsg(input,option,'face x coordinate',hint)
281) call InputReadDouble(input,option,temp_real_array(4,iconn))
282) call InputErrorMsg(input,option,'face y coordinate',hint)
283) call InputReadDouble(input,option,temp_real_array(5,iconn))
284) call InputErrorMsg(input,option,'face z coordinate',hint)
285) call InputReadDouble(input,option,temp_real_array(6,iconn))
286) call InputErrorMsg(input,option,'face area',hint)
287) enddo
288)
289) ! if the cells reside on io_rank
290) if (irank == option%io_rank) then
291) #if UGRID_DEBUG
292) write(string,*) num_connections_local
293) string = trim(adjustl(string)) // ' connections stored on p0'
294) print *, trim(string)
295) #endif
296) do iconn = 1, num_connections_local
297) explicit_grid%connections(1,iconn) = int(temp_real_array(1,iconn))
298) explicit_grid%connections(2,iconn) = int(temp_real_array(2,iconn))
299) explicit_grid%face_centroids(iconn)%x = temp_real_array(3,iconn)
300) explicit_grid%face_centroids(iconn)%y = temp_real_array(4,iconn)
301) explicit_grid%face_centroids(iconn)%z = temp_real_array(5,iconn)
302) explicit_grid%face_areas(iconn) = temp_real_array(6,iconn)
303) enddo
304) else
305) ! otherwise communicate to other ranks
306) #if UGRID_DEBUG
307) write(string,*) num_to_read
308) write(word,*) irank
309) string = trim(adjustl(string)) // ' connections sent from p0 to p' // &
310) trim(adjustl(word))
311) print *, trim(string)
312) #endif
313) int_mpi = num_to_read*6
314) call MPI_Send(temp_real_array,int_mpi,MPI_DOUBLE_PRECISION,irank, &
315) num_to_read,option%mycomm,ierr)
316) endif
317) enddo
318) else
319) ! other ranks post the recv
320) #if UGRID_DEBUG
321) write(string,*) num_connections_local
322) write(word,*) option%myrank
323) string = trim(adjustl(string)) // ' connections received from p0 at p' // &
324) trim(adjustl(word))
325) print *, trim(string)
326) #endif
327) allocate(temp_real_array(6,num_connections_local))
328) int_mpi = num_connections_local*6
329) call MPI_Recv(temp_real_array,int_mpi, &
330) MPI_DOUBLE_PRECISION,option%io_rank, &
331) MPI_ANY_TAG,option%mycomm,status_mpi,ierr)
332) do iconn = 1, num_connections_local
333) explicit_grid%connections(1,iconn) = int(temp_real_array(1,iconn))
334) explicit_grid%connections(2,iconn) = int(temp_real_array(2,iconn))
335) explicit_grid%face_centroids(iconn)%x = temp_real_array(3,iconn)
336) explicit_grid%face_centroids(iconn)%y = temp_real_array(4,iconn)
337) explicit_grid%face_centroids(iconn)%z = temp_real_array(5,iconn)
338) explicit_grid%face_areas(iconn) = temp_real_array(6,iconn)
339) enddo
340)
341) endif
342) deallocate(temp_real_array)
343)
344) if (option%myrank == option%io_rank) then
345) call InputReadPflotranString(input,option)
346) ! read ELEMENTS card, we only use this for tecplot output
347) ! not used while solving the PDEs
348) call InputReadWord(input,option,card,PETSC_TRUE)
349) word = 'ELEMENTS'
350) if (.not.StringCompare(word,card)) return
351) card = 'Explicit Unstruct. Grid ELEMENTS'
352) call InputReadInt(input,option,num_elems)
353) call InputErrorMsg(input,option,'number of elements',card)
354) explicit_grid%num_elems = num_elems
355) unstructured_grid%max_nvert_per_cell = 8 ! Initial guess
356) allocate(explicit_grid%cell_connectivity(0:unstructured_grid% &
357) max_nvert_per_cell,num_elems))
358) do iconn = 1, num_elems
359) call InputReadPflotranString(input,option)
360) call InputReadStringErrorMsg(input,option,card)
361) call InputReadWord(input,option,word,PETSC_TRUE)
362) call InputErrorMsg(input,option,'element_type',card)
363) call StringtoUpper(word)
364) select case (word)
365) case('H')
366) num_vertices = 8
367) case('W')
368) num_vertices = 6
369) case('P')
370) num_vertices = 5
371) case('T')
372) num_vertices = 4
373) case('Q')
374) num_vertices = 4
375) case('TRI')
376) num_vertices = 3
377) end select
378) explicit_grid%cell_connectivity(0,iconn) = num_vertices
379) do ivertex = 1, num_vertices
380) call InputReadInt(input,option,explicit_grid%cell_connectivity(ivertex,iconn))
381) call InputErrorMsg(input,option,'vertex id',hint)
382) enddo
383) enddo
384) call InputReadPflotranString(input,option)
385) ! read VERTICES card, not used for calcuations, only tecplot output
386) call InputReadWord(input,option,card,PETSC_TRUE)
387) word = 'VERTICES'
388) call InputErrorMsg(input,option,word,card)
389) if (.not.StringCompare(word,card)) then
390) option%io_buffer = 'Unrecognized keyword "' // trim(card) // &
391) '" in explicit grid file.'
392) call printErrMsgByRank(option)
393) endif
394) allocate(explicit_grid%vertex_coordinates(explicit_grid%num_cells_global))
395) do icell = 1, explicit_grid%num_cells_global
396) explicit_grid%vertex_coordinates(icell)%x = 0.d0
397) explicit_grid%vertex_coordinates(icell)%y = 0.d0
398) explicit_grid%vertex_coordinates(icell)%z = 0.d0
399) enddo
400) do icell = 1, explicit_grid%num_cells_global
401) call InputReadPflotranString(input,option)
402) call InputReadStringErrorMsg(input,option,card)
403) call InputReadDouble(input,option, &
404) explicit_grid%vertex_coordinates(icell)%x)
405) call InputErrorMsg(input,option,'vertex 1',card)
406) call InputReadDouble(input,option, &
407) explicit_grid%vertex_coordinates(icell)%y)
408) call InputErrorMsg(input,option,'vertex 2',card)
409) call InputReadDouble(input,option, &
410) explicit_grid%vertex_coordinates(icell)%z)
411) call InputErrorMsg(input,option,'vertex 3',card)
412) enddo
413) endif
414)
415) if (option%myrank == option%io_rank) then
416) call InputDestroy(input)
417) endif
418)
419) end subroutine UGridExplicitRead
420)
421) ! ************************************************************************** !
422)
423) subroutine UGridExplicitDecompose(ugrid,option)
424) !
425) ! Decomposes an explicit unstructured grid across
426) ! ranks
427) !
428) ! Author: Glenn Hammond
429) ! Date: 05/17/12
430) !
431)
432) use Option_module
433) use Utility_module, only: reallocateIntArray, SearchOrderedArray
434)
435) implicit none
436)
437) #include "petsc/finclude/petscvec.h"
438) #include "petsc/finclude/petscvec.h90"
439) #include "petsc/finclude/petscmat.h"
440) #include "petsc/finclude/petscmat.h90"
441) #include "petsc/finclude/petscdm.h"
442) #include "petsc/finclude/petscdm.h90"
443) #include "petsc/finclude/petscis.h"
444) #include "petsc/finclude/petscis.h90"
445) #include "petsc/finclude/petscviewer.h"
446)
447) type(grid_unstructured_type) :: ugrid
448) type(option_type) :: option
449)
450) type(unstructured_explicit_type), pointer :: explicit_grid
451) PetscViewer :: viewer
452)
453) Mat :: M_mat,M_mat_loc
454) Vec :: M_vec
455) Mat :: Adj_mat
456) Mat :: Dual_mat
457) MatPartitioning :: Part
458) IS :: is_new
459) IS :: is_scatter
460) IS :: is_gather
461) PetscInt :: num_cells_local_new, num_cells_local_old
462) Vec :: cells_old, cells_local
463) Vec :: connections_old, connections_local
464) VecScatter :: vec_scatter
465)
466) PetscInt :: global_offset_old
467) PetscInt :: global_offset_new
468) PetscInt :: ghosted_id
469) PetscInt, allocatable :: local_connections(:), local_connection_offsets(:)
470) PetscInt, allocatable :: local_connections2(:), local_connection_offsets2(:)
471) PetscInt, allocatable :: int_array(:), int_array2(:), int_array3(:)
472) PetscInt, allocatable :: int_array4(:)
473) PetscInt, allocatable :: int_array2d(:,:)
474) PetscInt :: num_connections_local_old, num_connections_local
475) PetscInt :: num_connections_total
476) PetscInt :: num_connections_global, global_connection_offset
477) PetscInt :: id_up, id_dn, iconn, icell, count, offset
478) PetscInt :: conn_id, dual_id
479) PetscBool :: found
480) PetscInt :: i, temp_int, idual
481) PetscReal :: temp_real
482)
483) PetscInt :: iflag
484) PetscBool :: success
485) PetscInt, pointer :: ia_ptr(:), ja_ptr(:)
486) PetscInt, pointer :: ia_ptr2(:), ja_ptr2(:)
487) PetscReal, pointer :: vec_ptr(:), vec_ptr2(:)
488) PetscInt :: num_rows, num_cols, istart, iend, icol
489) PetscInt :: cell_stride, dual_offset, connection_offset, connection_stride
490) PetscInt :: natural_id_offset
491) PetscErrorCode :: ierr
492) PetscInt :: icell_up,icell_dn
493)
494) character(len=MAXSTRINGLENGTH) :: string
495)
496) explicit_grid => ugrid%explicit_grid
497)
498) #if UGRID_DEBUG
499) call printMsg(option,'Adjacency matrix')
500) #endif
501)
502) num_cells_local_old = size(explicit_grid%cell_ids)
503)
504) call MPI_Allreduce(num_cells_local_old,ugrid%nmax, &
505) ONE_INTEGER_MPI,MPIU_INTEGER, &
506) MPI_SUM,option%mycomm,ierr)
507)
508) ! determine the global offset from 0 for cells on this rank
509) global_offset_old = 0
510) call MPI_Exscan(num_cells_local_old,global_offset_old, &
511) ONE_INTEGER_MPI,MPIU_INTEGER,MPI_SUM,option%mycomm,ierr)
512)
513) num_connections_local_old = size(explicit_grid%connections,2)
514)
515) call MPI_Allreduce(num_connections_local_old,num_connections_global, &
516) ONE_INTEGER_MPI,MPIU_INTEGER,MPI_SUM,option%mycomm,ierr)
517)
518) global_connection_offset = 0
519) call MPI_Exscan(num_connections_local_old,global_connection_offset, &
520) ONE_INTEGER_MPI,MPIU_INTEGER,MPI_SUM,option%mycomm,ierr)
521)
522) call VecCreateMPI(option%mycomm,num_cells_local_old,ugrid%nmax, &
523) M_vec,ierr);CHKERRQ(ierr)
524) call VecZeroEntries(M_vec,ierr);CHKERRQ(ierr)
525) do iconn = 1, num_connections_local_old
526) do i = 1, 2
527) icell = explicit_grid%connections(i,iconn)-1
528) call VecSetValue(M_vec,icell,1.d0,ADD_VALUES,ierr);CHKERRQ(ierr)
529) enddo
530) enddo
531) call VecAssemblyBegin(M_vec,ierr);CHKERRQ(ierr)
532) call VecAssemblyEnd(M_vec,ierr);CHKERRQ(ierr)
533)
534) #if UGRID_DEBUG
535) call PetscViewerASCIIOpen(option%mycomm,'M_vec.out',viewer, &
536) ierr);CHKERRQ(ierr)
537) call VecView(M_vec,viewer,ierr);CHKERRQ(ierr)
538) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
539) #endif
540)
541) call VecMax(M_vec,PETSC_NULL_INTEGER,temp_real,ierr);CHKERRQ(ierr)
542) call VecDestroy(M_vec,ierr);CHKERRQ(ierr)
543) ugrid%max_ndual_per_cell = int(temp_real+0.1d0)
544) call MatCreateAIJ(option%mycomm,num_cells_local_old,PETSC_DECIDE, &
545) ugrid%nmax,num_connections_global, &
546) ugrid%max_ndual_per_cell,PETSC_NULL_INTEGER, &
547) ugrid%max_ndual_per_cell,PETSC_NULL_INTEGER, &
548) M_mat,ierr);CHKERRQ(ierr)
549) call MatZeroEntries(M_mat,ierr);CHKERRQ(ierr)
550) do iconn = 1, num_connections_local_old
551) temp_int = iconn+global_connection_offset-1
552) do i = 1, 2
553) icell = explicit_grid%connections(i,iconn)-1
554) call MatSetValue(M_mat,icell,temp_int,1.d0,INSERT_VALUES, &
555) ierr);CHKERRQ(ierr)
556) enddo
557) enddo
558) call MatAssemblyBegin(M_mat,MAT_FINAL_ASSEMBLY,ierr);CHKERRQ(ierr)
559) call MatAssemblyEnd(M_mat,MAT_FINAL_ASSEMBLY,ierr);CHKERRQ(ierr)
560)
561) #if UGRID_DEBUG
562) call PetscViewerASCIIOpen(option%mycomm,'M_mat.out',viewer, &
563) ierr);CHKERRQ(ierr)
564) call MatView(M_mat,viewer,ierr);CHKERRQ(ierr)
565) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
566) #endif
567)
568) ! GB: When MatConvert() is used, the diagonal entries are lost in Adj_mat
569) !call MatConvert(M_mat,MATMPIADJ,MAT_INITIAL_MATRIX,Adj_mat,ierr)
570) !call MatDestroy(M_mat,ierr)
571)
572) ! Alternate method of creating Adj_mat
573) if (option%mycommsize>1) then
574) call MatMPIAIJGetLocalMat(M_mat,MAT_INITIAL_MATRIX,M_mat_loc, &
575) ierr);CHKERRQ(ierr)
576) call MatGetRowIJF90(M_mat_loc,ZERO_INTEGER,PETSC_FALSE,PETSC_FALSE,num_rows, &
577) ia_ptr,ja_ptr,success,ierr);CHKERRQ(ierr)
578) else
579) call MatGetRowIJF90(M_mat,ZERO_INTEGER,PETSC_FALSE,PETSC_FALSE,num_rows, &
580) ia_ptr,ja_ptr,success,ierr);CHKERRQ(ierr)
581) endif
582)
583) count=0
584) do icell = 1,num_rows
585) istart = ia_ptr(icell)
586) iend = ia_ptr(icell+1)-1
587) num_cols = iend-istart+1
588) count = count+num_cols
589) enddo
590) allocate(local_connections(count))
591) allocate(local_connection_offsets(num_rows+1))
592) local_connection_offsets(1:num_rows+1) = ia_ptr(1:num_rows+1)
593) local_connections(1:count) = ja_ptr(1:count)
594)
595) call MatCreateMPIAdj(option%mycomm,num_cells_local_old, &
596) num_connections_global, &
597) local_connection_offsets, &
598) local_connections,PETSC_NULL_INTEGER,Adj_mat, &
599) ierr);CHKERRQ(ierr)
600)
601) if (option%mycommsize>1) then
602) call MatRestoreRowIJF90(M_mat_loc,ZERO_INTEGER,PETSC_FALSE,PETSC_FALSE,num_rows, &
603) ia_ptr,ja_ptr,success,ierr);CHKERRQ(ierr)
604) else
605) call MatRestoreRowIJF90(M_mat,ZERO_INTEGER,PETSC_FALSE,PETSC_FALSE,num_rows, &
606) ia_ptr,ja_ptr,success,ierr);CHKERRQ(ierr)
607) endif
608) call MatDestroy(M_mat,ierr);CHKERRQ(ierr)
609)
610) #if UGRID_DEBUG
611) call PetscViewerASCIIOpen(option%mycomm,'Adj.out',viewer,ierr);CHKERRQ(ierr)
612) call MatView(Adj_mat,viewer,ierr);CHKERRQ(ierr)
613) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
614) #endif
615) ! call printErrMsg(option,'debugg')
616)
617) ! Create the Dual matrix.
618) call MatCreateAIJ(option%mycomm,num_cells_local_old,PETSC_DECIDE, &
619) ugrid%nmax,ugrid%nmax, &
620) ugrid%max_ndual_per_cell,PETSC_NULL_INTEGER, &
621) ugrid%max_ndual_per_cell,PETSC_NULL_INTEGER, &
622) M_mat,ierr);CHKERRQ(ierr)
623) do iconn = 1, num_connections_local_old
624) icell_up = explicit_grid%connections(1,iconn)-1
625) icell_dn = explicit_grid%connections(2,iconn)-1
626) call MatSetValue(M_mat,icell_up,icell_dn,1.d0,INSERT_VALUES, &
627) ierr);CHKERRQ(ierr)
628) call MatSetValue(M_mat,icell_dn,icell_up,1.d0,INSERT_VALUES, &
629) ierr);CHKERRQ(ierr)
630) enddo
631)
632) call MatAssemblyBegin(M_mat,MAT_FINAL_ASSEMBLY,ierr);CHKERRQ(ierr)
633) call MatAssemblyEnd(M_mat,MAT_FINAL_ASSEMBLY,ierr);CHKERRQ(ierr)
634)
635) !call MatConvert(M_mat,MATMPIADJ,MAT_INITIAL_MATRIX,Dual_mat,ierr)
636) !call MatDestroy(M_mat,ierr)
637)
638) ! Alternate method of creating Dual_mat
639) if (option%mycommsize>1) then
640) call MatMPIAIJGetLocalMat(M_mat,MAT_INITIAL_MATRIX,M_mat_loc, &
641) ierr);CHKERRQ(ierr)
642) call MatGetRowIJF90(M_mat_loc,ZERO_INTEGER,PETSC_FALSE,PETSC_FALSE,num_rows, &
643) ia_ptr,ja_ptr,success,ierr);CHKERRQ(ierr)
644) else
645) call MatGetRowIJF90(M_mat,ZERO_INTEGER,PETSC_FALSE,PETSC_FALSE,num_rows, &
646) ia_ptr,ja_ptr,success,ierr);CHKERRQ(ierr)
647) endif
648)
649) count=0
650) do icell = 1,num_rows
651) istart = ia_ptr(icell)
652) iend = ia_ptr(icell+1)-1
653) num_cols = iend-istart+1
654) count = count+num_cols
655) enddo
656) allocate(local_connections2(count))
657) allocate(local_connection_offsets2(num_rows+1))
658) local_connection_offsets2(1:num_rows+1) = ia_ptr(1:num_rows+1)
659) local_connections2(1:count) = ja_ptr(1:count)
660)
661) call MatCreateMPIAdj(option%mycomm,num_cells_local_old, &
662) ugrid%nmax, &
663) local_connection_offsets2, &
664) local_connections2,PETSC_NULL_INTEGER,Dual_mat, &
665) ierr);CHKERRQ(ierr)
666)
667) if (option%mycommsize>1) then
668) call MatRestoreRowIJF90(M_mat_loc,ZERO_INTEGER,PETSC_FALSE,PETSC_FALSE,num_rows, &
669) ia_ptr,ja_ptr,success,ierr);CHKERRQ(ierr)
670) else
671) call MatRestoreRowIJF90(M_mat,ZERO_INTEGER,PETSC_FALSE,PETSC_FALSE,num_rows, &
672) ia_ptr,ja_ptr,success,ierr);CHKERRQ(ierr)
673) endif
674) call MatDestroy(M_mat,ierr);CHKERRQ(ierr)
675)
676) #if UGRID_DEBUG
677) call PetscViewerASCIIOpen(option%mycomm,'Dual_mat.out',viewer, &
678) ierr);CHKERRQ(ierr)
679) call MatView(Dual_mat,viewer,ierr);CHKERRQ(ierr)
680) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
681) #endif
682)
683) call UGridPartition(ugrid,option,Dual_mat,is_new, &
684) num_cells_local_new)
685)
686) ! second argument of ZERO_INTEGER means to use 0-based indexing
687) ! MagGetRowIJF90 returns row and column pointers for compressed matrix data
688) call MatGetRowIJF90(Dual_mat,ZERO_INTEGER,PETSC_FALSE,PETSC_FALSE,num_rows, &
689) ia_ptr,ja_ptr,success,ierr);CHKERRQ(ierr)
690)
691) if (.not.success .or. num_rows /= num_cells_local_old) then
692) print *, option%myrank, num_rows, success, num_cells_local_old
693) option%io_buffer = 'Error getting IJ row indices from dual matrix'
694) call printErrMsg(option)
695) endif
696)
697) call MatRestoreRowIJF90(Dual_mat,ZERO_INTEGER,PETSC_FALSE,PETSC_FALSE, &
698) num_rows,ia_ptr,ja_ptr,success,ierr);CHKERRQ(ierr)
699)
700) ! in order to redistributed cell/connection data among ranks, I package it
701) ! in a crude way within a strided petsc vec and pass it. The stride
702) ! determines the size of each cells "packaged" data
703) connection_offset = 6 + 1 ! +1 for -777
704) dual_offset = connection_offset + ugrid%max_ndual_per_cell + 1 ! +1 for -888
705) cell_stride = dual_offset + ugrid%max_ndual_per_cell + 1 ! +1 for -999999
706) natural_id_offset = 2
707)
708) ! Information for each cell is packed in a strided petsc vec
709) ! The information is ordered within each stride as follows:
710) ! -cell_N ! global cell id (negative indicates 1-based)
711) ! natural cell id
712) ! cell x coordinate
713) ! cell y coordinate
714) ! cell z coordinate
715) ! cell volume
716) ! -777 ! separator between cell info and connection info
717) ! conn1 ! connection ids between cell_N and others
718) ! conn1
719) ! ...
720) ! connN
721) ! -888 ! separator between connection info and dual ids
722) ! dual1 ! dual ids between cell_N and others
723) ! dual2
724) ! ...
725) ! dualN
726) ! -999999 ! separator indicating end of information for cell_N
727)
728) ! the purpose of -888, and -999999 is to allow one to use cells of
729) ! various geometry.
730)
731) call UGridCreateOldVec(ugrid,option,cells_old, &
732) num_cells_local_old, &
733) is_new,is_scatter,cell_stride)
734)
735) ! 0 = 0-based indexing
736) ! MagGetRowIJF90 returns row and column pointers for compressed matrix data
737) ! pointers to Dual mat
738) call MatGetRowIJF90(Dual_mat,ZERO_INTEGER,PETSC_FALSE,PETSC_FALSE,num_rows, &
739) ia_ptr,ja_ptr,success,ierr);CHKERRQ(ierr)
740) ! pointers to Adj mat
741) call MatGetRowIJF90(Adj_mat,ZERO_INTEGER,PETSC_FALSE,PETSC_FALSE,temp_int, &
742) ia_ptr2,ja_ptr2,success,ierr);CHKERRQ(ierr)
743)
744) if (num_rows /= temp_int) then
745) write(string,*) num_rows, temp_int
746) option%io_buffer = 'Number of rows in Adj and Dual matrices inconsistent:'
747) option%io_buffer = trim(option%io_buffer) // trim(adjustl(string))
748) call printErrMsgByRank(option)
749) endif
750)
751) call VecGetArrayF90(cells_old,vec_ptr,ierr);CHKERRQ(ierr)
752) count = 0
753) do icell = 1, num_cells_local_old
754) count = count + 1
755) ! set global cell id
756) ! negate to indicate cell id with 1-based numbering (-0 = 0)
757) vec_ptr(count) = -(global_offset_old+icell)
758) count = count + 1
759) vec_ptr(count) = explicit_grid%cell_ids(icell)
760) count = count + 1
761) vec_ptr(count) = explicit_grid%cell_centroids(icell)%x
762) count = count + 1
763) vec_ptr(count) = explicit_grid%cell_centroids(icell)%y
764) count = count + 1
765) vec_ptr(count) = explicit_grid%cell_centroids(icell)%z
766) count = count + 1
767) vec_ptr(count) = explicit_grid%cell_volumes(icell)
768) ! add the separator
769) count = count + 1
770) vec_ptr(count) = -777 ! help differentiate
771)
772) ! add the connection ids
773) istart = ia_ptr2(icell)
774) iend = ia_ptr2(icell+1)-1
775) num_cols = iend-istart+1
776) if (num_cols > ugrid%max_ndual_per_cell) then
777) option%io_buffer = &
778) 'Number of columns in Adj matrix is larger then max_ndual_per_cell.'
779) call printErrMsgByRank(option)
780) endif
781) do icol = 1, ugrid%max_ndual_per_cell
782) count = count + 1
783) if (icol <= num_cols) then
784) ! increment for 1-based ordering
785) vec_ptr(count) = ja_ptr2(icol+istart) + 1
786) else
787) vec_ptr(count) = 0
788) endif
789) enddo
790)
791) ! add the separator
792) count = count + 1
793) vec_ptr(count) = -888 ! help differentiate
794)
795) ! add the dual ids
796) istart = ia_ptr(icell)
797) iend = ia_ptr(icell+1)-1
798) num_cols = iend-istart+1
799) if (num_cols > ugrid%max_ndual_per_cell) then
800) option%io_buffer = &
801) 'Number of columns in Dual matrix is larger then max_ndual_per_cell.'
802) call printErrMsgByRank(option)
803) endif
804) do icol = 1, ugrid%max_ndual_per_cell
805) count = count + 1
806) if (icol <= num_cols) then
807) ! increment for 1-based ordering
808) vec_ptr(count) = ja_ptr(icol+istart) + 1
809) else
810) vec_ptr(count) = 0
811) endif
812) enddo
813)
814) ! final separator
815) count = count + 1
816) vec_ptr(count) = -999999 ! help differentiate
817) enddo
818) call VecRestoreArrayF90(cells_old,vec_ptr,ierr);CHKERRQ(ierr)
819)
820) ! pointers to Dual mat
821) call MatRestoreRowIJF90(Dual_mat,ZERO_INTEGER,PETSC_FALSE,PETSC_FALSE, &
822) num_rows,ia_ptr,ja_ptr,success,ierr);CHKERRQ(ierr)
823) ! pointers to Adj mat
824) call MatRestoreRowIJF90(Adj_mat,ZERO_INTEGER,PETSC_FALSE,PETSC_FALSE, &
825) temp_int,ia_ptr2,ja_ptr2,success,ierr);CHKERRQ(ierr)
826) call MatDestroy(Dual_mat,ierr);CHKERRQ(ierr)
827) call MatDestroy(Adj_mat,ierr);CHKERRQ(ierr)
828) deallocate(local_connections)
829) deallocate(local_connection_offsets)
830) deallocate(local_connections2)
831) deallocate(local_connection_offsets2)
832)
833)
834) ! is_scatter is destroyed within UGridNaturalToPetsc
835) call UGridNaturalToPetsc(ugrid,option, &
836) cells_old,cells_local, &
837) num_cells_local_new,cell_stride,dual_offset, &
838) natural_id_offset,is_scatter)
839)
840) call VecDestroy(cells_old,ierr);CHKERRQ(ierr)
841)
842) ! set up connections
843) connection_stride = 8
844) ! create strided vector with the old connection distribution
845) call VecCreate(option%mycomm,connections_old,ierr);CHKERRQ(ierr)
846) call VecSetSizes(connections_old, &
847) connection_stride*num_connections_local_old, &
848) PETSC_DECIDE,ierr);CHKERRQ(ierr)
849) call VecSetFromOptions(connections_old,ierr);CHKERRQ(ierr)
850)
851) call VecGetArrayF90(connections_old,vec_ptr,ierr);CHKERRQ(ierr)
852) do iconn = 1, num_connections_local_old
853) offset = (iconn-1)*connection_stride
854) vec_ptr(offset+1) = explicit_grid%connections(1,iconn)
855) vec_ptr(offset+2) = explicit_grid%connections(2,iconn)
856) vec_ptr(offset+3) = explicit_grid%face_centroids(iconn)%x
857) vec_ptr(offset+4) = explicit_grid%face_centroids(iconn)%y
858) vec_ptr(offset+5) = explicit_grid%face_centroids(iconn)%z
859) vec_ptr(offset+6) = explicit_grid%face_areas(iconn)
860) vec_ptr(offset+7) = 1.d0 ! flag for local connections
861) vec_ptr(offset+8) = -888.d0
862) enddo
863) call VecRestoreArrayF90(connections_old,vec_ptr,ierr);CHKERRQ(ierr)
864)
865)
866) #if UGRID_DEBUG
867) call PetscViewerASCIIOpen(option%mycomm,'connections_old.out',viewer, &
868) ierr);CHKERRQ(ierr)
869) call VecView(connections_old,viewer,ierr);CHKERRQ(ierr)
870) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
871) #endif
872)
873) ! count the number of cells and their duals
874) call VecGetArrayF90(cells_local,vec_ptr,ierr);CHKERRQ(ierr)
875) count = 0
876) do ghosted_id=1, ugrid%ngmax
877) do iconn = 1, ugrid%max_ndual_per_cell
878) conn_id = int(vec_ptr(iconn + connection_offset + (ghosted_id-1)*cell_stride))
879) if (conn_id < 1) exit ! here we hit the 0 at the end of last connection
880) ! yes, we will be counting them twice
881) count = count + 1
882) enddo
883) enddo
884) num_connections_total = count ! many of these are redundant and will be removed
885) ! allocate and fill an array with the natural cell and dual ids
886) allocate(int_array(num_connections_total))
887) count = 0
888) do ghosted_id=1, ugrid%ngmax
889) do iconn = 1, ugrid%max_ndual_per_cell
890) conn_id = int(vec_ptr(iconn + connection_offset + (ghosted_id-1)*cell_stride))
891) if (conn_id < 1) exit ! again we hit the 0
892) count = count + 1
893) int_array(count) = conn_id
894) enddo
895) enddo
896) call VecRestoreArrayF90(cells_local,vec_ptr,ierr);CHKERRQ(ierr)
897)
898) allocate(int_array2(num_connections_total))
899) do iconn = 1, num_connections_total
900) int_array2(iconn) = iconn
901) enddo
902)
903) ! sort connections - int_array2 will return the reordering while int_array
904) ! remains the same.
905) int_array2 = int_array2 - 1
906) call PetscSortIntWithPermutation(num_connections_total,int_array, &
907) int_array2,ierr);CHKERRQ(ierr)
908) int_array2 = int_array2 + 1
909)
910) ! permute, remove duplicate connections, and renumber to local ordering
911) allocate(int_array3(num_connections_total))
912) allocate(int_array4(num_connections_total))
913) int_array3 = UNINITIALIZED_INTEGER
914) int_array4 = UNINITIALIZED_INTEGER
915) int_array3(1) = int_array(int_array2(1))
916) count = 1
917) do iconn = 1, num_connections_total
918) if (int_array(int_array2(iconn)) > int_array3(count)) then
919) count = count + 1
920) int_array3(count) = int_array(int_array2(iconn))
921) endif
922) int_array4(int_array2(iconn)) = count
923) enddo
924) deallocate(int_array)
925) deallocate(int_array2)
926)
927) num_connections_local = count ! new compressed count
928) allocate(int_array(num_connections_local))
929) int_array = int_array3(1:num_connections_local)
930) deallocate(int_array3)
931)
932) ! replace original connections ids (naturally numbered) with locally
933) ! numbered connection ids (int_array4)
934) call VecGetArrayF90(cells_local,vec_ptr,ierr);CHKERRQ(ierr)
935) count = 0
936) do ghosted_id=1, ugrid%ngmax
937) do iconn = 1, ugrid%max_ndual_per_cell
938) conn_id = int(vec_ptr(iconn + connection_offset + (ghosted_id-1)*cell_stride))
939) if (conn_id < 1) exit ! again we hit the 0
940) count = count + 1
941) vec_ptr(iconn + connection_offset + (ghosted_id-1)*cell_stride) = &
942) int_array4(count)
943) enddo
944) enddo
945) deallocate(int_array4)
946) call VecRestoreArrayF90(cells_local,vec_ptr,ierr);CHKERRQ(ierr)
947)
948) ! check to ensure that the number before/after are consistent
949) if (count /= num_connections_total) then
950) write(string,'(2i6)') count, num_connections_total
951) option%io_buffer = 'Inconsistent values for num_connections_total: ' // &
952) trim(adjustl(string))
953) call printErrMsgByRank(option)
954) endif
955) num_connections_total = UNINITIALIZED_INTEGER ! set to uninitialized value to catch bugs
956)
957) call VecCreate(PETSC_COMM_SELF,connections_local,ierr);CHKERRQ(ierr)
958) call VecSetSizes(connections_local,num_connections_local*connection_stride, &
959) PETSC_DECIDE,ierr);CHKERRQ(ierr)
960) call VecSetBlockSize(connections_local,connection_stride,ierr);CHKERRQ(ierr)
961) call VecSetFromOptions(connections_local,ierr);CHKERRQ(ierr)
962)
963) int_array = int_array-1
964) call ISCreateBlock(option%mycomm,connection_stride,num_connections_local, &
965) int_array,PETSC_COPY_VALUES,is_scatter, &
966) ierr);CHKERRQ(ierr)
967) do iconn = 1, num_connections_local
968) int_array(iconn) = iconn-1
969) enddo
970) call ISCreateBlock(option%mycomm,connection_stride,num_connections_local, &
971) int_array,PETSC_COPY_VALUES,is_gather,ierr);CHKERRQ(ierr)
972) deallocate(int_array)
973)
974) #if UGRID_DEBUG
975) call PetscViewerASCIIOpen(option%mycomm,'is_scatter_conn_old_to_local.out',viewer, &
976) ierr);CHKERRQ(ierr)
977) call ISView(is_scatter,viewer,ierr);CHKERRQ(ierr)
978) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
979) call PetscViewerASCIIOpen(option%mycomm,'is_gather_conn_old_to_local.out',viewer, &
980) ierr);CHKERRQ(ierr)
981) call ISView(is_gather,viewer,ierr);CHKERRQ(ierr)
982) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
983) call printMsg(option,'Scatter/gathering local connection info')
984) #endif
985)
986) ! scatter all the connection data from the old to local
987) call VecScatterCreate(connections_old,is_scatter,connections_local, &
988) is_gather,vec_scatter,ierr);CHKERRQ(ierr)
989) call ISDestroy(is_gather,ierr);CHKERRQ(ierr)
990) call ISDestroy(is_scatter,ierr);CHKERRQ(ierr)
991) call VecScatterBegin(vec_scatter,connections_old,connections_local, &
992) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
993) call VecScatterEnd(vec_scatter,connections_old,connections_local, &
994) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
995) call VecScatterDestroy(vec_scatter,ierr);CHKERRQ(ierr)
996)
997)
998) #if UGRID_DEBUG
999) write(string,*) option%myrank
1000) string = 'connections_local_nat' // trim(adjustl(string)) // '.out'
1001) call PetscViewerASCIIOpen(PETSC_COMM_SELF,trim(string),viewer, &
1002) ierr);CHKERRQ(ierr)
1003) call VecView(connections_local,viewer,ierr);CHKERRQ(ierr)
1004) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
1005) #endif
1006)
1007) ! loop over cells and change the natural ids in the duals to local ids
1008) allocate(int_array2d(2,num_connections_local))
1009) int_array2d = UNINITIALIZED_INTEGER
1010) call VecGetArrayF90(cells_local,vec_ptr,ierr);CHKERRQ(ierr)
1011) do ghosted_id=1, ugrid%ngmax
1012) do iconn = 1, ugrid%max_ndual_per_cell
1013) ! this connection id is now local
1014) conn_id = int(vec_ptr(iconn + connection_offset + (ghosted_id-1)*cell_stride))
1015) if (conn_id < 1) exit ! again we hit the 0
1016) do i = 1, 2
1017) if (int_array2d(i,conn_id) <= UNINITIALIZED_INTEGER) then
1018) int_array2d(i,conn_id) = ghosted_id
1019) exit
1020) endif
1021) if (i > 2) then
1022) write(string,'(2i5)') ghosted_id, conn_id
1023) option%io_buffer = 'Too many local cells match connection: ' // &
1024) trim(adjustl(string))
1025) call printErrMsgByRank(option)
1026) endif
1027) enddo
1028) enddo
1029) enddo
1030) call VecRestoreArrayF90(cells_local,vec_ptr,ierr);CHKERRQ(ierr)
1031)
1032) ! map natural ids in connections to local ids
1033) ! negate connection ids as a flag
1034) int_array2d = -1*int_array2d
1035) call VecGetArrayF90(connections_local,vec_ptr,ierr);CHKERRQ(ierr)
1036) do iconn = 1, num_connections_local
1037) offset = connection_stride*(iconn-1)
1038) ! all values should be negative at this point, unless uninitialized
1039) if (maxval(int_array2d(:,iconn)) >= 999) then
1040) ! connection is between two ghosted cells
1041) vec_ptr(offset+7) = 0.d0
1042) cycle
1043) endif
1044) id_up = int(vec_ptr(offset+1)) ! this is the natural id
1045) id_dn = int(vec_ptr(offset+2))
1046) count = 0
1047) found = PETSC_FALSE
1048) do i = 1, 2
1049) if (ugrid%cell_ids_natural(abs(int_array2d(i,iconn))) == id_up) then
1050) int_array2d(i,iconn) = abs(int_array2d(i,iconn))
1051) found = PETSC_TRUE
1052) count = count + 1
1053) endif
1054) if (ugrid%cell_ids_natural(abs(int_array2d(i,iconn))) == id_dn) then
1055) int_array2d(i,iconn) = abs(int_array2d(i,iconn))
1056) found = PETSC_TRUE
1057) count = count - 1
1058) endif
1059) enddo
1060) ! count should be zero, meaning it found the upwind and downwind cell
1061) ! ids
1062) if (count /= 0 .or. .not.found) then
1063) write(string,*) iconn, id_up, id_dn
1064) if (.not.found) then
1065) option%io_buffer = 'upwind/downwind cells not found: '
1066) else if (count < 0) then
1067) option%io_buffer = 'upwind cell not found: '
1068) else
1069) option%io_buffer = 'downwind cell not found: '
1070) endif
1071) option%io_buffer = trim(option%io_buffer) // trim(adjustl(string))
1072) call printErrMsgByRank(option)
1073) endif
1074) id_up = int_array2d(1,iconn)
1075) id_dn = int_array2d(2,iconn)
1076) if (id_up < id_dn) then
1077) vec_ptr(offset+1) = id_up ! now local ids
1078) vec_ptr(offset+2) = id_dn
1079) else
1080) vec_ptr(offset+1) = id_dn
1081) vec_ptr(offset+2) = id_up
1082) endif
1083) if (id_up > ugrid%nlmax .and. id_dn > ugrid%nlmax) then
1084) ! connection is between two ghosted cells
1085) vec_ptr(offset+7) = 0.d0
1086) endif
1087) enddo
1088) call VecRestoreArrayF90(connections_local,vec_ptr,ierr);CHKERRQ(ierr)
1089) deallocate(int_array2d)
1090)
1091) #if UGRID_DEBUG
1092) write(string,*) option%myrank
1093) string = 'connections_local_loc' // trim(adjustl(string)) // '.out'
1094) call PetscViewerASCIIOpen(PETSC_COMM_SELF,trim(string),viewer, &
1095) ierr);CHKERRQ(ierr)
1096) call VecView(connections_local,viewer,ierr);CHKERRQ(ierr)
1097) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
1098) #endif
1099)
1100) ! deallocate/allocate grid cell info locally
1101) deallocate(explicit_grid%cell_ids)
1102) deallocate(explicit_grid%cell_volumes)
1103) deallocate(explicit_grid%cell_centroids)
1104)
1105) allocate(explicit_grid%cell_ids(ugrid%ngmax))
1106) explicit_grid%cell_ids = UNINITIALIZED_INTEGER
1107) allocate(explicit_grid%cell_volumes(ugrid%ngmax))
1108) explicit_grid%cell_volumes = UNINITIALIZED_DOUBLE
1109) allocate(explicit_grid%cell_centroids(ugrid%ngmax))
1110) do icell = 1, ugrid%ngmax
1111) explicit_grid%cell_centroids(icell)%x = UNINITIALIZED_DOUBLE
1112) explicit_grid%cell_centroids(icell)%y = UNINITIALIZED_DOUBLE
1113) explicit_grid%cell_centroids(icell)%z = UNINITIALIZED_DOUBLE
1114) enddo
1115)
1116) call VecGetArrayF90(cells_local,vec_ptr,ierr);CHKERRQ(ierr)
1117) do ghosted_id=1, ugrid%ngmax
1118) offset = cell_stride*(ghosted_id-1)
1119) explicit_grid%cell_ids(ghosted_id) = int(vec_ptr(offset + 2))
1120) explicit_grid%cell_centroids(ghosted_id)%x = vec_ptr(offset + 3)
1121) explicit_grid%cell_centroids(ghosted_id)%y = vec_ptr(offset + 4)
1122) explicit_grid%cell_centroids(ghosted_id)%z = vec_ptr(offset + 5)
1123) explicit_grid%cell_volumes(ghosted_id) = vec_ptr(offset + 6)
1124) enddo
1125) call VecRestoreArrayF90(cells_local,vec_ptr,ierr);CHKERRQ(ierr)
1126)
1127) #if UGRID_DEBUG
1128) write(string,*) option%myrank
1129) string = 'cells_local_raw' // trim(adjustl(string)) // '.out'
1130) open(unit=86,file=trim(string))
1131) do ghosted_id = 1, ugrid%ngmax
1132) write(86,'(i5,4f10.3)') explicit_grid%cell_ids(ghosted_id), &
1133) explicit_grid%cell_centroids(ghosted_id)%x, &
1134) explicit_grid%cell_centroids(ghosted_id)%y, &
1135) explicit_grid%cell_centroids(ghosted_id)%z, &
1136) explicit_grid%cell_volumes(ghosted_id)
1137) enddo
1138) close(86)
1139) #endif
1140)
1141) ! deallocate/allocate connection info locally
1142) deallocate(explicit_grid%connections)
1143) deallocate(explicit_grid%face_areas)
1144) deallocate(explicit_grid%face_centroids)
1145)
1146) count = 0
1147) call VecGetArrayF90(connections_local,vec_ptr,ierr);CHKERRQ(ierr)
1148) do iconn = 1, num_connections_local
1149) offset = connection_stride*(iconn-1)
1150) if (vec_ptr(offset+7) > 0.1d0) count = count + 1
1151) enddo
1152) call VecRestoreArrayF90(connections_local,vec_ptr,ierr);CHKERRQ(ierr)
1153)
1154) allocate(explicit_grid%connections(2,count))
1155) explicit_grid%connections = 0
1156) allocate(explicit_grid%face_areas(count))
1157) explicit_grid%face_areas = 0
1158) allocate(explicit_grid%face_centroids(count))
1159) do iconn = 1, count
1160) explicit_grid%face_centroids(iconn)%x = 0.d0
1161) explicit_grid%face_centroids(iconn)%y = 0.d0
1162) explicit_grid%face_centroids(iconn)%z = 0.d0
1163) enddo
1164) call VecGetArrayF90(connections_local,vec_ptr,ierr);CHKERRQ(ierr)
1165) count = 0
1166) do iconn = 1, num_connections_local
1167) offset = connection_stride*(iconn-1)
1168) if (vec_ptr(offset+7) > 0.1d0) then
1169) count = count + 1
1170) explicit_grid%connections(1,count) = int(vec_ptr(offset+1))
1171) explicit_grid%connections(2,count) = int(vec_ptr(offset+2))
1172) explicit_grid%face_centroids(count)%x = vec_ptr(offset+3)
1173) explicit_grid%face_centroids(count)%y = vec_ptr(offset+4)
1174) explicit_grid%face_centroids(count)%z = vec_ptr(offset+5)
1175) explicit_grid%face_areas(count) = vec_ptr(offset+6)
1176) endif
1177) enddo
1178) call VecRestoreArrayF90(connections_local,vec_ptr,ierr);CHKERRQ(ierr)
1179) num_connections_local = count
1180)
1181) #if UGRID_DEBUG
1182) write(string,*) option%myrank
1183) string = 'connections_local_raw' // trim(adjustl(string)) // '.out'
1184) open(unit=86,file=trim(string))
1185) do iconn = 1, num_connections_local
1186) write(86,'(2i5,4f7.3)') explicit_grid%connections(1,iconn), &
1187) explicit_grid%connections(2,iconn), &
1188) explicit_grid%face_centroids(iconn)%x, &
1189) explicit_grid%face_centroids(iconn)%y, &
1190) explicit_grid%face_centroids(iconn)%z, &
1191) explicit_grid%face_areas(iconn)
1192) enddo
1193) close(86)
1194) #endif
1195)
1196) call VecDestroy(connections_old,ierr);CHKERRQ(ierr)
1197) call VecDestroy(connections_local,ierr);CHKERRQ(ierr)
1198) call VecDestroy(cells_local,ierr);CHKERRQ(ierr)
1199)
1200) end subroutine UGridExplicitDecompose
1201)
1202) ! ************************************************************************** !
1203)
1204) subroutine UGridExplicitSetCellCentroids(explicit_grid,x,y,z, &
1205) x_min,x_max,y_min,y_max,z_min,z_max)
1206) !
1207) ! Sets the centroid of each grid cell
1208) !
1209) ! Author: Glenn Hammond
1210) ! Date: 05/17/12
1211) !
1212)
1213) use Option_module
1214)
1215) implicit none
1216)
1217) type(unstructured_explicit_type) :: explicit_grid
1218) PetscReal :: x(:), y(:), z(:)
1219) PetscReal :: x_min, x_max, y_min, y_max, z_min, z_max
1220)
1221) PetscInt :: icell
1222)
1223) do icell = 1, size(explicit_grid%cell_centroids)
1224) x(icell) = explicit_grid%cell_centroids(icell)%x
1225) y(icell) = explicit_grid%cell_centroids(icell)%y
1226) z(icell) = explicit_grid%cell_centroids(icell)%z
1227) enddo
1228)
1229) x_min = minval(x)
1230) x_max = maxval(x)
1231) y_min = minval(y)
1232) y_max = maxval(y)
1233) z_min = minval(z)
1234) z_max = maxval(z)
1235)
1236) end subroutine UGridExplicitSetCellCentroids
1237)
1238) ! ************************************************************************** !
1239)
1240) function UGridExplicitSetInternConnect(explicit_grid,option)
1241) !
1242) ! Sets up the internal connectivity within
1243) ! the connectivity object
1244) !
1245) ! Author: Glenn Hammond
1246) ! Date: 05/17/12
1247) !
1248)
1249) use Utility_module
1250) use Connection_module
1251) use Option_module
1252)
1253) implicit none
1254)
1255) type(connection_set_type), pointer :: UGridExplicitSetInternConnect
1256)
1257) type(unstructured_explicit_type) :: explicit_grid
1258) type(option_type) :: option
1259)
1260) type(connection_set_type), pointer :: connections
1261) PetscInt :: num_connections
1262) PetscInt :: iconn
1263) PetscInt :: id_up, id_dn
1264) PetscReal :: v(3), v_up(3), v_dn(3)
1265) PetscReal :: distance
1266) character(len=MAXSTRINGLENGTH) :: string
1267) PetscBool :: error
1268)
1269) num_connections = size(explicit_grid%connections,2)
1270) connections => ConnectionCreate(num_connections,INTERNAL_CONNECTION_TYPE)
1271)
1272) error = PETSC_FALSE
1273) do iconn = 1, num_connections
1274) id_up = explicit_grid%connections(1,iconn)
1275) id_dn = explicit_grid%connections(2,iconn)
1276) connections%id_up(iconn) = id_up
1277) connections%id_dn(iconn) = id_dn
1278)
1279) v_up(1) = explicit_grid%face_centroids(iconn)%x - &
1280) explicit_grid%cell_centroids(id_up)%x
1281) v_up(2) = explicit_grid%face_centroids(iconn)%y - &
1282) explicit_grid%cell_centroids(id_up)%y
1283) v_up(3) = explicit_grid%face_centroids(iconn)%z - &
1284) explicit_grid%cell_centroids(id_up)%z
1285)
1286) v_dn(1) = explicit_grid%cell_centroids(id_dn)%x - &
1287) explicit_grid%face_centroids(iconn)%x
1288) v_dn(2) = explicit_grid%cell_centroids(id_dn)%y - &
1289) explicit_grid%face_centroids(iconn)%y
1290) v_dn(3) = explicit_grid%cell_centroids(id_dn)%z - &
1291) explicit_grid%face_centroids(iconn)%z
1292)
1293) v = v_up + v_dn
1294) distance = sqrt(DotProduct(v,v))
1295) if (dabs(distance) < 1.d-40) then
1296) write(string,'(2(es16.9,","),es16.9)') &
1297) explicit_grid%face_centroids(iconn)%x, &
1298) explicit_grid%face_centroids(iconn)%y, &
1299) explicit_grid%face_centroids(iconn)%z
1300) error = PETSC_TRUE
1301) option%io_buffer = 'Coincident cell and face centroids found at (' // &
1302) trim(adjustl(string)) // ') '
1303) call printMsgByRank(option)
1304) endif
1305) connections%dist(-1,iconn) = sqrt(DotProduct(v_up,v_up))/distance
1306) connections%dist(0,iconn) = distance
1307) connections%dist(1:3,iconn) = v/distance
1308) connections%area(iconn) = explicit_grid%face_areas(iconn)
1309) enddo
1310) if (error) then
1311) option%io_buffer = 'Coincident cell and face centroids found in ' // &
1312) 'UGridExplicitSetInternConnect(). See details above.'
1313) call printErrMsgByRank(option)
1314) endif
1315)
1316) UGridExplicitSetInternConnect => connections
1317)
1318) end function UGridExplicitSetInternConnect
1319)
1320) ! ************************************************************************** !
1321)
1322) subroutine UGridExplicitComputeVolumes(ugrid,option,volume)
1323) !
1324) ! Sets the volume of each grid cell
1325) !
1326) ! Author: Glenn Hammond
1327) ! Date: 05/17/12
1328) !
1329)
1330) use Option_module
1331)
1332) implicit none
1333)
1334) type(grid_unstructured_type) :: ugrid
1335) type(option_type) :: option
1336) Vec :: volume
1337)
1338) type(unstructured_explicit_type), pointer :: explicit_grid
1339)
1340) PetscInt :: icell
1341) PetscReal, pointer :: vec_ptr(:)
1342) PetscErrorCode :: ierr
1343)
1344) explicit_grid => ugrid%explicit_grid
1345)
1346) call VecGetArrayF90(volume,vec_ptr,ierr);CHKERRQ(ierr)
1347) do icell = 1, ugrid%nlmax
1348) vec_ptr(icell) = explicit_grid%cell_volumes(icell)
1349) enddo
1350) call VecRestoreArrayF90(volume,vec_ptr,ierr);CHKERRQ(ierr)
1351)
1352) end subroutine UGridExplicitComputeVolumes
1353)
1354) ! ************************************************************************** !
1355)
1356) function UGridExplicitSetBoundaryConnect(explicit_grid,cell_ids, &
1357) face_centroids,face_areas, &
1358) region_name,option)
1359) !
1360) ! Sets up the boundary connectivity within
1361) ! the connectivity object
1362) !
1363) ! Author: Glenn Hammond
1364) ! Date: 05/18/12
1365) !
1366)
1367) use Utility_module
1368) use Connection_module
1369) use Option_module
1370)
1371) implicit none
1372)
1373) type(connection_set_type), pointer :: UGridExplicitSetBoundaryConnect
1374)
1375) type(unstructured_explicit_type) :: explicit_grid
1376) PetscInt :: cell_ids(:)
1377) type(point3d_type) :: face_centroids(:)
1378) PetscReal :: face_areas(:)
1379) character(len=MAXWORDLENGTH) :: region_name
1380) type(option_type) :: option
1381)
1382) type(connection_set_type), pointer :: connections
1383) PetscInt :: num_connections
1384) PetscInt :: iconn
1385) PetscInt :: id
1386) PetscReal :: v(3)
1387) PetscReal :: distance
1388) character(len=MAXSTRINGLENGTH) :: string
1389) PetscBool :: error
1390)
1391) num_connections = size(cell_ids)
1392) connections => ConnectionCreate(num_connections,BOUNDARY_CONNECTION_TYPE)
1393)
1394) error = PETSC_FALSE
1395) do iconn = 1, num_connections
1396) id = cell_ids(iconn)
1397) connections%id_dn(iconn) = id
1398)
1399) v(1) = explicit_grid%cell_centroids(id)%x - &
1400) face_centroids(iconn)%x
1401) v(2) = explicit_grid%cell_centroids(id)%y - &
1402) face_centroids(iconn)%y
1403) v(3) = explicit_grid%cell_centroids(id)%z - &
1404) face_centroids(iconn)%z
1405)
1406) distance = sqrt(DotProduct(v,v))
1407) if (dabs(distance) < 1.d-40) then
1408) write(string,'(2(es16.9,","),es16.9)') &
1409) face_centroids(iconn)%x, face_centroids(iconn)%y, &
1410) face_centroids(iconn)%z
1411) error = PETSC_TRUE
1412) option%io_buffer = 'Coincident cell and face centroids found at (' // &
1413) trim(adjustl(string)) // ') '
1414) call printMsgByRank(option)
1415) endif
1416) connections%dist(-1,iconn) = 0.d0
1417) connections%dist(0,iconn) = distance
1418) connections%dist(1:3,iconn) = v/distance
1419) connections%area(iconn) = face_areas(iconn)
1420) enddo
1421) if (error) then
1422) option%io_buffer = 'Coincident cell and face centroids found in ' // &
1423) 'UGridExplicitSetBoundaryConnect() for region "' // trim(region_name) // &
1424) '". See details above.'
1425) call printErrMsgByRank(option)
1426) endif
1427)
1428) UGridExplicitSetBoundaryConnect => connections
1429)
1430) end function UGridExplicitSetBoundaryConnect
1431)
1432) ! ************************************************************************** !
1433)
1434) function UGridExplicitSetConnections(explicit_grid,cell_ids,connection_type, &
1435) option)
1436) !
1437) ! Sets up the connectivity for a region
1438) !
1439) ! Author: Glenn Hammond
1440) ! Date: 05/18/12
1441) !
1442)
1443) use Utility_module
1444) use Connection_module
1445) use Option_module
1446)
1447) implicit none
1448)
1449) type(connection_set_type), pointer :: UGridExplicitSetConnections
1450)
1451) type(unstructured_explicit_type) :: explicit_grid
1452) PetscInt, pointer :: cell_ids(:)
1453) PetscInt :: connection_type
1454) type(option_type) :: option
1455)
1456) type(connection_set_type), pointer :: connections
1457) PetscInt :: num_connections
1458) PetscInt :: iconn
1459) PetscInt :: id
1460)
1461) num_connections = 0
1462) if (associated(cell_ids)) then
1463) num_connections = size(cell_ids)
1464) endif
1465) connections => ConnectionCreate(num_connections,connection_type)
1466)
1467) do iconn = 1, num_connections
1468) id = cell_ids(iconn)
1469) connections%id_dn(iconn) = id
1470) enddo
1471)
1472) UGridExplicitSetConnections => connections
1473)
1474) end function UGridExplicitSetConnections
1475)
1476) end module Grid_Unstructured_Explicit_module