output_geomechanics.F90 coverage: 16.67 %func 1.92 %block
1) module Output_Geomechanics_module
2)
3) use Output_Aux_module
4) use Output_Tecplot_module
5) use Output_Common_module
6) use Output_HDF5_module
7) use PFLOTRAN_Constants_module
8)
9) implicit none
10)
11) private
12)
13) #include "petsc/finclude/petscsys.h"
14) #include "petsc/finclude/petscvec.h"
15) #include "petsc/finclude/petscvec.h90"
16) #include "petsc/finclude/petscdm.h"
17) #include "petsc/finclude/petscdm.h90"
18) #include "petsc/finclude/petsclog.h"
19)
20) PetscInt, save, public :: max_local_node_size_saved = -1
21) PetscBool :: geomech_hdf5_first
22)
23) public :: OutputGeomechanics, &
24) OutputGeomechInit, &
25) OutputGeomechGetVarFromArray
26)
27) contains
28)
29) ! ************************************************************************** !
30)
31) subroutine OutputGeomechInit(num_steps)
32) !
33) ! Initializes module variables for geomechanics variables
34) !
35) ! Author: Satish Karra, LANL
36) ! Date: 07/2/13
37) !
38)
39) use Option_module
40)
41) implicit none
42)
43) PetscInt :: num_steps
44)
45) if (num_steps == 0) then
46) geomech_hdf5_first = PETSC_TRUE
47) else
48) geomech_hdf5_first = PETSC_FALSE
49) endif
50)
51) end subroutine OutputGeomechInit
52)
53) ! ************************************************************************** !
54)
55) subroutine OutputGeomechanics(geomech_realization,snapshot_plot_flag, &
56) observation_plot_flag,massbal_plot_flag)
57) !
58) ! Main driver for all geomechanics output
59) !
60) ! Author: Satish Karra, LANL
61) ! Date: 07/2/13
62) !
63)
64) use Geomechanics_Realization_class
65) use Logging_module
66) use Option_module, only : OptionCheckTouch, option_type, &
67) printMsg, printErrMsg
68)
69) implicit none
70)
71) type(realization_geomech_type) :: geomech_realization
72) PetscBool :: snapshot_plot_flag
73) PetscBool :: observation_plot_flag
74) PetscBool :: massbal_plot_flag
75)
76) character(len=MAXSTRINGLENGTH) :: string
77) PetscErrorCode :: ierr
78) PetscLogDouble :: tstart, tend
79) type(option_type), pointer :: option
80)
81) option => geomech_realization%option
82)
83) call PetscLogStagePush(logging%stage(OUTPUT_STAGE),ierr);CHKERRQ(ierr)
84)
85) ! check for plot request from active directory
86) if (.not.snapshot_plot_flag) then
87)
88) if (option%use_touch_options) then
89) string = 'plot'
90) if (OptionCheckTouch(option,string)) then
91) geomech_realization%output_option%plot_name = 'plot'
92) snapshot_plot_flag = PETSC_TRUE
93) endif
94) endif
95) endif
96)
97) !.....................................
98) if (snapshot_plot_flag) then
99)
100) if (geomech_realization%output_option%print_hdf5) then
101) call OutputHDF5UGridXDMFGeomech(geomech_realization, &
102) INSTANTANEOUS_VARS)
103) endif
104)
105) if (geomech_realization%output_option%print_tecplot) then
106) call PetscTime(tstart,ierr);CHKERRQ(ierr)
107) call OutputTecplotGeomechanics(geomech_realization)
108) call PetscTime(tend,ierr);CHKERRQ(ierr)
109) endif
110)
111) endif
112)
113) !......................................
114) if (observation_plot_flag) then
115) endif
116)
117) !......................................
118) if (massbal_plot_flag) then
119) endif
120)
121) ! Increment the plot number
122) if (snapshot_plot_flag) then
123) geomech_realization%output_option%plot_number = &
124) geomech_realization%output_option%plot_number + 1
125) endif
126)
127) call PetscLogStagePop(ierr);CHKERRQ(ierr)
128)
129) end subroutine OutputGeomechanics
130)
131) ! ************************************************************************** !
132)
133) subroutine OutputTecplotGeomechanics(geomech_realization)
134) !
135) ! Tecplot output for geomechanics
136) !
137) ! Author: Satish Karra, LANL
138) ! Date: 07/2/13
139) !
140)
141) use Geomechanics_Realization_class
142) use Geomechanics_Discretization_module
143) use Geomechanics_Grid_module
144) use Geomechanics_Grid_Aux_module
145) use Geomechanics_Field_module
146) use Geomechanics_Patch_module
147) use Option_module
148)
149) implicit none
150)
151) type(realization_geomech_type) :: geomech_realization
152)
153) PetscInt, parameter :: icolumn = -1
154) character(len=MAXSTRINGLENGTH) :: filename, string, string2
155) character(len=MAXSTRINGLENGTH) :: tmp_global_prefix
156) character(len=MAXWORDLENGTH) :: word
157) type(geomech_grid_type), pointer :: grid
158) type(option_type), pointer :: option
159) type(geomech_discretization_type), pointer :: geomech_discretization
160) type(geomech_field_type), pointer :: geomech_field
161) type(geomech_patch_type), pointer :: patch
162) type(output_option_type), pointer :: output_option
163) type(output_variable_type), pointer :: cur_variable
164) PetscReal, pointer :: vec_ptr(:)
165) Vec :: global_vertex_vec
166) Vec :: global_cconn_vec
167) Vec :: global_vec
168) Vec :: natural_vec
169) PetscInt :: ivar, isubvar, var_type
170) PetscErrorCode :: ierr
171)
172) type(gmdm_type), pointer :: gmdm_element
173)
174) geomech_discretization => geomech_realization%geomech_discretization
175) patch => geomech_realization%geomech_patch
176) grid => patch%geomech_grid
177) option => geomech_realization%option
178) geomech_field => geomech_realization%geomech_field
179) output_option => geomech_realization%output_option
180)
181) tmp_global_prefix = option%global_prefix
182) option%global_prefix = trim(tmp_global_prefix) // '-geomech'
183) filename = OutputFilename(output_option,option,'tec','')
184) option%global_prefix = tmp_global_prefix
185)
186) if (option%myrank == option%io_rank) then
187) option%io_buffer = '--> write tecplot geomech output file: ' // &
188) trim(filename)
189) call printMsg(option)
190) open(unit=OUTPUT_UNIT,file=filename,action="write")
191) call OutputTecplotHeader(OUTPUT_UNIT,geomech_realization,icolumn)
192) endif
193)
194) ! write blocks
195) ! write out data sets
196) call GeomechDiscretizationCreateVector(geomech_discretization,ONEDOF, &
197) global_vec,GLOBAL,option)
198) call GeomechDiscretizationCreateVector(geomech_discretization,ONEDOF, &
199) natural_vec,NATURAL,option)
200)
201) ! write out coordinates
202) call WriteTecplotGeomechGridVertices(OUTPUT_UNIT,geomech_realization)
203)
204) ! loop over variables and write to file
205) cur_variable => output_option%output_snap_variable_list%first
206) do
207) if (.not.associated(cur_variable)) exit
208) call OutputGeomechGetVarFromArray(geomech_realization,global_vec, &
209) cur_variable%ivar, &
210) cur_variable%isubvar)
211) call GeomechDiscretizationGlobalToNatural(geomech_discretization, &
212) global_vec, &
213) natural_vec,ONEDOF)
214) if (cur_variable%iformat == 0) then
215) call WriteTecplotDataSetGeomechFromVec(OUTPUT_UNIT,geomech_realization, &
216) natural_vec, &
217) TECPLOT_REAL)
218) else
219) call WriteTecplotDataSetGeomechFromVec(OUTPUT_UNIT,geomech_realization, &
220) natural_vec, &
221) TECPLOT_INTEGER)
222) endif
223) cur_variable => cur_variable%next
224) enddo
225)
226) call VecDestroy(natural_vec,ierr);CHKERRQ(ierr)
227) call VecDestroy(global_vec,ierr);CHKERRQ(ierr)
228)
229) ! write vertices
230) call WriteTecplotGeomechGridElements(OUTPUT_UNIT,geomech_realization)
231)
232) if (option%myrank == option%io_rank) close(OUTPUT_UNIT)
233)
234) end subroutine OutputTecplotGeomechanics
235)
236) ! ************************************************************************** !
237)
238) subroutine WriteTecplotGeomechGridElements(fid,geomech_realization)
239) !
240) ! This subroutine writes unstructured grid elements
241) ! for geomechanics grid
242) !
243) ! Author: Satish Karra
244) ! Date: 07/03/2013
245) !
246)
247) use Geomechanics_Realization_class
248) use Geomechanics_Grid_module
249) use Geomechanics_Grid_Aux_module
250) use Option_module
251) use Geomechanics_Patch_module
252)
253) implicit none
254)
255) PetscInt :: fid
256) type(realization_geomech_type) :: geomech_realization
257)
258) type(geomech_grid_type), pointer :: grid
259) type(option_type), pointer :: option
260) type(geomech_patch_type), pointer :: patch
261) Vec :: global_cconn_vec
262) type(gmdm_type), pointer :: gmdm_element
263) PetscReal, pointer :: vec_ptr(:)
264) PetscErrorCode :: ierr
265)
266) Vec :: global_vec
267) Vec :: natural_vec
268)
269) patch => geomech_realization%geomech_patch
270) grid => patch%geomech_grid
271) option => geomech_realization%option
272)
273) call GMCreateGMDM(grid,gmdm_element,EIGHT_INTEGER,option)
274) call GMGridDMCreateVectorElem(grid,gmdm_element,global_vec, &
275) GLOBAL,option)
276) call GMGridDMCreateVectorElem(grid,gmdm_element,natural_vec, &
277) NATURAL,option)
278) call GetCellConnectionsGeomech(grid,global_vec)
279) call VecScatterBegin(gmdm_element%scatter_gton_elem,global_vec,natural_vec, &
280) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
281) call VecScatterEnd(gmdm_element%scatter_gton_elem,global_vec,natural_vec, &
282) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
283) call VecGetArrayF90(natural_vec,vec_ptr,ierr);CHKERRQ(ierr)
284) call WriteTecplotDataSetNumPerLineGeomech(fid,geomech_realization,vec_ptr, &
285) TECPLOT_INTEGER, &
286) grid%nlmax_elem*8, &
287) EIGHT_INTEGER)
288) call VecRestoreArrayF90(natural_vec,vec_ptr,ierr);CHKERRQ(ierr)
289) call VecDestroy(global_vec,ierr);CHKERRQ(ierr)
290) call VecDestroy(natural_vec,ierr);CHKERRQ(ierr)
291) call GMDMDestroy(gmdm_element)
292)
293) end subroutine WriteTecplotGeomechGridElements
294)
295) ! ************************************************************************** !
296)
297) subroutine GetCellConnectionsGeomech(grid,vec)
298) !
299) ! This routine returns a vector containing vertex ids
300) ! in natural order of local cells for geomech grid
301) !
302) ! Author: Satish Karra
303) ! Date: 07/03/2013
304) !
305)
306) use Geomechanics_Grid_Aux_module
307) use Geomechanics_Grid_module
308) use Grid_Unstructured_Cell_module
309)
310) implicit none
311)
312) type(geomech_grid_type) :: grid
313) Vec :: vec
314) PetscInt :: local_id
315) PetscInt :: ghosted_id
316) PetscInt :: offset
317) PetscInt :: ivertex
318) PetscReal, pointer :: vec_ptr(:)
319) PetscErrorCode :: ierr
320)
321) call GeomechGridVecGetArrayF90(grid,vec,vec_ptr,ierr)
322)
323) ! initialize
324) vec_ptr = UNINITIALIZED_DOUBLE
325) do local_id = 1, grid%nlmax_elem
326) ghosted_id = local_id
327) select case(grid%elem_type(ghosted_id))
328) case(HEX_TYPE)
329) offset = (local_id-1)*8
330) do ivertex = 1, 8
331) vec_ptr(offset + ivertex) = &
332) grid%node_ids_ghosted_natural(grid%elem_nodes(ivertex,local_id))
333) enddo
334) case(WEDGE_TYPE)
335) offset = (local_id-1)*8
336) vec_ptr(offset + 1) = &
337) grid%node_ids_ghosted_natural(grid%elem_nodes(1,local_id))
338) vec_ptr(offset + 2) = &
339) grid%node_ids_ghosted_natural(grid%elem_nodes(1,local_id))
340) vec_ptr(offset + 3) = &
341) grid%node_ids_ghosted_natural(grid%elem_nodes(4,local_id))
342) vec_ptr(offset + 4) = &
343) grid%node_ids_ghosted_natural(grid%elem_nodes(4,local_id))
344) vec_ptr(offset + 5) = &
345) grid%node_ids_ghosted_natural(grid%elem_nodes(3,local_id))
346) vec_ptr(offset + 6) = &
347) grid%node_ids_ghosted_natural(grid%elem_nodes(2,local_id))
348) vec_ptr(offset + 7) = &
349) grid%node_ids_ghosted_natural(grid%elem_nodes(5,local_id))
350) vec_ptr(offset + 8) = &
351) grid%node_ids_ghosted_natural(grid%elem_nodes(6,local_id))
352) case (PYR_TYPE)
353) offset = (local_id-1)*8
354) ! from Tecplot 360 Data Format Guide
355) ! n1=vert1,n2=vert2,n3=vert3,n4=vert4,n5=n6=n7=n8=vert5
356) do ivertex = 1, 4
357) vec_ptr(offset + ivertex) = &
358) grid%node_ids_ghosted_natural(grid%elem_nodes(ivertex,local_id))
359) enddo
360) do ivertex = 5, 8
361) vec_ptr(offset + ivertex) = &
362) grid%node_ids_ghosted_natural(grid%elem_nodes(5,local_id))
363) enddo
364) case (TET_TYPE)
365) offset = (local_id-1)*8
366) ! from Tecplot 360 Data Format Guide
367) ! n1=vert1,n2=vert2,n3=n4=vert3,n5=vert5=n6=n7=n8=vert4
368) do ivertex = 1, 3
369) vec_ptr(offset + ivertex) = &
370) grid%node_ids_ghosted_natural(grid%elem_nodes(ivertex,local_id))
371) enddo
372) vec_ptr(offset + 4) = &
373) grid%node_ids_ghosted_natural(grid%elem_nodes(3,local_id))
374) do ivertex = 5, 8
375) vec_ptr(offset + ivertex) = &
376) grid%node_ids_ghosted_natural(grid%elem_nodes(4,local_id))
377) enddo
378) case (QUAD_TYPE)
379) offset = (local_id-1)*4
380) do ivertex = 1, 4
381) vec_ptr(offset + ivertex) = &
382) grid%node_ids_ghosted_natural(grid%elem_nodes(ivertex,local_id))
383) enddo
384) case (TRI_TYPE)
385) offset = (local_id-1)*4
386) do ivertex = 1, 3
387) vec_ptr(offset + ivertex) = &
388) grid%node_ids_ghosted_natural(grid%elem_nodes(ivertex,local_id))
389) enddo
390) ivertex = 4
391) vec_ptr(offset + ivertex) = &
392) grid%node_ids_ghosted_natural(grid%elem_nodes(3,local_id))
393) end select
394) enddo
395)
396) call GeomechGridVecRestoreArrayF90(grid,vec,vec_ptr,ierr)
397)
398) end subroutine GetCellConnectionsGeomech
399)
400) ! ************************************************************************** !
401)
402) subroutine OutputTecplotHeader(fid,geomech_realization,icolumn)
403) !
404) ! Prints Tecplot header for geomechanics
405) !
406) ! Author: Satish Karra, LANL
407) ! Date: 07/2/13
408) !
409)
410) use Geomechanics_Realization_class
411) use Geomechanics_Grid_module
412) use Geomechanics_Grid_Aux_module
413) use Option_module
414) use Geomechanics_Patch_module
415)
416) implicit none
417)
418) PetscInt :: fid
419) type(realization_geomech_type) :: geomech_realization
420) PetscInt :: icolumn
421)
422) character(len=MAXSTRINGLENGTH) :: string, string2
423) character(len=MAXWORDLENGTH) :: word
424) type(geomech_grid_type), pointer :: grid
425) type(option_type), pointer :: option
426) type(geomech_patch_type), pointer :: patch
427) type(output_option_type), pointer :: output_option
428) PetscInt :: variable_count
429)
430) patch => geomech_realization%geomech_patch
431) grid => patch%geomech_grid
432) option => geomech_realization%option
433) output_option => geomech_realization%output_option
434)
435) ! write header
436) ! write title
437) write(fid,'(''TITLE = "'',1es13.5," [",a1,'']"'')') &
438) option%geomech_time/output_option%tconv,output_option%tunit
439)
440) ! initial portion of header
441) string = 'VARIABLES=' // &
442) '"X [m]",' // &
443) '"Y [m]",' // &
444) '"Z [m]"'
445) write(fid,'(a)',advance="no") trim(string)
446)
447) call OutputWriteVariableListToHeader(fid, &
448) output_option%output_snap_variable_list, &
449) '',icolumn,PETSC_TRUE,variable_count)
450) ! need to terminate line
451) write(fid,'(a)') ''
452) ! add x, y, z variables to count
453) variable_count = variable_count + 3
454)
455) !geh: due to pgi bug, cannot embed functions with calls to write() within
456) ! write statement
457) call OutputWriteTecplotZoneHeader(fid,geomech_realization,variable_count, &
458) output_option%tecplot_format)
459)
460) end subroutine OutputTecplotHeader
461)
462) ! ************************************************************************** !
463)
464) subroutine OutputWriteTecplotZoneHeader(fid,geomech_realization, &
465) variable_count,tecplot_format)
466) !
467) ! Prints zone header to a tecplot file
468) !
469) ! Author: Satish Karra, LANL
470) ! Date: 07/2/13
471) !
472)
473) use Geomechanics_Realization_class
474) use Geomechanics_Grid_module
475) use Geomechanics_Grid_Aux_module
476) use Option_module
477) use String_module
478)
479) implicit none
480)
481) PetscInt :: fid
482) type(realization_geomech_type) :: geomech_realization
483) PetscInt :: variable_count
484) PetscInt :: tecplot_format
485)
486) character(len=MAXSTRINGLENGTH) :: string, string2, string3
487) type(geomech_grid_type), pointer :: grid
488) type(option_type), pointer :: option
489) type(output_option_type), pointer :: output_option
490)
491) grid => geomech_realization%geomech_patch%geomech_grid
492) option => geomech_realization%option
493) output_option => geomech_realization%output_option
494)
495)
496) string = 'ZONE T="' // &
497) trim(StringFormatDouble(option%time/output_option%tconv)) // &
498) '"'
499) string2 = ''
500) select case(tecplot_format)
501) case (TECPLOT_POINT_FORMAT)
502) string2 = 'POINT format not supported for geomechanics ' // &
503) 'unstructured'
504) string2 = trim(string2) // &
505) ', DATAPACKING=POINT'
506) case default !(TECPLOT_BLOCK_FORMAT,TECPLOT_FEBRICK_FORMAT)
507) string2 = ', N=' // &
508) trim(StringFormatInt(grid%nmax_node)) // &
509) ', ELEMENTS=' // &
510) trim(StringFormatInt(grid%nmax_elem))
511) string2 = trim(string2) // ', ZONETYPE=FEBRICK'
512)
513) string3 = ', VARLOCATION=(NODAL)'
514)
515) string2 = trim(string2) // trim(string3) // ', DATAPACKING=BLOCK'
516) end select
517)
518) write(fid,'(a)') trim(string) // trim(string2)
519)
520) end subroutine OutputWriteTecplotZoneHeader
521)
522) ! ************************************************************************** !
523)
524) subroutine WriteTecplotGeomechGridVertices(fid,geomech_realization)
525) !
526) ! Prints zone header to a tecplot file
527) !
528) ! Author: Satish Karra, LANL
529) ! Date: 07/2/13
530) !
531)
532) use Geomechanics_Realization_class
533) use Geomechanics_Grid_Aux_module
534) use Geomechanics_Grid_module
535) use Option_module
536) use Geomechanics_Patch_module
537) use Variables_module
538)
539) implicit none
540)
541) PetscInt :: fid
542) type(realization_geomech_type) :: geomech_realization
543)
544) type(geomech_grid_type), pointer :: grid
545) type(option_type), pointer :: option
546) type(geomech_patch_type), pointer :: patch
547) PetscReal, pointer :: vec_ptr(:)
548) Vec :: global_vertex_vec
549) PetscInt :: local_size
550) PetscErrorCode :: ierr
551)
552) patch => geomech_realization%geomech_patch
553) grid => patch%geomech_grid
554) option => geomech_realization%option
555)
556) call VecCreateMPI(option%mycomm,PETSC_DECIDE, &
557) grid%nmax_node, &
558) global_vertex_vec,ierr);CHKERRQ(ierr)
559) call VecGetLocalSize(global_vertex_vec,local_size,ierr);CHKERRQ(ierr)
560) call GetVertexCoordinatesGeomech(grid,global_vertex_vec,X_COORDINATE,option)
561) call VecGetArrayF90(global_vertex_vec,vec_ptr,ierr);CHKERRQ(ierr)
562) call WriteTecplotDataSetGeomech(fid,geomech_realization,vec_ptr, &
563) TECPLOT_REAL,local_size)
564) call VecRestoreArrayF90(global_vertex_vec,vec_ptr,ierr);CHKERRQ(ierr)
565)
566) call GetVertexCoordinatesGeomech(grid,global_vertex_vec,Y_COORDINATE,option)
567) call VecGetArrayF90(global_vertex_vec,vec_ptr,ierr);CHKERRQ(ierr)
568) call WriteTecplotDataSetGeomech(fid,geomech_realization,vec_ptr, &
569) TECPLOT_REAL,local_size)
570) call VecRestoreArrayF90(global_vertex_vec,vec_ptr,ierr);CHKERRQ(ierr)
571)
572) call GetVertexCoordinatesGeomech(grid,global_vertex_vec,Z_COORDINATE,option)
573) call VecGetArrayF90(global_vertex_vec,vec_ptr,ierr);CHKERRQ(ierr)
574) call WriteTecplotDataSetGeomech(fid,geomech_realization,vec_ptr, &
575) TECPLOT_REAL,local_size)
576) call VecRestoreArrayF90(global_vertex_vec,vec_ptr,ierr);CHKERRQ(ierr)
577)
578) call VecDestroy(global_vertex_vec,ierr);CHKERRQ(ierr)
579)
580) end subroutine WriteTecplotGeomechGridVertices
581)
582) ! ************************************************************************** !
583)
584) subroutine GetVertexCoordinatesGeomech(grid,vec,direction,option)
585) !
586) ! Extracts vertex coordinates of cells into
587) ! a PetscVec
588) !
589) ! Author: Satish Karra, LANL
590) ! Date: 07/02/2013
591) !
592)
593) use Geomechanics_Grid_module
594) use Geomechanics_Grid_Aux_module
595) use Option_module
596) use Variables_module, only : X_COORDINATE, Y_COORDINATE, Z_COORDINATE
597)
598) implicit none
599)
600) #include "petsc/finclude/petscvec.h"
601) #include "petsc/finclude/petscvec.h90"
602)
603) type(geomech_grid_type) :: grid
604) Vec :: vec
605) PetscInt :: direction
606) type(option_type) :: option
607)
608) PetscInt :: ivertex
609) PetscReal, pointer :: vec_ptr(:)
610) PetscInt, allocatable :: indices(:)
611) PetscReal, allocatable :: values(:)
612) PetscErrorCode :: ierr
613)
614) if (option%mycommsize == 1) then
615) call VecGetArrayF90(vec,vec_ptr,ierr);CHKERRQ(ierr)
616) select case(direction)
617) case(X_COORDINATE)
618) do ivertex = 1,grid%nlmax_node
619) vec_ptr(ivertex) = grid%nodes(ivertex)%x
620) enddo
621) case(Y_COORDINATE)
622) do ivertex = 1,grid%nlmax_node
623) vec_ptr(ivertex) = grid%nodes(ivertex)%y
624) enddo
625) case(Z_COORDINATE)
626) do ivertex = 1,grid%nlmax_node
627) vec_ptr(ivertex) = grid%nodes(ivertex)%z
628) enddo
629) end select
630) call VecRestoreArrayF90(vec,vec_ptr,ierr);CHKERRQ(ierr)
631) else
632) ! initialize to UNINITIALIZED_DOUBLE to catch bugs
633) call VecSet(vec,UNINITIALIZED_DOUBLE,ierr);CHKERRQ(ierr)
634) allocate(values(grid%nlmax_node))
635) allocate(indices(grid%nlmax_node))
636) select case(direction)
637) case(X_COORDINATE)
638) do ivertex = 1,grid%nlmax_node
639) values(ivertex) = grid%nodes(ivertex)%x
640) enddo
641) case(Y_COORDINATE)
642) do ivertex = 1,grid%nlmax_node
643) values(ivertex) = grid%nodes(ivertex)%y
644) enddo
645) case(Z_COORDINATE)
646) do ivertex = 1,grid%nlmax_node
647) values(ivertex) = grid%nodes(ivertex)%z
648) enddo
649) end select
650) indices(:) = grid%node_ids_local_natural(:)-1
651) call VecSetValues(vec,grid%nlmax_node, &
652) indices,values,INSERT_VALUES,ierr);CHKERRQ(ierr)
653) call VecAssemblyBegin(vec,ierr);CHKERRQ(ierr)
654) deallocate(values)
655) deallocate(indices)
656) call VecAssemblyEnd(vec,ierr);CHKERRQ(ierr)
657) endif
658)
659) end subroutine GetVertexCoordinatesGeomech
660)
661) ! ************************************************************************** !
662)
663) subroutine OutputGeomechGetVarFromArray(geomech_realization,vec,ivar,isubvar, &
664) isubvar1)
665) !
666) ! Gets variables from an array
667) !
668) ! Author: Satish Karra, LANL
669) ! Date: 07/3/13
670) !
671)
672) use Geomechanics_Realization_class
673) use Geomechanics_Grid_Aux_module
674) use Option_module
675) use Geomechanics_Field_module
676)
677) implicit none
678)
679) #include "petsc/finclude/petscvec.h"
680) #include "petsc/finclude/petscvec.h90"
681) #include "petsc/finclude/petsclog.h"
682)
683) class(realization_geomech_type) :: geomech_realization
684) Vec :: vec
685) PetscInt :: ivar
686) PetscInt :: isubvar
687) PetscInt, optional :: isubvar1
688)
689) PetscErrorCode :: ierr
690)
691) call GeomechRealizGetDataset(geomech_realization,vec,ivar,isubvar,isubvar1)
692)
693) end subroutine OutputGeomechGetVarFromArray
694)
695) ! ************************************************************************** !
696)
697) subroutine WriteTecplotDataSetGeomechFromVec(fid,geomech_realization,vec, &
698) datatype)
699) !
700) ! Writes data from a Petsc Vec within a block
701) ! of a Tecplot file
702) !
703) ! Author: Satish Karra
704) ! Date: 07/03//13
705) !
706)
707) use Geomechanics_Realization_class
708)
709) implicit none
710)
711) PetscInt :: fid
712) type(realization_geomech_type) :: geomech_realization
713) Vec :: vec
714) PetscInt :: datatype
715) PetscErrorCode :: ierr
716)
717) PetscReal, pointer :: vec_ptr(:)
718)
719) call VecGetArrayF90(vec,vec_ptr,ierr);CHKERRQ(ierr)
720) call WriteTecplotDataSetGeomech(fid,geomech_realization,vec_ptr, &
721) datatype,ZERO_INTEGER)
722) call VecRestoreArrayF90(vec,vec_ptr,ierr);CHKERRQ(ierr)
723)
724) end subroutine WriteTecplotDataSetGeomechFromVec
725)
726) ! ************************************************************************** !
727)
728) subroutine WriteTecplotDataSetGeomech(fid,geomech_realization,array,datatype, &
729) size_flag)
730) !
731) ! Writes data from an array within a block
732) ! of a Tecplot file
733) !
734) ! Author: Satish Karra
735) ! Date: 07/02//13
736) !
737)
738) use Geomechanics_Realization_class
739) use Geomechanics_Grid_Aux_module
740) use Option_module
741) use Geomechanics_Patch_module
742)
743) implicit none
744)
745) PetscInt :: fid
746) type(realization_geomech_type) :: geomech_realization
747) PetscReal :: array(:)
748) PetscInt :: datatype
749) PetscInt :: size_flag ! if size_flag /= 0, use size_flag as the local size
750)
751) PetscInt, parameter :: num_per_line = 10
752)
753) call WriteTecplotDataSetNumPerLineGeomech(fid,geomech_realization,array, &
754) datatype,size_flag,num_per_line)
755)
756) end subroutine WriteTecplotDataSetGeomech
757)
758) ! ************************************************************************** !
759)
760) subroutine WriteTecplotDataSetNumPerLineGeomech(fid,geomech_realization, &
761) array,datatype, &
762) size_flag,num_per_line)
763) !
764) ! WriteTecplotDataSetNumPerLine: Writes data from an array within a block
765) ! of a Tecplot file with a specified number
766) ! of values per line
767) !
768) ! Author: Glenn Hammond
769) ! Date: 10/25/07, 12/02/11, Satish Karra 07/02/13
770) !
771)
772) use Geomechanics_Realization_class
773) use Geomechanics_Grid_module
774) use Geomechanics_Grid_Aux_module
775) use Option_module
776) use Geomechanics_Patch_module
777)
778) implicit none
779)
780) PetscInt :: fid
781) type(realization_geomech_type) :: geomech_realization
782) PetscReal :: array(:)
783) PetscInt :: datatype
784) PetscInt :: size_flag ! if size_flag /= 0, use size_flag as the local size
785) PetscInt :: num_per_line
786)
787) type(geomech_grid_type), pointer :: grid
788) type(option_type), pointer :: option
789) type(geomech_patch_type), pointer :: patch
790) PetscInt :: i
791) PetscInt :: max_proc, max_proc_prefetch
792) PetscMPIInt :: iproc_mpi, recv_size_mpi
793) PetscInt :: max_local_size
794) PetscMPIInt :: local_size_mpi
795) PetscInt :: istart, iend, num_in_array
796) PetscMPIInt :: status_mpi(MPI_STATUS_SIZE)
797) PetscInt, allocatable :: integer_data(:), integer_data_recv(:)
798) PetscReal, allocatable :: real_data(:), real_data_recv(:)
799) PetscErrorCode :: ierr
800)
801) 1000 format(100(i2,1x))
802) 1001 format(100(i4,1x))
803) 1002 format(100(i6,1x))
804) 1003 format(100(i8,1x))
805) 1004 format(100(i10,1x))
806) 1010 format(100(es13.6,1x))
807)
808) patch => geomech_realization%geomech_patch
809) grid => patch%geomech_grid
810) option => geomech_realization%option
811)
812) ! if num_per_line exceeds 100, need to change the format statement below
813) if (num_per_line > 100) then
814) option%io_buffer = 'Number of values to be written to line in ' // &
815) 'WriteTecplotDataSetNumPerLine() exceeds 100. ' // &
816) 'Must fix format statements.'
817) call printErrMsg(option)
818) endif
819)
820) ! maximum number of initial messages
821) #define HANDSHAKE
822) max_proc = option%io_handshake_buffer_size
823) max_proc_prefetch = option%io_handshake_buffer_size / 10
824)
825) if (size_flag /= 0) then
826) call MPI_Allreduce(size_flag,max_local_size,ONE_INTEGER_MPI,MPIU_INTEGER, &
827) MPI_MAX,option%mycomm,ierr)
828) local_size_mpi = size_flag
829) else
830) ! if first time, determine the maximum size of any local array across
831) ! all procs
832) if (max_local_node_size_saved < 0) then
833) call MPI_Allreduce(grid%nlmax_node,max_local_size,ONE_INTEGER_MPI, &
834) MPIU_INTEGER,MPI_MAX,option%mycomm,ierr)
835) max_local_node_size_saved = max_local_size
836) write(option%io_buffer,'("max_local_node_size_saved: ",i9)') &
837) max_local_size
838) call printMsg(option)
839) endif
840) max_local_size = max_local_node_size_saved
841) local_size_mpi = grid%nlmax_node
842) endif
843)
844) ! transfer the data to an integer or real array
845) if (datatype == TECPLOT_INTEGER) then
846) allocate(integer_data(max_local_size+10))
847) allocate(integer_data_recv(max_local_size))
848) do i=1,local_size_mpi
849) integer_data(i) = int(array(i))
850) enddo
851) else
852) allocate(real_data(max_local_size+10))
853) allocate(real_data_recv(max_local_size))
854) do i=1,local_size_mpi
855) real_data(i) = array(i)
856) enddo
857) endif
858)
859) ! communicate data to processor 0, round robin style
860) if (option%myrank == option%io_rank) then
861) if (datatype == TECPLOT_INTEGER) then
862) ! This approach makes output files identical, regardless of processor
863) ! distribution. It is necessary when diffing files.
864) iend = 0
865) do
866) istart = iend+1
867) if (iend+num_per_line > local_size_mpi) exit
868) iend = istart+(num_per_line-1)
869) i = abs(maxval(integer_data(istart:iend)))
870) if (i < 10) then
871) write(fid,1000) integer_data(istart:iend)
872) else if (i < 1000) then
873) write(fid,1001) integer_data(istart:iend)
874) else if (i < 100000) then
875) write(fid,1002) integer_data(istart:iend)
876) else if (i < 10000000) then
877) write(fid,1003) integer_data(istart:iend)
878) else
879) write(fid,1004) integer_data(istart:iend)
880) endif
881) enddo
882) ! shift remaining data to front of array
883) integer_data(1:local_size_mpi-iend) = integer_data(iend+1:local_size_mpi)
884) num_in_array = local_size_mpi-iend
885) else
886) iend = 0
887) do
888) istart = iend+1
889) if (iend+num_per_line > local_size_mpi) exit
890) iend = istart+(num_per_line-1)
891) ! if num_per_line exceeds 100, need to change the format statement below
892) write(fid,1010) real_data(istart:iend)
893) enddo
894) ! shift remaining data to front of array
895) real_data(1:local_size_mpi-iend) = real_data(iend+1:local_size_mpi)
896) num_in_array = local_size_mpi-iend
897) endif
898) do iproc_mpi=1,option%mycommsize-1
899) #ifdef HANDSHAKE
900) if (option%io_handshake_buffer_size > 0 .and. &
901) iproc_mpi+max_proc_prefetch >= max_proc) then
902) max_proc = max_proc + option%io_handshake_buffer_size
903) call MPI_Bcast(max_proc,ONE_INTEGER_MPI,MPIU_INTEGER,option%io_rank, &
904) option%mycomm,ierr)
905) endif
906) #endif
907) call MPI_Probe(iproc_mpi,MPI_ANY_TAG,option%mycomm,status_mpi,ierr)
908) recv_size_mpi = status_mpi(MPI_TAG)
909) if (datatype == TECPLOT_INTEGER) then
910) call MPI_Recv(integer_data_recv,recv_size_mpi,MPIU_INTEGER,iproc_mpi, &
911) MPI_ANY_TAG,option%mycomm,status_mpi,ierr)
912) if (recv_size_mpi > 0) then
913) integer_data(num_in_array+1:num_in_array+recv_size_mpi) = &
914) integer_data_recv(1:recv_size_mpi)
915) num_in_array = num_in_array+recv_size_mpi
916) endif
917) iend = 0
918) do
919) istart = iend+1
920) if (iend+num_per_line > num_in_array) exit
921) iend = istart+(num_per_line-1)
922) i = abs(maxval(integer_data(istart:iend)))
923) if (i < 10) then
924) write(fid,1000) integer_data(istart:iend)
925) else if (i < 1000) then
926) write(fid,1001) integer_data(istart:iend)
927) else if (i < 100000) then
928) write(fid,1002) integer_data(istart:iend)
929) else if (i < 10000000) then
930) write(fid,1003) integer_data(istart:iend)
931) else
932) write(fid,1004) integer_data(istart:iend)
933) endif
934) enddo
935) if (iend > 0) then
936) integer_data(1:num_in_array-iend) = integer_data(iend+1:num_in_array)
937) num_in_array = num_in_array-iend
938) endif
939) else
940) call MPI_Recv(real_data_recv,recv_size_mpi, &
941) MPI_DOUBLE_PRECISION,iproc_mpi, &
942) MPI_ANY_TAG,option%mycomm,status_mpi,ierr)
943) if (recv_size_mpi > 0) then
944) real_data(num_in_array+1:num_in_array+recv_size_mpi) = &
945) real_data_recv(1:recv_size_mpi)
946) num_in_array = num_in_array+recv_size_mpi
947) endif
948) iend = 0
949) do
950) istart = iend+1
951) if (iend+num_per_line > num_in_array) exit
952) iend = istart+(num_per_line-1)
953) ! if num_per_line exceeds 100, need to change the format statement below
954) write(fid,1010) real_data(istart:iend)
955) enddo
956) if (iend > 0) then
957) real_data(1:num_in_array-iend) = real_data(iend+1:num_in_array)
958) num_in_array = num_in_array-iend
959) endif
960) endif
961) enddo
962) #ifdef HANDSHAKE
963) if (option%io_handshake_buffer_size > 0) then
964) max_proc = -1
965) call MPI_Bcast(max_proc,ONE_INTEGER_MPI,MPIU_INTEGER,option%io_rank, &
966) option%mycomm,ierr)
967) endif
968) #endif
969) ! Print the remaining values, if they exist
970) if (datatype == TECPLOT_INTEGER) then
971) if (num_in_array > 0) then
972) i = abs(maxval(integer_data(1:num_in_array)))
973) if (i < 10) then
974) write(fid,1000) integer_data(1:num_in_array)
975) else if (i < 1000) then
976) write(fid,1001) integer_data(1:num_in_array)
977) else if (i < 100000) then
978) write(fid,1002) integer_data(1:num_in_array)
979) else if (i < 10000000) then
980) write(fid,1003) integer_data(1:num_in_array)
981) else
982) write(fid,1004) integer_data(1:num_in_array)
983) endif
984) endif
985) else
986) if (num_in_array > 0) &
987) write(fid,1010) real_data(1:num_in_array)
988) endif
989) else
990) #ifdef HANDSHAKE
991) if (option%io_handshake_buffer_size > 0) then
992) do
993) if (option%myrank < max_proc) exit
994) call MPI_Bcast(max_proc,1,MPIU_INTEGER,option%io_rank,option%mycomm, &
995) ierr)
996) enddo
997) endif
998) #endif
999) if (datatype == TECPLOT_INTEGER) then
1000) call MPI_Send(integer_data,local_size_mpi,MPIU_INTEGER,option%io_rank, &
1001) local_size_mpi,option%mycomm,ierr)
1002) else
1003) call MPI_Send(real_data,local_size_mpi, &
1004) MPI_DOUBLE_PRECISION,option%io_rank, &
1005) local_size_mpi,option%mycomm,ierr)
1006) endif
1007) #ifdef HANDSHAKE
1008) if (option%io_handshake_buffer_size > 0) then
1009) do
1010) call MPI_Bcast(max_proc,1,MPIU_INTEGER,option%io_rank,option%mycomm, &
1011) ierr)
1012) if (max_proc < 0) exit
1013) enddo
1014) endif
1015) #endif
1016) #undef HANDSHAKE
1017) endif
1018)
1019) if (datatype == TECPLOT_INTEGER) then
1020) deallocate(integer_data)
1021) else
1022) deallocate(real_data)
1023) endif
1024)
1025) end subroutine WriteTecplotDataSetNumPerLineGeomech
1026)
1027) ! ************************************************************************** !
1028)
1029) subroutine OutputXMFHeaderGeomech(fid,time,nmax,xmf_vert_len,ngvert,filename)
1030) !
1031) ! This subroutine writes header to a .xmf file
1032) !
1033) ! Author: Satish Karra, LANL
1034) ! Date: 07/3/13
1035) !
1036)
1037) implicit none
1038)
1039) PetscInt :: fid, vert_count
1040) PetscReal :: time
1041) PetscInt :: nmax,xmf_vert_len,ngvert
1042) character(len=MAXSTRINGLENGTH) :: filename
1043)
1044) character(len=MAXSTRINGLENGTH) :: string, string2
1045) character(len=MAXWORDLENGTH) :: word
1046)
1047) string="<?xml version=""1.0"" ?>"
1048) write(fid,'(a)') trim(string)
1049)
1050) string="<!DOCTYPE Xdmf SYSTEM ""Xdmf.dtd"" []>"
1051) write(fid,'(a)') trim(string)
1052)
1053) string="<Xdmf>"
1054) write(fid,'(a)') trim(string)
1055)
1056) string=" <Domain>"
1057) write(fid,'(a)') trim(string)
1058)
1059) string=" <Grid Name=""Mesh"">"
1060) write(fid,'(a)') trim(string)
1061)
1062) write(string2,'(es13.5)') time
1063) string=" <Time Value = """ // trim(adjustl(string2)) // """ />"
1064) write(fid,'(a)') trim(string)
1065)
1066) write(string2,*) nmax
1067) string=" <Topology Type=""Mixed"" NumberOfElements=""" // &
1068) trim(adjustl(string2)) // """ >"
1069) write(fid,'(a)') trim(string)
1070)
1071) write(string2,*) xmf_vert_len
1072) string=" <DataItem Format=""HDF"" DataType=""Int"" Dimensions=""" // &
1073) trim(adjustl(string2)) // """>"
1074) write(fid,'(a)') trim(string)
1075)
1076) string=" "//trim(filename) //":/Domain/Cells"
1077) write(fid,'(a)') trim(string)
1078)
1079) string=" </DataItem>"
1080) write(fid,'(a)') trim(string)
1081)
1082) string=" </Topology>"
1083) write(fid,'(a)') trim(string)
1084)
1085) string=" <Geometry GeometryType=""XYZ"">"
1086) write(fid,'(a)') trim(string)
1087)
1088) write(string2,*) ngvert
1089) string=" <DataItem Format=""HDF"" Dimensions=""" // &
1090) trim(adjustl(string2)) // " 3"">"
1091) write(fid,'(a)') trim(string)
1092)
1093) string=" "//trim(filename) //":/Domain/Vertices"
1094) write(fid,'(a)') trim(string)
1095)
1096) string=" </DataItem>"
1097) write(fid,'(a)') trim(string)
1098)
1099) string=" </Geometry>"
1100) write(fid,'(a)') trim(string)
1101)
1102) #if 0
1103) string=" <Attribute Name=""X"" AttributeType=""Scalar"" Center=""Node"">"
1104) write(fid,'(a)') trim(string)
1105)
1106) write(string2,*) ngvert
1107) string=" <DataItem Dimensions=""" // &
1108) trim(adjustl(string2)) // " 1"" Format=""HDF""> "
1109) write(fid,'(a)') trim(string)
1110)
1111) string=" " // trim(filename) //":/Domain/X"
1112) write(fid,'(a)') trim(string)
1113)
1114) string=" </DataItem> "
1115) write(fid,'(a)') trim(string)
1116)
1117) string=" </Attribute>"
1118) write(fid,'(a)') trim(string)
1119) #endif
1120)
1121) end subroutine OutputXMFHeaderGeomech
1122)
1123) ! ************************************************************************** !
1124)
1125) subroutine OutputXMFFooterGeomech(fid)
1126) !
1127) ! This subroutine writes footer to a .xmf file
1128) !
1129) ! Author: Satish Karra, LANL
1130) ! Date: 07/3/13
1131) !
1132)
1133) implicit none
1134)
1135) PetscInt :: fid
1136)
1137) character(len=MAXSTRINGLENGTH) :: string
1138)
1139) string=" </Grid>"
1140) write(fid,'(a)') trim(string)
1141)
1142) string=" </Domain>"
1143) write(fid,'(a)') trim(string)
1144)
1145) string="</Xdmf>"
1146) write(fid,'(a)') trim(string)
1147)
1148) end subroutine OutputXMFFooterGeomech
1149)
1150) ! ************************************************************************** !
1151)
1152) subroutine OutputXMFAttributeGeomech(fid,nmax,attname,att_datasetname)
1153) !
1154) ! This subroutine writes an attribute to a .xmf file
1155) !
1156) ! Author: Satish Karra, LANL
1157) ! Date: 07/3/13
1158) !
1159)
1160) implicit none
1161)
1162) PetscInt :: fid,nmax
1163)
1164) character(len=MAXSTRINGLENGTH) :: attname, att_datasetname
1165) character(len=MAXSTRINGLENGTH) :: string,string2
1166) string=" <Attribute Name=""" // trim(attname) // &
1167) """ AttributeType=""Scalar"" Center=""Node"">"
1168) write(fid,'(a)') trim(string)
1169)
1170) write(string2,*) nmax
1171) string=" <DataItem Dimensions=""" // &
1172) trim(adjustl(string2)) // " 1"" Format=""HDF""> "
1173) write(fid,'(a)') trim(string)
1174)
1175) string=" " // trim(att_datasetname)
1176) write(fid,'(a)') trim(string)
1177)
1178) string=" </DataItem> "
1179) write(fid,'(a)') trim(string)
1180)
1181) string=" </Attribute>"
1182) write(fid,'(a)') trim(string)
1183)
1184) end subroutine OutputXMFAttributeGeomech
1185)
1186) ! ************************************************************************** !
1187)
1188) subroutine OutputHDF5UGridXDMFGeomech(geomech_realization,var_list_type)
1189) !
1190) ! This routine writes unstructured grid data
1191) ! in HDF5 XDMF format
1192) !
1193) ! Author: Satish Karra, LANL
1194) ! Date: 07/3/13
1195) !
1196)
1197) use Geomechanics_Realization_class
1198) use Geomechanics_Discretization_module
1199) use Option_module
1200) use Geomechanics_Grid_module
1201) use Geomechanics_Grid_Aux_module
1202) use Geomechanics_Field_module
1203) use Geomechanics_Patch_module
1204)
1205) #if !defined(PETSC_HAVE_HDF5)
1206)
1207) implicit none
1208)
1209) type(realization_geomech_type) :: geomech_realization
1210) PetscInt :: var_list_type
1211)
1212) call printMsg(geomech_realization%option,'')
1213) write(geomech_realization%option%io_buffer, &
1214) '("PFLOTRAN must be compiled with HDF5 to &
1215) &write HDF5 formatted structured grids Darn.")')
1216) call printErrMsg(geomech_realization%option)
1217)
1218) #else
1219)
1220) ! 64-bit stuff
1221) #ifdef PETSC_USE_64BIT_INDICES
1222) !#define HDF_NATIVE_INTEGER H5T_STD_I64LE
1223) #define HDF_NATIVE_INTEGER H5T_NATIVE_INTEGER
1224) #else
1225) #define HDF_NATIVE_INTEGER H5T_NATIVE_INTEGER
1226) #endif
1227)
1228) use hdf5
1229) use HDF5_module, only : HDF5WriteDataSetFromVec
1230) use HDF5_Aux_module
1231)
1232) implicit none
1233)
1234) #include "petsc/finclude/petscvec.h"
1235) #include "petsc/finclude/petscvec.h90"
1236) #include "petsc/finclude/petsclog.h"
1237)
1238) type(realization_geomech_type) :: geomech_realization
1239) PetscInt :: var_list_type
1240)
1241) #if defined(SCORPIO_WRITE)
1242) integer :: file_id
1243) integer :: data_type
1244) integer :: grp_id
1245) integer :: file_space_id
1246) integer :: memory_space_id
1247) integer :: data_set_id
1248) integer :: realization_set_id
1249) integer :: prop_id
1250) PetscMPIInt :: rank
1251) integer :: rank_mpi,file_space_rank_mpi
1252) integer :: dims(3)
1253) integer :: start(3), length(3), stride(3)
1254) #else
1255) integer(HID_T) :: file_id
1256) integer(HID_T) :: data_type
1257) integer(HID_T) :: grp_id
1258) integer(HID_T) :: file_space_id
1259) integer(HID_T) :: realization_set_id
1260) integer(HID_T) :: memory_space_id
1261) integer(HID_T) :: data_set_id
1262) integer(HID_T) :: prop_id
1263) PetscMPIInt :: rank
1264) PetscMPIInt :: rank_mpi,file_space_rank_mpi
1265) integer(HSIZE_T) :: dims(3)
1266) integer(HSIZE_T) :: start(3), length(3), stride(3)
1267) #endif
1268)
1269) type(geomech_grid_type), pointer :: grid
1270) type(geomech_discretization_type), pointer :: geomech_discretization
1271) type(geomech_field_type), pointer :: field
1272) type(geomech_patch_type), pointer :: patch
1273) type(output_option_type), pointer :: output_option
1274) type(option_type), pointer :: option
1275) type(output_variable_type), pointer :: cur_variable
1276)
1277) Vec :: global_vec
1278) Vec :: natural_vec
1279) PetscReal, pointer :: v_ptr
1280)
1281) character(len=MAXSTRINGLENGTH) :: filename
1282) character(len=MAXSTRINGLENGTH) :: xmf_filename, att_datasetname, group_name
1283) character(len=MAXSTRINGLENGTH) :: string, string2,string3
1284) character(len=MAXWORDLENGTH) :: word
1285) character(len=2) :: free_mol_char, tot_mol_char, sec_mol_char
1286) PetscReal, pointer :: array(:)
1287) PetscInt :: istart
1288) PetscInt :: i
1289) PetscInt :: nviz_flow, nviz_tran, nviz_dof
1290) PetscInt :: current_component
1291) PetscMPIInt, parameter :: ON=1, OFF=0
1292) PetscFortranAddr :: app_ptr
1293) PetscMPIInt :: hdf5_err
1294) PetscBool :: first
1295) PetscInt :: ivar, isubvar, var_type
1296) PetscInt :: vert_count
1297) PetscErrorCode :: ierr
1298)
1299) geomech_discretization => geomech_realization%geomech_discretization
1300) patch => geomech_realization%geomech_patch
1301) option => geomech_realization%option
1302) field => geomech_realization%geomech_field
1303) output_option => geomech_realization%output_option
1304)
1305) select case (var_list_type)
1306) case (INSTANTANEOUS_VARS)
1307) string2=''
1308) write(string3,'(i4)') output_option%plot_number
1309) xmf_filename = OutputFilename(output_option,option,'xmf','geomech')
1310) case (AVERAGED_VARS)
1311) string2='-aveg'
1312) write(string3,'(i4)') &
1313) int(option%time/output_option%periodic_snap_output_time_incr)
1314) xmf_filename = OutputFilename(output_option,option,'xmf','geomech_aveg')
1315) end select
1316) if (output_option%print_single_h5_file) then
1317) first = geomech_hdf5_first
1318) filename = trim(option%global_prefix) // trim(string2) // &
1319) trim(option%group_prefix) // '-geomech.h5'
1320) else
1321) string = OutputHDF5FilenameID(output_option,option,var_list_type)
1322) select case (var_list_type)
1323) case (INSTANTANEOUS_VARS)
1324) if (mod(output_option%plot_number, &
1325) output_option%times_per_h5_file)==0) then
1326) first = PETSC_TRUE
1327) else
1328) first = PETSC_FALSE
1329) endif
1330) case (AVERAGED_VARS)
1331) if (mod((option%time-output_option%periodic_snap_output_time_incr)/ &
1332) output_option%periodic_snap_output_time_incr, &
1333) dble(output_option%times_per_h5_file))==0) then
1334) first = PETSC_TRUE
1335) else
1336) first = PETSC_FALSE
1337) endif
1338) end select
1339)
1340) filename = trim(option%global_prefix) // trim(option%group_prefix) // &
1341) trim(string2) // '-' // trim(string) // '-geomech.h5'
1342) endif
1343)
1344) grid => patch%geomech_grid
1345)
1346) #ifdef SCORPIO_WRITE
1347) option%io_buffer='OutputHDF5UGridXDMF not supported with SCORPIO_WRITE'
1348) call printErrMsg(option)
1349) #endif
1350)
1351) ! initialize fortran interface
1352) call h5open_f(hdf5_err)
1353)
1354) call h5pcreate_f(H5P_FILE_ACCESS_F,prop_id,hdf5_err)
1355) #ifndef SERIAL_HDF5
1356) call h5pset_fapl_mpio_f(prop_id,option%mycomm,MPI_INFO_NULL,hdf5_err)
1357) #endif
1358)
1359) if (.not.first) then
1360) call h5eset_auto_f(OFF,hdf5_err)
1361) call h5fopen_f(filename,H5F_ACC_RDWR_F,file_id,hdf5_err,prop_id)
1362) if (hdf5_err /= 0) first = PETSC_TRUE
1363) call h5eset_auto_f(ON,hdf5_err)
1364) endif
1365) if (first) then
1366) call h5fcreate_f(filename,H5F_ACC_TRUNC_F,file_id,hdf5_err, &
1367) H5P_DEFAULT_F,prop_id)
1368) endif
1369) call h5pclose_f(prop_id,hdf5_err)
1370)
1371) if (first) then
1372) option%io_buffer = '--> creating hdf5 geomech output file: ' // &
1373) trim(filename)
1374) else
1375) option%io_buffer = '--> appending to hdf5 geomech output file: ' // &
1376) trim(filename)
1377) endif
1378) call printMsg(option)
1379)
1380) if (first) then
1381) ! create a group for the coordinates data set
1382) string = "Domain"
1383) call h5gcreate_f(file_id,string,grp_id,hdf5_err,OBJECT_NAMELEN_DEFAULT_F)
1384) call WriteHDF5CoordinatesXDMFGeomech(geomech_realization,option,grp_id)
1385) call h5gclose_f(grp_id,hdf5_err)
1386) endif
1387)
1388) if (option%myrank == option%io_rank) then
1389) option%io_buffer = '--> write xmf geomech output file: ' // &
1390) trim(xmf_filename)
1391) call printMsg(option)
1392) open(unit=OUTPUT_UNIT,file=xmf_filename,action="write")
1393) call OutputXMFHeaderGeomech(OUTPUT_UNIT, &
1394) option%time/output_option%tconv, &
1395) grid%nmax_elem, &
1396) geomech_realization%output_option%xmf_vert_len, &
1397) grid%nmax_node,filename)
1398) endif
1399)
1400) ! create a group for the data set
1401) write(string,'(''Time'',es13.5,x,a1)') &
1402) option%time/output_option%tconv,output_option%tunit
1403) if (len_trim(output_option%plot_name) > 2) then
1404) string = trim(string) // ' ' // output_option%plot_name
1405) endif
1406) string = trim(string3) // ' ' // trim(string)
1407)
1408) call h5eset_auto_f(OFF,hdf5_err)
1409) call h5gopen_f(file_id,string,grp_id,hdf5_err)
1410) group_name=string
1411) if (hdf5_err /= 0) then
1412) call h5gcreate_f(file_id,string,grp_id,hdf5_err,OBJECT_NAMELEN_DEFAULT_F)
1413) endif
1414) call h5eset_auto_f(ON,hdf5_err)
1415)
1416) ! write out data sets
1417) call GeomechDiscretizationCreateVector(geomech_discretization,ONEDOF, &
1418) global_vec, &
1419) GLOBAL,option)
1420) call GeomechDiscretizationCreateVector(geomech_discretization,ONEDOF, &
1421) natural_vec, &
1422) NATURAL,option)
1423)
1424) select case (var_list_type)
1425)
1426) case (INSTANTANEOUS_VARS)
1427) ! loop over variables and write to file
1428) cur_variable => output_option%output_snap_variable_list%first
1429) do
1430) if (.not.associated(cur_variable)) exit
1431) call OutputGeomechGetVarFromArray(geomech_realization,global_vec, &
1432) cur_variable%ivar, &
1433) cur_variable%isubvar)
1434) call GeomechDiscretizationGlobalToNatural(geomech_discretization, &
1435) global_vec, &
1436) natural_vec,ONEDOF)
1437) string = cur_variable%name
1438) if (len_trim(cur_variable%units) > 0) then
1439) word = cur_variable%units
1440) call HDF5MakeStringCompatible(word)
1441) string = trim(string) // ' [' // trim(word) // ']'
1442) endif
1443) if (cur_variable%iformat == 0) then
1444) call HDF5WriteDataSetFromVec(string,option, &
1445) natural_vec,grp_id,H5T_NATIVE_DOUBLE)
1446) else
1447) call HDF5WriteDataSetFromVec(string,option, &
1448) natural_vec,grp_id,H5T_NATIVE_INTEGER)
1449) endif
1450) att_datasetname = trim(filename) // ":/" // &
1451) trim(group_name) // "/" // trim(string)
1452) if (option%myrank == option%io_rank) then
1453) call OutputXMFAttributeGeomech(OUTPUT_UNIT,grid%nmax_node,string, &
1454) att_datasetname)
1455) endif
1456) cur_variable => cur_variable%next
1457) enddo
1458)
1459) #if 0
1460) case (AVERAGED_VARS)
1461) if (associated(output_option%aveg_output_variable_list%first)) then
1462) cur_variable => output_option%aveg_output_variable_list%first
1463) do ivar = 1,output_option%aveg_output_variable_list%nvars
1464) string = 'Aveg. ' // cur_variable%name
1465) if (len_trim(cur_variable%units) > 0) then
1466) word = cur_variable%units
1467) call HDF5MakeStringCompatible(word)
1468) string = trim(string) // ' [' // trim(word) // ']'
1469) endif
1470)
1471) call GeomechDiscretizationGlobalToNatural(geomech_discretization, &
1472) field%avg_vars_vec(ivar), &
1473) natural_vec,ONEDOF)
1474) call HDF5WriteDataSetFromVec(string,option, &
1475) natural_vec,grp_id,H5T_NATIVE_DOUBLE)
1476) att_datasetname = trim(filename) // ":/" // &
1477) trim(group_name) // "/" // trim(string)
1478) if (option%myrank == option%io_rank) then
1479) call OutputXMFAttributeGeomech(OUTPUT_UNIT,grid%nlmax_node,string, &
1480) att_datasetname)
1481) endif
1482) cur_variable => cur_variable%next
1483) enddo
1484) endif
1485) #endif
1486) ! 0
1487)
1488) end select
1489)
1490) call VecDestroy(global_vec,ierr);CHKERRQ(ierr)
1491) call VecDestroy(natural_vec,ierr);CHKERRQ(ierr)
1492) call h5gclose_f(grp_id,hdf5_err)
1493)
1494) call h5fclose_f(file_id,hdf5_err)
1495) call h5close_f(hdf5_err)
1496)
1497) if (option%myrank == option%io_rank) then
1498) call OutputXMFFooterGeomech(OUTPUT_UNIT)
1499) close(OUTPUT_UNIT)
1500) endif
1501)
1502) geomech_hdf5_first = PETSC_FALSE
1503)
1504) #endif
1505) ! !defined(PETSC_HAVE_HDF5)
1506)
1507) end subroutine OutputHDF5UGridXDMFGeomech
1508)
1509) #if defined(PETSC_HAVE_HDF5)
1510)
1511) ! ************************************************************************** !
1512)
1513) subroutine WriteHDF5CoordinatesXDMFGeomech(geomech_realization, &
1514) option,file_id)
1515) !
1516) ! Writes the geomech coordinates in HDF5 file
1517) !
1518) ! Author: Satish Karra, LANL
1519) ! Date: 07/3/13
1520) !
1521)
1522) use hdf5
1523) use HDF5_module, only : HDF5WriteDataSetFromVec
1524) use Geomechanics_Realization_class
1525) use Geomechanics_Grid_module
1526) use Geomechanics_Grid_Aux_module
1527) use Option_module
1528) use Variables_module
1529)
1530) implicit none
1531)
1532) #include "petsc/finclude/petscvec.h"
1533) #include "petsc/finclude/petscvec.h90"
1534) #include "petsc/finclude/petsclog.h"
1535)
1536) type(realization_geomech_type) :: geomech_realization
1537) type(option_type), pointer :: option
1538)
1539) #if defined(SCORPIO_WRITE)
1540) integer :: file_id
1541) integer :: data_type
1542) integer :: grp_id
1543) integer :: file_space_id
1544) integer :: memory_space_id
1545) integer :: data_set_id
1546) integer :: realization_set_id
1547) integer :: prop_id
1548) integer :: dims(3)
1549) integer :: start(3), length(3), stride(3)
1550) integer :: rank_mpi,file_space_rank_mpi
1551) integer :: hdf5_flag
1552) integer, parameter :: ON=1, OFF=0
1553) #else
1554) integer(HID_T) :: file_id
1555) integer(HID_T) :: data_type
1556) integer(HID_T) :: grp_id
1557) integer(HID_T) :: file_space_id
1558) integer(HID_T) :: realization_set_id
1559) integer(HID_T) :: memory_space_id
1560) integer(HID_T) :: data_set_id
1561) integer(HID_T) :: prop_id
1562) integer(HSIZE_T) :: dims(3)
1563) integer(HSIZE_T) :: start(3), length(3), stride(3)
1564) PetscMPIInt :: rank_mpi,file_space_rank_mpi
1565) PetscMPIInt :: hdf5_flag
1566) PetscMPIInt, parameter :: ON=1, OFF=0
1567) #endif
1568)
1569) PetscInt :: istart
1570) type(geomech_grid_type), pointer :: grid
1571) character(len=MAXSTRINGLENGTH) :: string
1572) PetscMPIInt :: hdf5_err
1573)
1574) PetscInt :: local_size,vert_count,nverts
1575) PetscInt :: i,j
1576) PetscReal, pointer :: vec_x_ptr(:),vec_y_ptr(:),vec_z_ptr(:)
1577) PetscReal, pointer :: double_array(:)
1578) Vec :: global_x_vertex_vec,global_y_vertex_vec,global_z_vertex_vec
1579) Vec :: natural_x_vertex_vec,natural_y_vertex_vec,natural_z_vertex_vec
1580)
1581) PetscReal, pointer :: vec_ptr(:)
1582) Vec :: global_vec, natural_vec
1583) PetscInt, pointer :: int_array(:)
1584) type(gmdm_type),pointer :: gmdm_element
1585) PetscErrorCode :: ierr
1586)
1587) PetscInt :: TET_ID_XDMF = 6
1588) PetscInt :: PYR_ID_XDMF = 7
1589) PetscInt :: WED_ID_XDMF = 8
1590) PetscInt :: HEX_ID_XDMF = 9
1591)
1592) grid => geomech_realization%geomech_patch%geomech_grid
1593)
1594) call VecCreateMPI(option%mycomm,PETSC_DECIDE, &
1595) grid%nmax_node, &
1596) global_x_vertex_vec,ierr);CHKERRQ(ierr)
1597) call VecCreateMPI(option%mycomm,PETSC_DECIDE, &
1598) grid%nmax_node, &
1599) global_y_vertex_vec,ierr);CHKERRQ(ierr)
1600) call VecCreateMPI(option%mycomm,PETSC_DECIDE, &
1601) grid%nmax_node, &
1602) global_z_vertex_vec,ierr);CHKERRQ(ierr)
1603)
1604) call VecGetLocalSize(global_x_vertex_vec,local_size,ierr);CHKERRQ(ierr)
1605) call VecGetLocalSize(global_y_vertex_vec,local_size,ierr);CHKERRQ(ierr)
1606) call VecGetLocalSize(global_z_vertex_vec,local_size,ierr);CHKERRQ(ierr)
1607)
1608) call GetVertexCoordinatesGeomech(grid,global_x_vertex_vec, &
1609) X_COORDINATE,option)
1610) call GetVertexCoordinatesGeomech(grid,global_y_vertex_vec, &
1611) Y_COORDINATE,option)
1612) call GetVertexCoordinatesGeomech(grid,global_z_vertex_vec, &
1613) Z_COORDINATE,option)
1614)
1615) call VecGetArrayF90(global_x_vertex_vec,vec_x_ptr,ierr);CHKERRQ(ierr)
1616) call VecGetArrayF90(global_y_vertex_vec,vec_y_ptr,ierr);CHKERRQ(ierr)
1617) call VecGetArrayF90(global_z_vertex_vec,vec_z_ptr,ierr);CHKERRQ(ierr)
1618)
1619) #if defined(SCORPIO_WRITE)
1620) write(*,*) 'SCORPIO_WRITE'
1621) option%io_buffer = 'WriteHDF5CoordinatesUGrid not supported for SCORPIO_WRITE'
1622) call printErrMsg(option)
1623) #else
1624)
1625) !
1626) ! not(SCORPIO_WRITE)
1627) !
1628)
1629) ! memory space which is a 1D vector
1630) rank_mpi = 1
1631) dims = 0
1632) dims(1) = local_size * 3
1633) call h5screate_simple_f(rank_mpi,dims,memory_space_id,hdf5_err,dims)
1634)
1635) ! file space which is a 2D block
1636) rank_mpi = 2
1637) dims = 0
1638) dims(2) = grid%nmax_node
1639) dims(1) = 3
1640) call h5pcreate_f(H5P_DATASET_CREATE_F,prop_id,hdf5_err)
1641)
1642) string = "Vertices" // CHAR(0)
1643)
1644) call h5eset_auto_f(OFF,hdf5_err)
1645) call h5dopen_f(file_id,string,data_set_id,hdf5_err)
1646) hdf5_flag = hdf5_err
1647) call h5eset_auto_f(ON,hdf5_err)
1648) if (hdf5_flag < 0) then
1649) call h5screate_simple_f(rank_mpi,dims,file_space_id,hdf5_err,dims)
1650) call h5dcreate_f(file_id,string,H5T_NATIVE_DOUBLE,file_space_id, &
1651) data_set_id,hdf5_err,prop_id)
1652) else
1653) call h5dget_space_f(data_set_id,file_space_id,hdf5_err)
1654) endif
1655)
1656) call h5pclose_f(prop_id,hdf5_err)
1657)
1658) istart = 0
1659) call MPI_Exscan(local_size, istart, ONE_INTEGER_MPI, &
1660) MPIU_INTEGER, MPI_SUM, option%mycomm, ierr)
1661)
1662) start(2) = istart
1663) start(1) = 0
1664)
1665) length(2) = local_size
1666) length(1) = 3
1667)
1668) stride = 1
1669) call h5sselect_hyperslab_f(file_space_id,H5S_SELECT_SET_F,start,length, &
1670) hdf5_err,stride,stride)
1671) ! write the data
1672) call h5pcreate_f(H5P_DATASET_XFER_F,prop_id,hdf5_err)
1673) #ifndef SERIAL_HDF5
1674) call h5pset_dxpl_mpio_f(prop_id,H5FD_MPIO_INDEPENDENT_F, &
1675) hdf5_err)
1676) #endif
1677)
1678) allocate(double_array(local_size*3))
1679) do i=1,local_size
1680) double_array((i-1)*3+1) = vec_x_ptr(i)
1681) double_array((i-1)*3+2) = vec_y_ptr(i)
1682) double_array((i-1)*3+3) = vec_z_ptr(i)
1683) enddo
1684)
1685) call h5dwrite_f(data_set_id,H5T_NATIVE_DOUBLE,double_array,dims, &
1686) hdf5_err,memory_space_id,file_space_id,prop_id)
1687)
1688) deallocate(double_array)
1689) call h5pclose_f(prop_id,hdf5_err)
1690)
1691) call h5dclose_f(data_set_id,hdf5_err)
1692) call h5sclose_f(file_space_id,hdf5_err)
1693)
1694) call VecRestoreArrayF90(global_x_vertex_vec,vec_x_ptr,ierr);CHKERRQ(ierr)
1695) call VecRestoreArrayF90(global_y_vertex_vec,vec_y_ptr,ierr);CHKERRQ(ierr)
1696) call VecRestoreArrayF90(global_z_vertex_vec,vec_z_ptr,ierr);CHKERRQ(ierr)
1697)
1698) ! Vertex X/Y/Z
1699) ! X coord
1700) string = "X" // CHAR(0)
1701)
1702) call HDF5WriteDataSetFromVec(string,option, &
1703) global_x_vertex_vec, &
1704) file_id,H5T_NATIVE_DOUBLE)
1705)
1706) ! Y coord
1707) string = "Y" // CHAR(0)
1708)
1709) call HDF5WriteDataSetFromVec(string,option, &
1710) global_y_vertex_vec, &
1711) file_id,H5T_NATIVE_DOUBLE)
1712)
1713) ! Z coord
1714) string = "Z" // CHAR(0)
1715)
1716) call HDF5WriteDataSetFromVec(string,option, &
1717) global_z_vertex_vec, &
1718) file_id,H5T_NATIVE_DOUBLE)
1719)
1720)
1721) call VecDestroy(global_x_vertex_vec,ierr);CHKERRQ(ierr)
1722) call VecDestroy(global_y_vertex_vec,ierr);CHKERRQ(ierr)
1723) call VecDestroy(global_z_vertex_vec,ierr);CHKERRQ(ierr)
1724)
1725) !
1726) ! Write elements
1727) !
1728)
1729) call GMCreateGMDM(grid,gmdm_element,EIGHT_INTEGER,option)
1730) call GMGridDMCreateVectorElem(grid,gmdm_element,global_vec, &
1731) GLOBAL,option)
1732) call GMGridDMCreateVectorElem(grid,gmdm_element,natural_vec, &
1733) NATURAL,option)
1734) call GetCellConnectionsGeomech(grid,global_vec)
1735) call VecScatterBegin(gmdm_element%scatter_gton_elem,global_vec,natural_vec, &
1736) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1737) call VecScatterEnd(gmdm_element%scatter_gton_elem,global_vec,natural_vec, &
1738) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1739) call VecGetArrayF90(natural_vec,vec_ptr,ierr);CHKERRQ(ierr)
1740)
1741) local_size = grid%nlmax_elem
1742)
1743) vert_count=0
1744) do i=1,local_size*EIGHT_INTEGER
1745) if (int(vec_ptr(i)) >0 ) vert_count=vert_count+1
1746) enddo
1747) vert_count=vert_count+grid%nlmax_elem
1748)
1749) ! memory space which is a 1D vector
1750) rank_mpi = 1
1751) dims(1) = vert_count
1752) call h5screate_simple_f(rank_mpi,dims,memory_space_id,hdf5_err,dims)
1753)
1754) call MPI_Allreduce(vert_count,dims(1),ONE_INTEGER_MPI, &
1755) MPIU_INTEGER,MPI_SUM,option%mycomm,ierr)
1756) geomech_realization%output_option%xmf_vert_len=int(dims(1))
1757)
1758) ! file space which is a 2D block
1759) rank_mpi = 1
1760) call h5pcreate_f(H5P_DATASET_CREATE_F,prop_id,hdf5_err)
1761)
1762) string = "Cells" // CHAR(0)
1763)
1764) call h5eset_auto_f(OFF,hdf5_err)
1765) call h5dopen_f(file_id,string,data_set_id,hdf5_err)
1766) hdf5_flag = hdf5_err
1767) call h5eset_auto_f(ON,hdf5_err)
1768) if (hdf5_flag < 0) then
1769) call h5screate_simple_f(rank_mpi,dims,file_space_id,hdf5_err,dims)
1770) call h5dcreate_f(file_id,string,H5T_NATIVE_INTEGER,file_space_id, &
1771) data_set_id,hdf5_err,prop_id)
1772) else
1773) call h5dget_space_f(data_set_id,file_space_id,hdf5_err)
1774) endif
1775)
1776) call h5pclose_f(prop_id,hdf5_err)
1777)
1778) istart = 0
1779) call MPI_Exscan(vert_count, istart, ONE_INTEGER_MPI, &
1780) MPIU_INTEGER, MPI_SUM, option%mycomm, ierr)
1781)
1782) start(1) = istart
1783) length(1) = vert_count
1784) stride = 1
1785) call h5sselect_hyperslab_f(file_space_id,H5S_SELECT_SET_F,start,length, &
1786) hdf5_err,stride,stride)
1787)
1788) ! write the data
1789) call h5pcreate_f(H5P_DATASET_XFER_F,prop_id,hdf5_err)
1790) #ifndef SERIAL_HDF5
1791) call h5pset_dxpl_mpio_f(prop_id,H5FD_MPIO_INDEPENDENT_F, &
1792) hdf5_err)
1793) #endif
1794)
1795) allocate(int_array(vert_count))
1796)
1797) vert_count=0
1798) do i=1,local_size
1799) nverts=0
1800) do j=1,8
1801) if (vec_ptr((i-1)*8+j)>0) nverts=nverts+1
1802) enddo
1803) vert_count=vert_count+1
1804) select case (nverts)
1805) case (4) ! Tetrahedron
1806) int_array(vert_count) = TET_ID_XDMF
1807) case (5) ! Pyramid
1808) int_array(vert_count) = PYR_ID_XDMF
1809) case (6) ! Wedge
1810) int_array(vert_count) = WED_ID_XDMF
1811) case (8) ! Hexahedron
1812) int_array(vert_count) = HEX_ID_XDMF
1813) end select
1814)
1815) do j=1,8
1816) if (vec_ptr((i-1)*8+j)>0) then
1817) vert_count=vert_count+1
1818) int_array(vert_count) = INT(vec_ptr((i-1)*8+j))-1
1819) endif
1820) enddo
1821) enddo
1822)
1823) call h5dwrite_f(data_set_id,H5T_NATIVE_INTEGER,int_array,dims, &
1824) hdf5_err,memory_space_id,file_space_id,prop_id)
1825)
1826) deallocate(int_array)
1827) call h5pclose_f(prop_id,hdf5_err)
1828)
1829) call h5dclose_f(data_set_id,hdf5_err)
1830) call h5sclose_f(file_space_id,hdf5_err)
1831)
1832) call VecRestoreArrayF90(natural_vec,vec_ptr,ierr);CHKERRQ(ierr)
1833) call VecDestroy(global_vec,ierr);CHKERRQ(ierr)
1834) call VecDestroy(natural_vec,ierr);CHKERRQ(ierr)
1835) call GMDMDestroy(gmdm_element)
1836)
1837) #endif
1838) !if defined(SCORPIO_WRITE)
1839)
1840) end subroutine WriteHDF5CoordinatesXDMFGeomech
1841) #endif
1842) ! defined(PETSC_HAVE_HDF5)
1843)
1844) end module Output_Geomechanics_module