geomech_grid.F90 coverage: 100.00 %func 80.61 %block
1) module Geomechanics_Grid_module
2)
3) use Geomechanics_Grid_Aux_module
4) use Grid_Unstructured_Cell_module
5) use PFLOTRAN_Constants_module
6)
7) implicit none
8)
9) private
10)
11) #include "petsc/finclude/petscsys.h"
12) #include "petsc/finclude/petscvec.h"
13) #include "petsc/finclude/petscvec.h90"
14) #include "petsc/finclude/petscis.h"
15) #include "petsc/finclude/petscis.h90"
16) #if defined(SCORPIO)
17) include "scorpiof.h"
18) #endif
19)
20) ! PetscInt, parameter :: HEX_TYPE = 1
21) ! PetscInt, parameter :: TET_TYPE = 2
22) ! PetscInt, parameter :: WEDGE_TYPE = 3
23) ! PetscInt, parameter :: PYR_TYPE = 4
24) ! PetscInt, parameter :: TRI_FACE_TYPE = 1
25) ! PetscInt, parameter :: QUAD_FACE_TYPE = 2
26) ! PetscInt, parameter :: MAX_VERT_PER_FACE = 4
27)
28) public :: CopySubsurfaceGridtoGeomechGrid, &
29) GeomechGridLocalizeRegions, &
30) GeomechGridVecGetArrayF90, &
31) GeomechGridVecRestoreArrayF90, &
32) GeomechGridCopyIntegerArrayToVec, &
33) GeomechGridCopyVecToIntegerArray, &
34) GeomechSubsurfMapFromFilename
35)
36) contains
37)
38) ! ************************************************************************** !
39) !
40) ! CopySubsurfaceGridtoGeomechGrid: Subroutine to copy subsurface grid info.
41) ! to geomechanics grid
42) ! author: Satish Karra, LANL
43) ! date: 05/30/13
44) !
45) ! ************************************************************************** !
46) subroutine CopySubsurfaceGridtoGeomechGrid(ugrid,geomech_grid,option)
47)
48) use Grid_Unstructured_Aux_module
49) use Geomechanics_Grid_Aux_module
50) use Option_module
51) use Gauss_module
52)
53) implicit none
54)
55) #include "petsc/finclude/petscvec.h"
56) #include "petsc/finclude/petscvec.h90"
57) #include "petsc/finclude/petscmat.h"
58) #include "petsc/finclude/petscmat.h90"
59) #include "petsc/finclude/petscdm.h"
60) #include "petsc/finclude/petscdm.h90"
61) #include "petsc/finclude/petscis.h"
62) #include "petsc/finclude/petscis.h90"
63) #include "petsc/finclude/petscviewer.h"
64)
65) type(grid_unstructured_type), pointer :: ugrid
66) type(geomech_grid_type), pointer :: geomech_grid
67) type(option_type), pointer :: option
68) PetscInt :: local_id
69) PetscInt :: ghosted_id
70) PetscInt :: vertex_count
71) PetscInt :: ivertex
72) PetscInt :: vertex_id
73) PetscInt :: count
74) PetscInt, allocatable :: int_array(:)
75) PetscInt, allocatable :: int_array2(:)
76) PetscInt, allocatable :: int_array3(:)
77) PetscInt, allocatable :: int_array4(:)
78) PetscErrorCode :: ierr
79) character(len=MAXSTRINGLENGTH) :: string, string1
80) PetscInt :: global_offset_old
81) PetscInt :: global_offset
82) Mat :: Rank_Mat
83) PetscReal :: rank
84) PetscViewer :: viewer
85) PetscReal, pointer :: vec_ptr(:)
86) PetscInt :: istart,iend
87) PetscBool :: vertex_found
88) PetscInt :: int_rank
89) PetscInt :: vertex_count2
90) IS :: is_rank
91) IS :: is_rank_new
92) IS :: is_natural
93) IS :: is_ghost_petsc
94) PetscReal :: max_val
95) PetscInt :: row
96) PetscScalar, allocatable :: val(:)
97) PetscInt :: ncols
98) PetscInt, allocatable :: cols(:)
99) AO :: ao_natural_to_petsc_nodes
100) PetscInt :: nlmax_node
101) PetscInt, pointer :: int_ptr(:)
102) type(point_type), pointer :: vertices(:)
103) PetscInt, allocatable :: vertex_count_array(:)
104) PetscInt, allocatable :: vertex_count_array2(:)
105)
106)
107) #ifdef GEOMECH_DEBUG
108) call printMsg(option,'Copying unstructured grid to geomechanics grid')
109) #endif
110)
111) geomech_grid%global_offset_elem = ugrid%global_offset
112) geomech_grid%nmax_elem = ugrid%nmax
113) geomech_grid%nlmax_elem = ugrid%nlmax
114) geomech_grid%nmax_node = ugrid%num_vertices_global
115)
116) ! Element natural ids
117) allocate(geomech_grid%elem_ids_natural(geomech_grid%nlmax_elem))
118) do local_id = 1, geomech_grid%nlmax_elem
119) geomech_grid%elem_ids_natural(local_id) = ugrid%cell_ids_natural(local_id)
120) enddo
121)
122) allocate(geomech_grid%elem_ids_petsc(size(ugrid%cell_ids_petsc)))
123) geomech_grid%elem_ids_petsc = ugrid%cell_ids_petsc
124) geomech_grid%ao_natural_to_petsc = ugrid%ao_natural_to_petsc
125) geomech_grid%max_ndual_per_elem = ugrid%max_ndual_per_cell
126) geomech_grid%max_nnode_per_elem = ugrid%max_nvert_per_cell
127) geomech_grid%max_elem_sharing_a_node = ugrid%max_cells_sharing_a_vertex
128) geomech_grid%nlmax_node = ugrid%num_vertices_natural
129)
130) #ifdef GEOMECH_DEBUG
131) call printMsg(option,'Removing ghosted elements (cells)')
132) #endif
133)
134) ! Type of element
135) allocate(geomech_grid%elem_type(geomech_grid%nlmax_elem))
136) do local_id = 1, geomech_grid%nlmax_elem
137) geomech_grid%elem_type(local_id) = ugrid%cell_type(local_id)
138) enddo
139)
140) #ifdef GEOMECH_DEBUG
141) call printMsg(option,'Reordering nodes (vertices)')
142) #endif
143)
144) ! First calculate number of elements on local domain (without ghosted elements)
145) vertex_count = 0
146) do local_id = 1, geomech_grid%nlmax_elem
147) vertex_count = vertex_count + ugrid%cell_vertices(0,local_id)
148) enddo
149)
150) ! Store all the vertices in int_array
151) count = 0
152) allocate(int_array(vertex_count))
153) do local_id = 1, geomech_grid%nlmax_elem
154) do ivertex = 1, ugrid%cell_vertices(0,local_id)
155) count = count + 1
156) int_array(count) = ugrid%cell_vertices(ivertex,local_id)
157) enddo
158) enddo
159)
160) ! Sort the vertex ids
161) allocate(int_array2(vertex_count))
162) do ivertex = 1, vertex_count
163) int_array2(ivertex) = ivertex
164) enddo
165) int_array2 = int_array2 - 1
166) call PetscSortIntWithPermutation(vertex_count,int_array,int_array2, &
167) ierr);CHKERRQ(ierr)
168) int_array2 = int_array2+1
169)
170) ! Remove duplicates
171) allocate(int_array3(vertex_count))
172) allocate(int_array4(vertex_count))
173) int_array3 = 0
174) int_array4 = 0
175) int_array3(1) = int_array(int_array2(1))
176) count = 1
177) int_array4(int_array2(1)) = count
178) do ivertex = 2, vertex_count
179) vertex_id = int_array(int_array2(ivertex))
180) if (vertex_id > int_array3(count)) then
181) count = count + 1
182) int_array3(count) = vertex_id
183) endif
184) int_array4(int_array2(ivertex)) = count
185) enddo
186) vertex_count = count
187) deallocate(int_array)
188)
189) ! Store the vertices including ghosted ones in natural order
190) allocate(geomech_grid%node_ids_ghosted_natural(vertex_count))
191) do ivertex = 1, vertex_count
192) geomech_grid%node_ids_ghosted_natural(ivertex) = ugrid% &
193) vertex_ids_natural(int_array3(ivertex))
194) enddo
195)
196) #ifdef GEOMECH_DEBUG
197) write(string,*) option%myrank
198) string = 'geomech_node_ids_ghosted_natural' &
199) // trim(adjustl(string)) // '.out'
200) open(unit=86,file=trim(string))
201) do local_id = 1, vertex_count
202) write(86,'(i5)') geomech_grid%node_ids_ghosted_natural(local_id)
203) enddo
204) close(86)
205) #endif
206)
207) allocate(geomech_grid%elem_nodes( &
208) 0:geomech_grid%max_nnode_per_elem,geomech_grid%nlmax_elem))
209) geomech_grid%elem_nodes = 0
210)
211)
212) ! Store the vertices (natural ordering) of each element
213) count = 0
214) do local_id = 1, geomech_grid%nlmax_elem
215) geomech_grid%elem_nodes(0,local_id) = ugrid%cell_vertices(0,local_id)
216) do ivertex = 1, geomech_grid%elem_nodes(0,local_id)
217) count = count + 1
218) geomech_grid%elem_nodes(ivertex,local_id) = int_array4(count)
219) enddo
220) enddo
221)
222) #ifdef GEOMECH_DEBUG
223) write(string,*) option%myrank
224) string = 'geomech_elem_nodes' // trim(adjustl(string)) // '.out'
225) open(unit=86,file=trim(string))
226) do local_id = 1, geomech_grid%nlmax_elem
227) write(86,'(i5)') geomech_grid%elem_nodes(0,local_id)
228) do ivertex = 1, geomech_grid%max_nnode_per_elem
229) write(86,'(i5)') geomech_grid%elem_nodes(ivertex,local_id)
230) enddo
231) enddo
232) close(86)
233) #endif
234)
235) deallocate(int_array2)
236) deallocate(int_array4)
237)
238) ! Store the coordinates of the vertices on each process
239) allocate(geomech_grid%nodes(vertex_count))
240) do ivertex = 1, vertex_count
241) geomech_grid%nodes(ivertex)%id = &
242) geomech_grid%node_ids_ghosted_natural(ivertex)
243) geomech_grid%nodes(ivertex)%x = ugrid%vertices(int_array3(ivertex))%x
244) geomech_grid%nodes(ivertex)%y = ugrid%vertices(int_array3(ivertex))%y
245) geomech_grid%nodes(ivertex)%z = ugrid%vertices(int_array3(ivertex))%z
246) enddo
247)
248) #ifdef GEOMECH_DEBUG
249) write(string,*) option%myrank
250) string = 'geomech_node_coordinates' // trim(adjustl(string)) // '.out'
251) open(unit=86,file=trim(string))
252) do ivertex = 1, vertex_count
253) write(86,'(i5)') geomech_grid%nodes(ivertex)%id
254) write(86,'(1pe12.4)') geomech_grid%nodes(ivertex)%x
255) write(86,'(1pe12.4)') geomech_grid%nodes(ivertex)%y
256) write(86,'(1pe12.4)') geomech_grid%nodes(ivertex)%z
257) enddo
258) close(86)
259) #endif
260)
261) geomech_grid%ngmax_node = vertex_count
262)
263) deallocate(int_array3)
264)
265) ! So far we have stored the ghosted vertices on each process
266) ! Now we will assign the local vertices on each process
267) ! For this, we first calculate the maximum rank of all the processes
268) ! that share a vertex. This rank will have the vertex as its local
269) ! vertex.
270)
271)
272) ! Create a nmax_node X mycommsize matrix
273) ! For each row numbered by vertex (-1, zero-base),
274) ! the ranks of the processes that possess the vertex are stored
275) call MatCreateAIJ(option%mycomm,PETSC_DECIDE,ONE_INTEGER, &
276) geomech_grid%nmax_node,option%mycommsize, &
277) option%mycommsize,PETSC_NULL_INTEGER, &
278) option%mycommsize,PETSC_NULL_INTEGER,Rank_Mat, &
279) ierr);CHKERRQ(ierr)
280)
281) call MatZeroEntries(Rank_Mat,ierr);CHKERRQ(ierr)
282)
283) rank = option%myrank + 1
284) do ivertex = 1, geomech_grid%ngmax_node
285) call MatSetValue(Rank_Mat, &
286) geomech_grid%node_ids_ghosted_natural(ivertex)-1, &
287) option%myrank,rank,INSERT_VALUES,ierr);CHKERRQ(ierr)
288) enddo
289) call MatAssemblyBegin(Rank_Mat,MAT_FINAL_ASSEMBLY,ierr);CHKERRQ(ierr)
290) call MatAssemblyEnd(Rank_Mat,MAT_FINAL_ASSEMBLY,ierr);CHKERRQ(ierr)
291)
292) #ifdef GEOMECH_DEBUG
293) call PetscViewerASCIIOpen(option%mycomm,'Rank_Mat.out', &
294) viewer,ierr);CHKERRQ(ierr)
295) call MatView(Rank_Mat,viewer,ierr);CHKERRQ(ierr)
296) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
297) #endif
298)
299) ! Now find the maximum of all the ranks for each vertex
300) allocate(val(option%mycommsize))
301) allocate(cols(option%mycommsize))
302) call MatGetOwnershipRange(Rank_Mat,istart,iend,ierr);CHKERRQ(ierr)
303) allocate(int_array(iend-istart))
304) count = 0
305) do row = istart, iend-1
306) call MatGetRow(Rank_Mat,row,ncols,cols,val,ierr);CHKERRQ(ierr)
307) max_val = 0.d0
308) do local_id = 1, ncols
309) max_val = max(max_val,val(local_id))
310) enddo
311) count = count + 1
312) int_array(count) = int(max_val)
313) call MatRestoreRow(Rank_Mat,row,ncols,cols,val,ierr);CHKERRQ(ierr)
314) enddo
315) deallocate(val)
316) deallocate(cols)
317) call MatDestroy(Rank_Mat,ierr);CHKERRQ(ierr)
318)
319) ! Change rank to start from 0
320) int_array = int_array - 1
321)
322) #ifdef GEOMECH_DEBUG
323) write(string,*) option%myrank
324) string = 'geomech_node_ranks' // trim(adjustl(string)) // '.out'
325) open(unit=86,file=trim(string))
326) do row = 1, count
327) write(86,'(i5)') int_array(row)
328) enddo
329) close(86)
330) #endif
331)
332) ! Now count the number of vertices that are local to each rank
333) allocate(vertex_count_array(option%mycommsize))
334) allocate(vertex_count_array2(option%mycommsize))
335) vertex_count_array = 0
336) vertex_count_array2 = 0
337)
338) do int_rank = 0, option%mycommsize
339) do local_id = 1, count
340) if (int_array(local_id) == int_rank) then
341) vertex_count_array(int_rank+1) = vertex_count_array(int_rank+1) + 1
342) endif
343) enddo
344) enddo
345) call MPI_Allreduce(vertex_count_array,vertex_count_array2, &
346) option%mycommsize,MPIU_INTEGER,MPI_SUM, &
347) option%mycomm,ierr)
348)
349) do int_rank = 0, option%mycommsize
350) if (option%myrank == int_rank) geomech_grid%nlmax_node = &
351) vertex_count_array2(int_rank+1)
352) if (geomech_grid%nlmax_node > geomech_grid%ngmax_node) then
353) option%io_buffer = 'Error: nlmax_node cannot be greater than' // &
354) ' ngmax_node.'
355) call printErrMsg(option)
356) endif
357) enddo
358)
359)
360) if (allocated(vertex_count_array)) deallocate(vertex_count_array)
361) if (allocated(vertex_count_array2)) deallocate(vertex_count_array2)
362)
363) ! Add a check on nlmax_node to see if there are too many processes
364) call MPI_Allreduce(geomech_grid%nlmax_node,nlmax_node, &
365) ONE_INTEGER_MPI,MPIU_INTEGER,MPI_MIN, &
366) option%mycomm,ierr)
367)
368) if (nlmax_node < 1) then
369) option%io_buffer = 'Error: Too many processes for the size of the domain.'
370) call printErrMsg(option)
371) endif
372)
373) ! Create an index set of the ranks of each vertex
374) call ISCreateGeneral(option%mycomm,count,int_array,PETSC_COPY_VALUES, &
375) is_rank,ierr);CHKERRQ(ierr)
376)
377) #if GEOMECH_DEBUG
378) call PetscViewerASCIIOpen(option%mycomm,'geomech_is_rank_nodes.out', &
379) viewer,ierr);CHKERRQ(ierr)
380) call ISView(is_rank,viewer,ierr);CHKERRQ(ierr)
381) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
382) #endif
383) deallocate(int_array)
384)
385) allocate(int_array(count))
386) global_offset_old = 0
387) call MPI_Exscan(count,global_offset_old, &
388) ONE_INTEGER_MPI,MPIU_INTEGER,MPI_SUM,option%mycomm,ierr)
389) do local_id = 1, count
390) int_array(local_id) = (local_id-1) + global_offset_old
391) enddo
392) call ISCreateGeneral(option%mycomm,count,int_array,PETSC_COPY_VALUES, &
393) is_natural,ierr);CHKERRQ(ierr)
394) deallocate(int_array)
395)
396) #if GEOMECH_DEBUG
397) call PetscViewerASCIIOpen(option%mycomm,'geomech_is_natural.out', &
398) viewer,ierr);CHKERRQ(ierr)
399) call ISView(is_natural,viewer,ierr);CHKERRQ(ierr)
400) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
401) #endif
402)
403) ! Find the global_offset for vertices on this rank
404) global_offset_old = 0
405) call MPI_Exscan(geomech_grid%nlmax_node,global_offset_old, &
406) ONE_INTEGER_MPI,MPIU_INTEGER,MPI_SUM,option%mycomm,ierr)
407)
408) geomech_grid%global_offset = global_offset_old
409)
410) #ifdef GEOMECH_DEBUG
411) write(string,*) option%myrank
412) print *, 'Number of local vertices on rank ' // &
413) trim(adjustl(string)) // ' is:', geomech_grid%nlmax_node
414) print *, 'Global offset of vertices on rank ' // &
415) trim(adjustl(string)) // ' is:', global_offset_old
416) print *, 'Number of ghosted vertices on rank ' // &
417) trim(adjustl(string)) // ' is:', geomech_grid%ngmax_node
418) #endif
419)
420) call ISPartitioningToNumbering(is_rank,is_rank_new,ierr);CHKERRQ(ierr)
421) call ISDestroy(is_rank,ierr);CHKERRQ(ierr)
422)
423) #if GEOMECH_DEBUG
424) call PetscViewerASCIIOpen(option%mycomm,'geomech_is_rank_nodes_new.out', &
425) viewer,ierr);CHKERRQ(ierr)
426) call ISView(is_rank_new,viewer,ierr);CHKERRQ(ierr)
427) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
428) #endif
429)
430) !Create an application ordering
431) call AOCreateBasicIS(is_natural,is_rank_new,ao_natural_to_petsc_nodes, &
432) ierr);CHKERRQ(ierr)
433) call ISDestroy(is_rank_new,ierr);CHKERRQ(ierr)
434) call ISDestroy(is_natural,ierr);CHKERRQ(ierr)
435)
436) #if GEOMECH_DEBUG
437) call PetscViewerASCIIOpen(option%mycomm, &
438) 'geomech_ao_natural_to_petsc_nodes.out', &
439) viewer,ierr);CHKERRQ(ierr)
440) call AOView(ao_natural_to_petsc_nodes,viewer,ierr);CHKERRQ(ierr)
441) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
442) #endif
443)
444) ! Create a new IS with local PETSc numbering
445) allocate(int_array(geomech_grid%nlmax_node))
446) do local_id = 1, geomech_grid%nlmax_node
447) int_array(local_id) = (local_id-1) + geomech_grid%global_offset
448) enddo
449) call ISCreateGeneral(option%mycomm,geomech_grid%nlmax_node, &
450) int_array,PETSC_COPY_VALUES,is_natural, &
451) ierr);CHKERRQ(ierr)
452) deallocate(int_array)
453)
454) #if GEOMECH_DEBUG
455) call PetscViewerASCIIOpen(option%mycomm,'geomech_is_petsc.out', &
456) viewer,ierr);CHKERRQ(ierr)
457) call ISView(is_natural,viewer,ierr);CHKERRQ(ierr)
458) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
459) #endif
460)
461) ! Rewrite the IS with natural numbering
462) call AOPetscToApplicationIS(ao_natural_to_petsc_nodes,is_natural, &
463) ierr);CHKERRQ(ierr)
464)
465) ! These are the natural ids of the vertices local to each process
466) #if GEOMECH_DEBUG
467) call PetscViewerASCIIOpen(option%mycomm, &
468) 'geomech_is_natural_after_ordering.out', &
469) viewer,ierr);CHKERRQ(ierr)
470) call ISView(is_natural,viewer,ierr);CHKERRQ(ierr)
471) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
472) #endif
473)
474) ! Get the local indices (natural) and store them
475) allocate(geomech_grid%node_ids_local_natural(geomech_grid%nlmax_node))
476) call ISGetIndicesF90(is_natural,int_ptr,ierr);CHKERRQ(ierr)
477) do local_id = 1, geomech_grid%nlmax_node
478) geomech_grid%node_ids_local_natural(local_id) = int_ptr(local_id)
479) enddo
480) call ISRestoreIndicesF90(is_natural,int_ptr,ierr);CHKERRQ(ierr)
481) call ISDestroy(is_natural,ierr);CHKERRQ(ierr)
482)
483) ! Changing to 1-based
484) geomech_grid%node_ids_local_natural = geomech_grid%node_ids_local_natural + 1
485)
486) ! Find the natural ids of ghost nodes (vertices)
487) vertex_count = 0
488) if (geomech_grid%ngmax_node - geomech_grid%nlmax_node > 0) then
489) allocate(int_array2(geomech_grid%ngmax_node-geomech_grid%nlmax_node))
490) do ivertex = 1, geomech_grid%ngmax_node
491) do local_id = 1, geomech_grid%nlmax_node
492) vertex_found = PETSC_FALSE
493) if (geomech_grid%node_ids_ghosted_natural(ivertex) == &
494) geomech_grid%node_ids_local_natural(local_id)) then
495) vertex_found = PETSC_TRUE
496) exit
497) endif
498) enddo
499) if (.not.vertex_found) then
500) vertex_count = vertex_count + 1
501) int_array2(vertex_count) = &
502) geomech_grid%node_ids_ghosted_natural(ivertex)
503) endif
504) enddo
505) else
506) allocate(int_array2(1))
507) int_array2 = 0
508) endif
509)
510) if (vertex_count /= geomech_grid%ngmax_node - geomech_grid%nlmax_node) then
511) option%io_buffer = 'Error in number of ghost nodes!'
512) call printErrMsg(option)
513) endif
514)
515) #ifdef GEOMECH_DEBUG
516) write(string,*) option%myrank
517) string = 'geomech_node_ids_ghosts_natural' // trim(adjustl(string)) // '.out'
518) open(unit=86,file=trim(string))
519) if (vertex_count > 0) then
520) do local_id = 1, vertex_count
521) write(86,'(i5)') int_array2(local_id)
522) enddo
523) else
524) write(86,*) 'There are no ghost nodes (vertices) on this process.'
525) endif
526) close(86)
527) #endif
528)
529) geomech_grid%num_ghost_nodes = geomech_grid%ngmax_node - &
530) geomech_grid%nlmax_node
531)
532) ! Changing the index to 0 based
533) if (allocated(int_array2)) &
534) int_array2 = int_array2 - 1
535)
536) ! Create a new IS with local PETSc numbering
537) call ISCreateGeneral(option%mycomm,geomech_grid%num_ghost_nodes, &
538) int_array2,PETSC_COPY_VALUES,is_ghost_petsc, &
539) ierr);CHKERRQ(ierr)
540)
541) ! Natural ids of ghost vertices on each rank
542) #if GEOMECH_DEBUG
543) call PetscViewerASCIIOpen(option%mycomm,'geomech_is_ghost_natural.out', &
544) viewer,ierr);CHKERRQ(ierr)
545) call ISView(is_ghost_petsc,viewer,ierr);CHKERRQ(ierr)
546) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
547) #endif
548)
549) ! Rewrite the IS with natural numbering
550) call AOApplicationtoPetscIS(ao_natural_to_petsc_nodes,is_ghost_petsc, &
551) ierr);CHKERRQ(ierr)
552)
553) ! Petsc ids of ghost vertices on each rank
554) #if GEOMECH_DEBUG
555) call PetscViewerASCIIOpen(option%mycomm,'geomech_is_ghost_petsc.out', &
556) viewer,ierr);CHKERRQ(ierr)
557) call ISView(is_ghost_petsc,viewer,ierr);CHKERRQ(ierr)
558) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
559) #endif
560)
561) if (allocated(int_array2)) then
562) if (vertex_count > 0) then
563) allocate(geomech_grid%ghosted_node_ids_natural(vertex_count))
564) ! Change back to 1-based
565) geomech_grid%ghosted_node_ids_natural = int_array2 + 1
566) endif
567) deallocate(int_array2)
568) endif
569)
570) ! Store the petsc indices for ghost nodes
571) if (geomech_grid%num_ghost_nodes > 0) then
572) allocate(geomech_grid%ghosted_node_ids_petsc(geomech_grid%num_ghost_nodes))
573) call ISGetIndicesF90(is_ghost_petsc,int_ptr,ierr);CHKERRQ(ierr)
574) do ghosted_id = 1, geomech_grid%num_ghost_nodes
575) geomech_grid%ghosted_node_ids_petsc(ghosted_id) = int_ptr(ghosted_id)
576) enddo
577) call ISRestoreIndicesF90(is_ghost_petsc,int_ptr,ierr);CHKERRQ(ierr)
578) endif
579) call ISDestroy(is_ghost_petsc,ierr);CHKERRQ(ierr)
580)
581)
582) ! Changing back to 1-based
583) if (geomech_grid%num_ghost_nodes > 0) &
584) geomech_grid%ghosted_node_ids_petsc = geomech_grid%ghosted_node_ids_petsc + 1
585)
586) geomech_grid%ao_natural_to_petsc_nodes = ao_natural_to_petsc_nodes
587)
588) ! The following is for re-ordering of local ghosted numbering such that
589) ! the first nlmax_node values are local nodes and the rest are ghost nodes
590) allocate(int_array(geomech_grid%ngmax_node))
591) allocate(int_array2(geomech_grid%ngmax_node))
592) do ivertex = 1, geomech_grid%nlmax_node
593) int_array(ivertex) = geomech_grid%node_ids_local_natural(ivertex)
594) enddo
595) if (geomech_grid%num_ghost_nodes > 0) then
596) do ivertex = geomech_grid%nlmax_node+1, geomech_grid%ngmax_node
597) int_array(ivertex) = geomech_grid% &
598) ghosted_node_ids_natural(ivertex-geomech_grid%nlmax_node)
599) enddo
600) endif
601)
602) do ivertex = 1, geomech_grid%ngmax_node
603) int_array2(ivertex) = ivertex
604) enddo
605) int_array2 = int_array2 - 1
606) call PetscSortIntWithPermutation(geomech_grid%ngmax_node,int_array, &
607) int_array2,ierr);CHKERRQ(ierr)
608) int_array2 = int_array2+1
609)
610) ! Here the local ghosted ids in the elements need to be changed to reflect
611) ! the change in ordering to first nlmax_node being local nodes and
612) ! the rest being ghost nodes
613) do local_id = 1, geomech_grid%nlmax_elem
614) do ivertex = 1, geomech_grid%elem_nodes(0,local_id)
615) geomech_grid%elem_nodes(ivertex,local_id) = &
616) int_array2(geomech_grid%elem_nodes(ivertex,local_id))
617) enddo
618) enddo
619)
620) #ifdef GEOMECH_DEBUG
621) write(string,*) option%myrank
622) string = 'geomech_elem_nodes_reordered' // trim(adjustl(string)) // '.out'
623) open(unit=86,file=trim(string))
624) do local_id = 1, geomech_grid%nlmax_elem
625) write(86,'(i5)') geomech_grid%elem_nodes(0,local_id)
626) do ivertex = 1, geomech_grid%elem_nodes(0,local_id)
627) write(86,'(i5)') geomech_grid%elem_nodes(ivertex,local_id)
628) enddo
629) enddo
630) close(86)
631) #endif
632)
633) ! Now re-order the geomech_grid%nodes datastructure due to the re-ordering of
634) ! the vertices. Temporarily store in vertices
635) allocate(vertices(geomech_grid%ngmax_node))
636) do ghosted_id = 1, geomech_grid%ngmax_node
637) vertices(ghosted_id)%id = geomech_grid%nodes(ghosted_id)%id
638) vertices(ghosted_id)%x = geomech_grid%nodes(ghosted_id)%x
639) vertices(ghosted_id)%y = geomech_grid%nodes(ghosted_id)%y
640) vertices(ghosted_id)%z = geomech_grid%nodes(ghosted_id)%z
641) enddo
642)
643) do ghosted_id = 1, geomech_grid%ngmax_node
644) geomech_grid%nodes(int_array2(ghosted_id))%id = vertices(ghosted_id)%id
645) geomech_grid%nodes(int_array2(ghosted_id))%x = vertices(ghosted_id)%x
646) geomech_grid%nodes(int_array2(ghosted_id))%y = vertices(ghosted_id)%y
647) geomech_grid%nodes(int_array2(ghosted_id))%z = vertices(ghosted_id)%z
648) enddo
649)
650) #ifdef GEOMECH_DEBUG
651) write(string,*) option%myrank
652) string = 'geomech_node_coordinates_reordered' &
653) // trim(adjustl(string)) // '.out'
654) open(unit=86,file=trim(string))
655) do ivertex = 1, vertex_count
656) write(86,'(i5)') geomech_grid%nodes(ivertex)%id
657) write(86,'(1pe12.4)') geomech_grid%nodes(ivertex)%x
658) write(86,'(1pe12.4)') geomech_grid%nodes(ivertex)%y
659) write(86,'(1pe12.4)') geomech_grid%nodes(ivertex)%z
660) enddo
661) close(86)
662) #endif
663)
664) deallocate(int_array)
665) deallocate(vertices)
666)
667) allocate(int_array(geomech_grid%ngmax_node))
668)
669) int_array = geomech_grid%node_ids_ghosted_natural
670)
671) ! Store the natural ids of all the local and ghost nodes in the new
672) ! re-ordered system
673) do ghosted_id = 1, geomech_grid%ngmax_node
674) geomech_grid%node_ids_ghosted_natural(int_array2(ghosted_id)) = &
675) int_array(ghosted_id)
676) enddo
677)
678) #ifdef GEOMECH_DEBUG
679) write(string,*) option%myrank
680) string = 'geomech_node_ids_ghosted_natural_reordered' &
681) // trim(adjustl(string)) // '.out'
682) open(unit=86,file=trim(string))
683) do ghosted_id = 1, geomech_grid%ngmax_node
684) write(86,'(i5)') geomech_grid%node_ids_ghosted_natural(ghosted_id)
685) enddo
686) close(86)
687) #endif
688)
689) deallocate(int_array)
690) deallocate(int_array2)
691)
692) ! Initialize the Gauss data structure in each element
693) allocate(geomech_grid%gauss_node(geomech_grid%nlmax_elem))
694) do local_id = 1, geomech_grid%nlmax_elem
695) call GaussInitialize(geomech_grid%gauss_node(local_id))
696) geomech_grid%gauss_node(local_id)%Eletype = &
697) geomech_grid%Elem_type(local_id)
698) ! Set to 3D although we have gauss point calculations for 2D
699) geomech_grid%gauss_node(local_id)%dim = THREE_DIM_GRID
700) ! Three gauss points in each direction
701) geomech_grid%gauss_node(local_id)%NGPTS = THREE_INTEGER
702) if (geomech_grid%gauss_node(local_id)%Eletype == PYR_TYPE) &
703) geomech_grid%gauss_node(local_id)%NGPTS = FIVE_INTEGER
704) call GaussCalculatePoints(geomech_grid%gauss_node(local_id))
705) enddo
706)
707) ! Store petsc ids of the local and ghost nodes in the new re-ordered system on
708) ! each rank
709) allocate(int_array(geomech_grid%ngmax_node))
710) do local_id = 1, geomech_grid%nlmax_node
711) int_array(local_id) = local_id + geomech_grid%global_offset
712) enddo
713) do ghosted_id = geomech_grid%nlmax_node+1, geomech_grid%ngmax_node
714) int_array(ghosted_id) = &
715) geomech_grid%ghosted_node_ids_petsc(ghosted_id - geomech_grid%nlmax_node)
716) enddo
717) allocate(geomech_grid%node_ids_ghosted_petsc(geomech_grid%ngmax_node))
718) geomech_grid%node_ids_ghosted_petsc = int_array
719) deallocate(int_array)
720)
721) #ifdef GEOMECH_DEBUG
722) write(string,*) option%myrank
723) string = 'geomech_node_ids_ghosted_petsc' // trim(adjustl(string)) // '.out'
724) open(unit=86,file=trim(string))
725) do ghosted_id = 1, geomech_grid%ngmax_node
726) write(86,'(i5)') geomech_grid%node_ids_ghosted_petsc(ghosted_id)
727) enddo
728) close(86)
729) #endif
730)
731) ! Vector that stores for each vertex the number of elements it is shared by
732) ! locally on a rank
733) ! local vector
734) call VecCreate(PETSC_COMM_SELF,geomech_grid%no_elems_sharing_node_loc, &
735) ierr);CHKERRQ(ierr)
736) call VecSetSizes(geomech_grid%no_elems_sharing_node_loc, &
737) geomech_grid%ngmax_node, &
738) PETSC_DECIDE,ierr);CHKERRQ(ierr)
739) call VecSetFromOptions(geomech_grid%no_elems_sharing_node_loc, &
740) ierr);CHKERRQ(ierr)
741)
742) ! Vector that stores for each vertex the number of elements it is shared by
743) ! globally across all ranks
744) ! global vector
745) call VecCreate(option%mycomm,geomech_grid%no_elems_sharing_node, &
746) ierr);CHKERRQ(ierr)
747) call VecSetSizes(geomech_grid%no_elems_sharing_node, &
748) geomech_grid%nlmax_node, &
749) PETSC_DECIDE,ierr);CHKERRQ(ierr)
750) call VecSetFromOptions(geomech_grid%no_elems_sharing_node, &
751) ierr);CHKERRQ(ierr)
752)
753) call VecSet(geomech_grid%no_elems_sharing_node_loc,0.d0,ierr);CHKERRQ(ierr)
754) call VecSet(geomech_grid%no_elems_sharing_node,0.d0,ierr);CHKERRQ(ierr)
755)
756) end subroutine CopySubsurfaceGridtoGeomechGrid
757)
758) ! ************************************************************************** !
759) !
760) ! GeomechGridLocalizeRegions: Resticts regions to vertices local
761) ! to processor for geomech grid when
762) ! the region is defined by a list of
763) ! vertex ids
764) ! author: Satish Karra
765) ! date: 06/13/13
766) !
767) ! ************************************************************************** !
768) subroutine GeomechGridLocalizeRegions(grid,region_list,option)
769)
770) use Option_module
771) use Geomechanics_Region_module
772)
773) implicit none
774)
775) type(gm_region_list_type), pointer :: region_list
776) type(geomech_grid_type), pointer :: grid
777) type(option_type) :: option
778)
779) type(gm_region_type), pointer :: region
780) character(len=MAXSTRINGLENGTH) :: string
781)
782)
783)
784)
785) region => region_list%first
786) do
787) if (.not.(associated(region))) exit
788)
789) if (.not.(associated(region%vertex_ids))) then
790) option%io_buffer = 'GeomechGridLocalizeRegions: define region only ' // &
791) 'by list of vertices is currently implemented: ' // &
792) trim(region%name)
793) call printErrMsg(option)
794) else
795) call GeomechGridLocalizeRegFromVertIDs(grid,region,option)
796) endif
797)
798) if (region%num_verts == 0 .and. associated(region%vertex_ids)) then
799) deallocate(region%vertex_ids)
800) nullify(region%vertex_ids)
801) endif
802)
803) region => region%next
804)
805) enddo
806)
807)
808) end subroutine GeomechGridLocalizeRegions
809)
810) ! ************************************************************************** !
811) !
812) ! GeomechGridLocalizeRegFromVertIDs: Resticts regions to vertices local
813) ! to processor for geomech grid when
814) ! the region is defined by a list of
815) ! vertex ids
816) ! author: Satish Karra
817) ! date: 06/13/13
818) !
819) ! ************************************************************************** !
820) subroutine GeomechGridLocalizeRegFromVertIDs(geomech_grid,geomech_region, &
821) option)
822)
823)
824) use Option_module
825) use Geomechanics_Region_module
826)
827) implicit none
828)
829) #include "petsc/finclude/petsclog.h"
830) #include "petsc/finclude/petscviewer.h"
831) #include "petsc/finclude/petscvec.h"
832) #include "petsc/finclude/petscvec.h90"
833) #include "petsc/finclude/petscis.h"
834) #include "petsc/finclude/petscis.h90"
835) #include "petsc/finclude/petscmat.h"
836)
837)
838) type(geomech_grid_type) :: geomech_grid
839) type(gm_region_type) :: geomech_region
840) type(option_type) :: option
841)
842) Vec :: vec_vertex_ids,vec_vertex_ids_loc
843) IS :: is_from, is_to
844) VecScatter :: vec_scat
845) PetscErrorCode :: ierr
846) PetscViewer :: viewer
847) PetscInt :: ii,jj,kk,count
848) PetscInt :: istart,iend
849) PetscInt :: ghosted_id,local_id
850) PetscInt :: natural_id
851) PetscInt, pointer :: tmp_int_array(:)
852) PetscScalar, pointer :: v_loc_p(:)
853) PetscScalar, pointer :: tmp_scl_array(:)
854) character(len=MAXSTRINGLENGTH) :: string,string1
855)
856)
857)
858) if (associated(geomech_region%vertex_ids)) then
859) call VecCreateMPI(option%mycomm,geomech_grid%nlmax_node,PETSC_DECIDE, &
860) vec_vertex_ids,ierr);CHKERRQ(ierr)
861) call VecCreateMPI(option%mycomm,geomech_grid%nlmax_node,PETSC_DECIDE,&
862) vec_vertex_ids_loc,ierr);CHKERRQ(ierr)
863) call VecZeroEntries(vec_vertex_ids,ierr);CHKERRQ(ierr)
864)
865) allocate(tmp_int_array(geomech_region%num_verts))
866) allocate(tmp_scl_array(geomech_region%num_verts))
867)
868) count = 0
869) do ii = 1, geomech_region%num_verts
870) count = count + 1
871) ! Change to zero-based numbering
872) tmp_int_array(count) = geomech_region%vertex_ids(ii) - 1
873) tmp_scl_array(count) = 1.d0
874) enddo
875)
876) #ifdef GEOMECH_DEBUG
877) call PetscViewerASCIIOpen(option%mycomm,'vec_vertex_ids_bef.out', &
878) viewer,ierr);CHKERRQ(ierr)
879) call VecView(vec_vertex_ids,viewer,ierr);CHKERRQ(ierr)
880) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
881) #endif
882)
883) call VecSetValues(vec_vertex_ids,geomech_region%num_verts,tmp_int_array, &
884) tmp_scl_array,ADD_VALUES,ierr);CHKERRQ(ierr)
885)
886) deallocate(tmp_int_array)
887) deallocate(tmp_scl_array)
888)
889) call VecAssemblyBegin(vec_vertex_ids,ierr);CHKERRQ(ierr)
890) call VecAssemblyEnd(vec_vertex_ids,ierr);CHKERRQ(ierr)
891)
892) #ifdef GEOMECH_DEBUG
893) call PetscViewerASCIIOpen(option%mycomm,'vec_vertex_ids_aft.out', &
894) viewer,ierr);CHKERRQ(ierr)
895) call VecView(vec_vertex_ids,viewer,ierr);CHKERRQ(ierr)
896) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
897) #endif
898)
899) endif
900)
901) allocate(tmp_int_array(geomech_grid%nlmax_node))
902) count = 0
903) do ghosted_id = 1, geomech_grid%ngmax_node
904) local_id = geomech_grid%nG2L(ghosted_id)
905) if (local_id < 1) cycle
906) count = count + 1
907) natural_id = geomech_grid%nG2A(ghosted_id)
908) tmp_int_array(count) = natural_id
909) enddo
910)
911) tmp_int_array = tmp_int_array - 1
912) call ISCreateBlock(option%mycomm,1,geomech_grid%nlmax_node, &
913) tmp_int_array,PETSC_COPY_VALUES,is_from, &
914) ierr);CHKERRQ(ierr)
915)
916) call VecGetOwnershipRange(vec_vertex_ids_loc,istart,iend,ierr);CHKERRQ(ierr)
917) do ii = 1,geomech_grid%nlmax_node
918) tmp_int_array(ii) = ii + istart
919) enddo
920)
921) ! is_from is natural_numbering
922) ! is_to is PETSc_numbering
923)
924) tmp_int_array = tmp_int_array - 1
925) call ISCreateBlock(option%mycomm,1,geomech_grid%nlmax_node,&
926) tmp_int_array,PETSC_COPY_VALUES,is_to,ierr);CHKERRQ(ierr)
927)
928) deallocate(tmp_int_array)
929)
930) call VecScatterCreate(vec_vertex_ids,is_from,vec_vertex_ids_loc,is_to, &
931) vec_scat,ierr);CHKERRQ(ierr)
932)
933) call ISDestroy(is_from,ierr);CHKERRQ(ierr)
934) call ISDestroy(is_to,ierr);CHKERRQ(ierr)
935)
936) call VecScatterBegin(vec_scat,vec_vertex_ids,vec_vertex_ids_loc, &
937) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
938) call VecScatterEnd(vec_scat,vec_vertex_ids,vec_vertex_ids_loc, &
939) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
940) call VecScatterDestroy(vec_scat,ierr);CHKERRQ(ierr)
941)
942) #if GEOMECH_DEBUG
943) call PetscViewerASCIIOpen(option%mycomm,'vec_vertex_ids_loc.out', &
944) viewer,ierr);CHKERRQ(ierr)
945) call VecView(vec_vertex_ids_loc,viewer,ierr);CHKERRQ(ierr)
946) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
947) #endif
948)
949) call VecGetArrayF90(vec_vertex_ids_loc,v_loc_p,ierr);CHKERRQ(ierr)
950) count = 0
951) do ii = 1,geomech_grid%nlmax_node
952) if (v_loc_p(ii) == 1) count = count + 1
953) enddo
954)
955) geomech_region%num_verts = count
956) if (count > 0) then
957) allocate(tmp_int_array(count))
958) count = 0
959) do ii = 1,geomech_grid%nlmax_node
960) if (v_loc_p(ii) == 1) then
961) count = count + 1
962) tmp_int_array(count) = ii
963) endif
964) enddo
965)
966) deallocate(geomech_region%vertex_ids)
967) allocate(geomech_region%vertex_ids(geomech_region%num_verts))
968) geomech_region%vertex_ids = tmp_int_array
969) deallocate(tmp_int_array)
970) endif
971)
972) #ifdef GEOMECH_DEBUG
973) write(string,*) option%myrank
974) write(string1,*) geomech_region%name
975) string = 'vec_region_' // trim(adjustl(string1)) // &
976) '_mapped' // trim(adjustl(string)) // '.out'
977) open(unit=86,file=trim(string))
978) do ii = 1,geomech_region%num_verts
979) write(86,'(i5)') geomech_region%vertex_ids(ii)
980) enddo
981) close(86)
982) #endif
983)
984) call VecRestoreArrayF90(vec_vertex_ids_loc,v_loc_p,ierr);CHKERRQ(ierr)
985)
986) call VecDestroy(vec_vertex_ids,ierr);CHKERRQ(ierr)
987) call VecDestroy(vec_vertex_ids_loc,ierr);CHKERRQ(ierr)
988)
989) end subroutine GeomechGridLocalizeRegFromVertIDs
990)
991) ! ************************************************************************** !
992) !
993) ! GeomechGridCopyIntegerArrayToVec: Copies values from an integer array into a
994) ! PETSc Vec
995) ! author: Satish Karra, LANL
996) ! date: 06/17/13
997) !
998) ! ************************************************************************** !
999) subroutine GeomechGridCopyIntegerArrayToVec(grid,array,vector,num_values)
1000)
1001) implicit none
1002)
1003) #include "petsc/finclude/petscvec.h"
1004) #include "petsc/finclude/petscvec.h90"
1005)
1006) type(geomech_grid_type) :: grid
1007) PetscInt :: array(:)
1008) Vec :: vector
1009) PetscInt :: num_values
1010)
1011) PetscReal, pointer :: vec_ptr(:)
1012) PetscErrorCode :: ierr
1013)
1014) call GeomechGridVecGetArrayF90(grid,vector,vec_ptr,ierr)
1015) vec_ptr(1:num_values) = array(1:num_values)
1016) call GeomechGridVecRestoreArrayF90(grid, vector,vec_ptr,ierr)
1017)
1018) end subroutine GeomechGridCopyIntegerArrayToVec
1019)
1020) ! ************************************************************************** !
1021) !
1022) ! GeomechGridVecGetArrayF90: Returns pointer to geomech veretex-based vector
1023) ! values
1024) ! author: Satish Karra, LANL
1025) ! date: 06/17/13
1026) !
1027) ! ************************************************************************** !
1028) subroutine GeomechGridVecGetArrayF90(grid,vec,f90ptr,ierr)
1029)
1030) implicit none
1031)
1032) #include "petsc/finclude/petscvec.h"
1033) #include "petsc/finclude/petscvec.h90"
1034)
1035) type(geomech_grid_type) :: grid
1036) Vec :: vec
1037) PetscReal, pointer :: f90ptr(:)
1038) PetscErrorCode :: ierr
1039)
1040) call VecGetArrayF90(vec,f90ptr,ierr);CHKERRQ(ierr)
1041)
1042) end subroutine GeomechGridVecGetArrayF90
1043)
1044) ! ************************************************************************** !
1045) !
1046) ! GeomechGridVecRestoreArrayF90: Restores pointer to geomech vector values
1047) ! author: Satish Karra, LANL
1048) ! date: 06/17/13
1049) !
1050) ! ************************************************************************** !
1051) subroutine GeomechGridVecRestoreArrayF90(grid,vec,f90ptr,ierr)
1052)
1053) implicit none
1054)
1055) #include "petsc/finclude/petscvec.h"
1056) #include "petsc/finclude/petscvec.h90"
1057)
1058) type(geomech_grid_type) :: grid
1059) Vec :: vec
1060) PetscReal, pointer :: f90ptr(:)
1061) PetscErrorCode :: ierr
1062)
1063) call VecRestoreArrayF90(vec,f90ptr,ierr);CHKERRQ(ierr)
1064)
1065) end subroutine GeomechGridVecRestoreArrayF90
1066)
1067) ! ************************************************************************** !
1068) !
1069) ! GeomechGridCopyVecToIntegerArray: Copies values from a PETSc Vec to an
1070) ! integer array
1071) ! author: Satish Karra, LANL
1072) ! date: 06/17/13
1073) !
1074) ! ************************************************************************** !
1075) subroutine GeomechGridCopyVecToIntegerArray(grid,array,vector,num_values)
1076)
1077) implicit none
1078)
1079) #include "petsc/finclude/petscvec.h"
1080) #include "petsc/finclude/petscvec.h90"
1081)
1082) type(geomech_grid_type) :: grid
1083) PetscInt :: array(:)
1084) Vec :: vector
1085) PetscInt :: num_values
1086)
1087) PetscInt :: i
1088) PetscReal, pointer :: vec_ptr(:)
1089) PetscErrorCode :: ierr
1090)
1091) call GeomechGridVecGetArrayF90(grid,vector,vec_ptr,ierr)
1092) do i = 1,num_values
1093) if (vec_ptr(i) > 0.d0) then
1094) array(i) = int(vec_ptr(i)+1.d-4)
1095) else
1096) array(i) = int(vec_ptr(i)-1.d-4)
1097) endif
1098) enddo
1099) call GeomechGridVecRestoreArrayF90(grid,vector,vec_ptr,ierr)
1100)
1101) end subroutine GeomechGridCopyVecToIntegerArray
1102)
1103) ! ************************************************************************** !
1104) !
1105) ! GeomechSubsurfMapFromFilename: Reads a list of vertex ids from a file named
1106) ! filename
1107) ! author: Satish Karra, LANL
1108) ! date: 09/07/13
1109) !
1110) ! ************************************************************************** !
1111) subroutine GeomechSubsurfMapFromFilename(grid,filename,option)
1112)
1113) use Input_Aux_module
1114) use Option_module
1115) use Utility_module
1116)
1117) implicit none
1118)
1119) type(geomech_grid_type) :: grid
1120) type(option_type) :: option
1121) type(input_type), pointer :: input
1122) character(len=MAXSTRINGLENGTH) :: filename
1123)
1124) input => InputCreate(IUNIT_TEMP,filename,option)
1125) call GeomechSubsurfMapFromFileId(grid,input,option)
1126) call InputDestroy(input)
1127)
1128) end subroutine GeomechSubsurfMapFromFilename
1129)
1130) ! ************************************************************************** !
1131) !
1132) ! GeomechSubsurfMapFromFileId: Reads a list of vertex ids from an open file
1133) ! author: Satish Karra, LANL
1134) ! date: 09/07/13
1135) !
1136) ! ************************************************************************** !
1137) subroutine GeomechSubsurfMapFromFileId(grid,input,option)
1138)
1139) use Input_Aux_module
1140) use Option_module
1141) use Utility_module
1142) use Logging_module
1143) use Grid_Unstructured_Cell_module
1144)
1145) implicit none
1146)
1147) type(geomech_grid_type) :: grid
1148) type(option_type) :: option
1149) type(input_type), pointer :: input
1150)
1151) character(len=MAXWORDLENGTH) :: word
1152) character(len=1) :: backslash
1153) character(len=MAXSTRINGLENGTH) :: string, string1
1154)
1155) PetscInt, pointer :: temp_int_array(:)
1156) PetscInt, pointer :: vertex_ids_geomech(:)
1157) PetscInt, pointer :: cell_ids_flow(:)
1158) PetscInt :: max_size, max_size_old
1159) PetscInt :: count
1160) PetscInt :: temp_int
1161) PetscInt :: input_data_type
1162) PetscInt :: ii
1163) PetscInt :: istart
1164) PetscInt :: iend
1165) PetscInt :: remainder
1166) PetscErrorCode :: ierr
1167)
1168) max_size = 1000
1169) backslash = achar(92) ! 92 = "\" Some compilers choke on \" thinking it
1170) ! is a double quote as in c/c++
1171)
1172) allocate(temp_int_array(max_size))
1173) allocate(vertex_ids_geomech(max_size))
1174) allocate(cell_ids_flow(max_size))
1175)
1176) temp_int_array = 0
1177) vertex_ids_geomech = 0
1178) cell_ids_flow = 0
1179)
1180) count = 0
1181) call InputReadPflotranString(input, option)
1182) do
1183) call InputReadInt(input, option, temp_int)
1184) if (InputError(input)) exit
1185) count = count + 1
1186) temp_int_array(count) = temp_int
1187) enddo
1188)
1189) if (count == 2) then
1190) cell_ids_flow(1) = temp_int_array(1)
1191) vertex_ids_geomech(1) = temp_int_array(2)
1192) count = 1 ! reset the counter to represent the num of rows read
1193)
1194) ! Read the data
1195) do
1196) call InputReadPflotranString(input, option)
1197) if (InputError(input)) exit
1198) call InputReadInt(input, option, temp_int)
1199) if (InputError(input)) exit
1200) count = count + 1
1201) cell_ids_flow(count) = temp_int
1202)
1203) call InputReadInt(input,option,temp_int)
1204) if (InputError(input)) then
1205) option%io_buffer = 'ERROR while reading ' // &
1206) 'GEOMECHANICS_SUBSURFACE_COUPLING mapping file.'
1207) call printErrMsg(option)
1208) endif
1209) vertex_ids_geomech(count) = temp_int
1210) if (count+1 > max_size) then ! resize temporary array
1211) max_size_old = max_size
1212) call reallocateIntArray(cell_ids_flow, max_size_old)
1213) call reallocateIntArray(vertex_ids_geomech, max_size)
1214) endif
1215) enddo
1216)
1217) ! Depending on processor rank, save only a portion of data
1218) grid%mapping_num_cells = count/option%mycommsize
1219) remainder = count - grid%mapping_num_cells*option%mycommsize
1220) if (option%myrank < remainder) grid%mapping_num_cells = &
1221) grid%mapping_num_cells + 1
1222) istart = 0
1223) iend = 0
1224) call MPI_Exscan(grid%mapping_num_cells,istart,ONE_INTEGER_MPI, &
1225) MPIU_INTEGER, &
1226) MPI_SUM,option%mycomm,ierr)
1227) call MPI_Scan(grid%mapping_num_cells,iend,ONE_INTEGER_MPI,MPIU_INTEGER, &
1228) MPI_SUM,option%mycomm,ierr)
1229)
1230) ! Allocate memory and save the data
1231) allocate(grid%mapping_cell_ids_flow(grid%mapping_num_cells))
1232) allocate(grid%mapping_vertex_ids_geomech(grid%mapping_num_cells))
1233) grid%mapping_cell_ids_flow(1:grid%mapping_num_cells) = &
1234) cell_ids_flow(istart + 1:iend)
1235) grid%mapping_vertex_ids_geomech(1:grid%mapping_num_cells) = &
1236) vertex_ids_geomech(istart + 1:iend)
1237) deallocate(cell_ids_flow)
1238) deallocate(vertex_ids_geomech)
1239) else
1240) option%io_buffer = &
1241) 'Provide a flow cell_id and a geomech vertex_id per ' // &
1242) 'line in GEOMECHANICS_SUBSURFACE_COUPLING mapping file.'
1243) call printErrMsg(option)
1244) endif
1245)
1246) deallocate(temp_int_array)
1247)
1248) #ifdef GEOMECH_DEBUG
1249) write(string,*) option%myrank
1250) string = 'geomech_subsurf_mapping_vertex_ids_geomech' &
1251) // trim(adjustl(string)) // '.out'
1252) open(unit=86,file=trim(string))
1253) do ii = 1, grid%mapping_num_cells
1254) write(86,'(i5)') grid%mapping_vertex_ids_geomech(ii)
1255) enddo
1256) close(86)
1257)
1258) write(string,*) option%myrank
1259) string = 'geomech_subsurf_mapping_cell_ids_flow' &
1260) // trim(adjustl(string)) // '.out'
1261) open(unit=86,file=trim(string))
1262) do ii = 1, grid%mapping_num_cells
1263) write(86,'(i5)') grid%mapping_cell_ids_flow(ii)
1264) enddo
1265) close(86)
1266) #endif
1267)
1268) end subroutine GeomechSubsurfMapFromFileId
1269)
1270) end module Geomechanics_Grid_module