output_surface.F90 coverage: 12.50 %func 1.97 %block
1) module Output_Surface_module
2)
3) use Logging_module
4) use Output_Aux_module
5) use Output_Common_module
6) use Output_HDF5_module
7) use Output_Tecplot_module
8)
9) use PFLOTRAN_Constants_module
10)
11) implicit none
12)
13) private
14)
15) #include "petsc/finclude/petscsys.h"
16) #include "petsc/finclude/petscvec.h"
17) #include "petsc/finclude/petscvec.h90"
18) #include "petsc/finclude/petscdm.h"
19) #include "petsc/finclude/petscdm.h90"
20) #include "petsc/finclude/petsclog.h"
21)
22) #if defined(SCORPIO_WRITE)
23) include "scorpiof.h"
24) #endif
25)
26) ! flags signifying the first time a routine is called during a given
27) ! simulation
28) PetscBool :: hydrograph_first
29) PetscBool :: surf_hdf5_first
30)
31) public :: OutputSurfaceInit, &
32) OutputSurface, &
33) OutputSurfaceVariableRead
34)
35) contains
36)
37) ! ************************************************************************** !
38)
39) subroutine OutputSurfaceInit(num_steps)
40) !
41) ! Initializes module variables for surface variables
42) ! subroutine OutputSurfaceInit(realization_base,num_steps)
43) !
44) ! Author: Glenn Hammond
45) ! Date: 01/16/13
46) !
47)
48) !use Realization_Base_class, only : realization_base_type
49) use Option_module
50)
51) implicit none
52)
53) !class(realization_base_type) :: realization_base
54) PetscInt :: num_steps
55)
56) call OutputCommonInit()
57) !call OutputObservationInit(num_steps)
58) !call OutputHDF5Init(num_steps)
59) if (num_steps == 0) then
60) hydrograph_first = PETSC_TRUE
61) surf_hdf5_first = PETSC_TRUE
62) else
63) hydrograph_first = PETSC_FALSE
64) surf_hdf5_first = PETSC_FALSE
65) endif
66)
67) end subroutine OutputSurfaceInit
68)
69) ! ************************************************************************** !
70)
71) subroutine OutputSurface(surf_realization,realization,snapshot_plot_flag, &
72) observation_plot_flag,massbal_plot_flag)
73) !
74) ! This subroutine is main driver for all output subroutines related to
75) ! surface flows.
76) !
77) ! Author: Gautam Bisht, ORNL
78) ! Date: 05/29/12
79) !
80)
81) use Realization_Surface_class, only : realization_surface_type
82) use Realization_Subsurface_class, only : realization_subsurface_type
83) use Option_module, only : OptionCheckTouch, option_type, &
84) printMsg, printErrMsg
85) use PFLOTRAN_Constants_module
86)
87) implicit none
88)
89) class(realization_surface_type) :: surf_realization
90) class(realization_subsurface_type) :: realization
91) PetscBool :: snapshot_plot_flag
92) PetscBool :: observation_plot_flag
93) PetscBool :: massbal_plot_flag
94)
95) character(len=MAXSTRINGLENGTH) :: string
96) PetscErrorCode :: ierr
97) PetscLogDouble :: tstart, tend
98) type(option_type), pointer :: option
99)
100) option => surf_realization%option
101)
102) call PetscLogStagePush(logging%stage(OUTPUT_STAGE),ierr);CHKERRQ(ierr)
103)
104) ! check for plot request from active directory
105) if (.not.snapshot_plot_flag) then
106)
107) if (option%use_touch_options) then
108) string = 'plot'
109) if (OptionCheckTouch(option,string)) then
110) surf_realization%output_option%plot_name = 'plot'
111) snapshot_plot_flag = PETSC_TRUE
112) endif
113) endif
114) endif
115)
116) !......................................
117) if (snapshot_plot_flag) then
118) if (surf_realization%output_option%print_hdf5) then
119) call OutputSurfaceHDF5UGridXDMF(surf_realization,realization, &
120) INSTANTANEOUS_VARS)
121) endif
122)
123) if (surf_realization%output_option%print_tecplot) then
124) call PetscTime(tstart,ierr);CHKERRQ(ierr)
125) call PetscLogEventBegin(logging%event_output_tecplot,ierr);CHKERRQ(ierr)
126) select case(surf_realization%output_option%tecplot_format)
127) case (TECPLOT_FEQUADRILATERAL_FORMAT)
128) call OutputTecplotFEQUAD(surf_realization,realization)
129) end select
130) call PetscTime(tend,ierr);CHKERRQ(ierr)
131) call PetscLogEventEnd(logging%event_output_tecplot,ierr);CHKERRQ(ierr)
132) endif
133)
134) if (surf_realization%output_option%print_hydrograph) then
135) call PetscTime(tstart,ierr);CHKERRQ(ierr)
136) call PetscLogEventBegin(logging%event_output_hydrograph, &
137) ierr);CHKERRQ(ierr)
138) call OutputHydrograph(surf_realization)
139) call PetscTime(tend,ierr);CHKERRQ(ierr)
140) call PetscLogEventEnd(logging%event_output_hydrograph, &
141) ierr);CHKERRQ(ierr)
142) endif
143)
144) endif
145)
146) !......................................
147) if (observation_plot_flag) then
148) endif
149)
150) !......................................
151) if (massbal_plot_flag) then
152) endif
153)
154) ! Output temporally average variables
155) call OutputSurfaceAvegVars(surf_realization,realization)
156)
157) ! Increment the plot number
158) if (snapshot_plot_flag) then
159) surf_realization%output_option%plot_number = &
160) surf_realization%output_option%plot_number + 1
161) endif
162)
163) call PetscLogStagePop(ierr);CHKERRQ(ierr)
164)
165) end subroutine OutputSurface
166)
167) ! ************************************************************************** !
168)
169) subroutine OutputTecplotFEQUAD(surf_realization,realization)
170) !
171) ! This subroutine print to Tecplot file in FEQUADRILATERAL format for surface
172) ! flows.
173) !
174) ! Author: Gautam Bisht, ORNL
175) ! Date: 05/29/12
176) !
177)
178) use Realization_Surface_class, only : realization_surface_type
179) use Realization_Subsurface_class, only : realization_subsurface_type
180) use Discretization_module
181) use Grid_module
182) use Grid_Unstructured_Aux_module
183) use Option_module
184) use Surface_Field_module
185) use Patch_module
186)
187) implicit none
188)
189) class(realization_surface_type) :: surf_realization
190) class(realization_subsurface_type) :: realization
191)
192) PetscInt :: i
193) PetscInt, parameter :: icolumn = -1
194) character(len=MAXSTRINGLENGTH) :: filename, string, string2
195) character(len=MAXSTRINGLENGTH) :: tmp_global_prefix
196) character(len=MAXWORDLENGTH) :: word
197) type(grid_type), pointer :: grid
198) type(option_type), pointer :: option
199) type(discretization_type), pointer :: discretization
200) type(surface_field_type), pointer :: surf_field
201) type(patch_type), pointer :: patch
202) type(output_option_type), pointer :: output_option
203) type(output_variable_type), pointer :: cur_variable
204) PetscReal, pointer :: vec_ptr(:)
205) Vec :: global_vertex_vec
206) Vec :: global_cconn_vec
207) Vec :: global_vec
208) Vec :: natural_vec
209) PetscInt :: ivar, isubvar, var_type
210) PetscErrorCode :: ierr
211)
212) type(ugdm_type), pointer :: ugdm_element
213)
214) discretization => surf_realization%discretization
215) patch => surf_realization%patch
216) grid => patch%grid
217) option => surf_realization%option
218) surf_field => surf_realization%surf_field
219) output_option => surf_realization%output_option
220)
221) tmp_global_prefix = option%global_prefix
222) option%global_prefix = trim(tmp_global_prefix) // '-surf'
223) filename = OutputFilename(output_option,option,'tec','')
224) option%global_prefix = tmp_global_prefix
225)
226) if (option%myrank == option%io_rank) then
227) option%io_buffer = '--> write tecplot output file: ' // trim(filename)
228) call printMsg(option)
229) open(unit=OUTPUT_UNIT,file=filename,action="write")
230) call OutputTecplotHeader(OUTPUT_UNIT,surf_realization,icolumn)
231) endif
232)
233) ! write blocks
234) ! write out data sets
235) call DiscretizationCreateVector(discretization,ONEDOF,global_vec,GLOBAL, &
236) option)
237) call DiscretizationCreateVector(discretization,ONEDOF,natural_vec,NATURAL, &
238) option)
239)
240) ! write out coordinates
241) !TODO(gautam): check this. It was flipped (realization was uncommented) before.
242) call WriteTecplotUGridVertices(OUTPUT_UNIT,surf_realization)
243) !call WriteTecplotUGridVertices(OUTPUT_UNIT,realization)
244)
245) ! loop over snapshot variables and write to file
246) cur_variable => output_option%output_snap_variable_list%first
247) do
248) if (.not.associated(cur_variable)) exit
249) call OutputSurfaceGetVarFromArray(surf_realization,global_vec,cur_variable%ivar, &
250) cur_variable%isubvar)
251) call DiscretizationGlobalToNatural(discretization,global_vec, &
252) natural_vec,ONEDOF)
253) if (cur_variable%iformat == 0) then
254) call WriteTecplotDataSetFromVec(OUTPUT_UNIT,realization,natural_vec, &
255) TECPLOT_REAL)
256) else
257) call WriteTecplotDataSetFromVec(OUTPUT_UNIT,realization,natural_vec, &
258) TECPLOT_INTEGER)
259) endif
260) cur_variable => cur_variable%next
261) enddo
262)
263) call VecDestroy(natural_vec,ierr);CHKERRQ(ierr)
264) call VecDestroy(global_vec,ierr);CHKERRQ(ierr)
265)
266) ! write vertices
267) call WriteTecplotUGridElements(OUTPUT_UNIT,surf_realization)
268)
269) if (option%myrank == option%io_rank) close(OUTPUT_UNIT)
270)
271) end subroutine OutputTecplotFEQUAD
272)
273) ! ************************************************************************** !
274)
275) subroutine OutputTecplotHeader(fid,surf_realization,icolumn)
276) !
277) ! This subroutine prints header to Tecplot file for surface grid.
278) !
279) ! Author: Gautam Bisht, ORNL
280) ! Date: 05/30/12
281) !
282)
283) use Realization_Surface_class
284) use Grid_module
285) use Option_module
286) use Patch_module
287)
288) implicit none
289)
290) PetscInt :: fid
291) class(realization_surface_type) :: surf_realization
292) PetscInt :: icolumn
293)
294) character(len=MAXSTRINGLENGTH) :: string, string2
295) character(len=MAXWORDLENGTH) :: word
296) type(grid_type), pointer :: grid
297) type(option_type), pointer :: option
298) type(patch_type), pointer :: patch
299) type(output_option_type), pointer :: output_option
300) PetscInt :: variable_count
301)
302) patch => surf_realization%patch
303) grid => patch%grid
304) option => surf_realization%option
305) output_option => surf_realization%output_option
306)
307) ! write header
308) ! write title
309) write(fid,'(''TITLE = "'',1es13.5," [",a1,'']"'')') &
310) option%surf_flow_time/output_option%tconv,output_option%tunit
311)
312) ! initial portion of header
313) string = 'VARIABLES=' // &
314) '"X [m]",' // &
315) '"Y [m]",' // &
316) '"Z [m]"'
317) write(fid,'(a)',advance='no') trim(string)
318)
319) call OutputWriteVariableListToHeader(fid, &
320) output_option%output_snap_variable_list, &
321) '',icolumn,PETSC_TRUE,variable_count)
322) ! need to terminate line
323) write(fid,'(a)') ''
324) ! add x, y, z variables to count
325) variable_count = variable_count + 3
326)
327) !geh: due to pgi bug, cannot embed functions with calls to write() within
328) ! write statement
329) call OutputWriteTecplotZoneHeader(fid,surf_realization,variable_count, &
330) output_option%tecplot_format)
331)
332) end subroutine OutputTecplotHeader
333)
334) ! ************************************************************************** !
335)
336) subroutine OutputWriteTecplotZoneHeader(fid,surf_realization,variable_count, &
337) tecplot_format)
338) !
339) ! This subroutine prints zone header to Tecplot file.
340) !
341) ! Author: Gautam Bisht, ORNL
342) ! Date: 05/30/12
343) !
344)
345) use Realization_Surface_class
346) use Grid_module
347) use Option_module
348) use String_module
349)
350) implicit none
351)
352) PetscInt :: fid
353) class(realization_surface_type) :: surf_realization
354) PetscInt :: variable_count
355) PetscInt :: tecplot_format
356)
357) character(len=MAXSTRINGLENGTH) :: string, string2, string3
358) type(grid_type), pointer :: grid
359) type(option_type), pointer :: option
360) type(output_option_type), pointer :: output_option
361)
362) grid => surf_realization%patch%grid
363) option => surf_realization%option
364) output_option => surf_realization%output_option
365)
366)
367) string = 'ZONE T="' // &
368) trim(StringFormatDouble(option%time/output_option%tconv)) // &
369) '"'
370) string2 = ''
371) select case(tecplot_format)
372) case (TECPLOT_POINT_FORMAT)
373) if (surf_realization%discretization%itype == STRUCTURED_GRID) then
374) string2 = ', I=' // &
375) trim(StringFormatInt(grid%structured_grid%nx)) // &
376) ', J=' // &
377) trim(StringFormatInt(grid%structured_grid%ny)) // &
378) ', K=' // &
379) trim(StringFormatInt(grid%structured_grid%nz))
380) else
381) string2 = 'POINT format currently not supported for unstructured'
382) endif
383) string2 = trim(string2) // &
384) ', DATAPACKING=POINT'
385) case default !(TECPLOT_BLOCK_FORMAT,TECPLOT_FEBRICK_FORMAT)
386) if (surf_realization%discretization%itype == STRUCTURED_GRID) then
387) string2 = ', I=' // &
388) trim(StringFormatInt(grid%structured_grid%nx+1)) // &
389) ', J=' // &
390) trim(StringFormatInt(grid%structured_grid%ny+1)) // &
391) ', K=' // &
392) trim(StringFormatInt(grid%structured_grid%nz+1))
393) else
394) string2 = ', N=' // &
395) trim(StringFormatInt(grid%unstructured_grid%num_vertices_global)) // &
396) ', E=' // &
397) trim(StringFormatInt(grid%unstructured_grid%nmax))
398) string2 = trim(string2) // ', ZONETYPE=FEQUADRILATERAL'
399) endif
400)
401) if (variable_count > 4) then
402) string3 = ', VARLOCATION=([4-' // &
403) trim(StringFormatInt(variable_count)) // &
404) ']=CELLCENTERED)'
405) else
406) string3 = ', VARLOCATION=([4]=CELLCENTERED)'
407) endif
408) string2 = trim(string2) // trim(string3) // ', DATAPACKING=BLOCK'
409) end select
410)
411) write(fid,'(a)') trim(string) // trim(string2)
412)
413) end subroutine OutputWriteTecplotZoneHeader
414)
415) ! ************************************************************************** !
416)
417) subroutine WriteTecplotUGridElements(fid, &
418) surf_realization)
419) !
420) ! This subroutine writes unstructured grid elements for surface grid.
421) !
422) ! Author: Gautam Bisht, ORNL
423) ! Date: 05/30/12
424) !
425)
426) use Realization_Surface_class
427) use Grid_module
428) use Grid_Unstructured_Aux_module
429) use Option_module
430) use Patch_module
431)
432) implicit none
433)
434) PetscInt :: fid
435) class(realization_surface_type) :: surf_realization
436)
437) type(grid_type), pointer :: grid
438) type(option_type), pointer :: option
439) type(patch_type), pointer :: patch
440) Vec :: global_cconn_vec
441) type(ugdm_type), pointer :: ugdm_element
442) PetscReal, pointer :: vec_ptr(:),vec_ptr2(:)
443) PetscInt :: ii
444) PetscErrorCode :: ierr
445)
446) Vec :: global_vec
447) Vec :: natural_vec
448) Vec :: surface_natural_vec
449) PetscInt :: natural_vec_local_size, surface_natural_vec_local_size
450)
451) patch => surf_realization%patch
452) grid => patch%grid
453) option => surf_realization%option
454)
455) call UGridCreateUGDM(grid%unstructured_grid,ugdm_element,EIGHT_INTEGER,option)
456) call UGridDMCreateVector(grid%unstructured_grid,ugdm_element,global_vec, &
457) GLOBAL,option)
458) call UGridDMCreateVector(grid%unstructured_grid,ugdm_element,natural_vec, &
459) NATURAL,option)
460) call GetCellConnectionsTecplot(grid,global_vec)
461) call VecScatterBegin(ugdm_element%scatter_gton,global_vec,natural_vec, &
462) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
463) call VecScatterEnd(ugdm_element%scatter_gton,global_vec,natural_vec, &
464) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
465)
466) surface_natural_vec_local_size = 0
467) call VecGetLocalSize(natural_vec,natural_vec_local_size,ierr);CHKERRQ(ierr)
468) call VecGetArrayF90(natural_vec,vec_ptr,ierr);CHKERRQ(ierr)
469) do ii = 1,natural_vec_local_size
470) if (Initialized(vec_ptr(ii))) &
471) surface_natural_vec_local_size = surface_natural_vec_local_size + 1
472) enddo
473) call VecRestoreArrayF90(natural_vec,vec_ptr,ierr);CHKERRQ(ierr)
474)
475) call VecCreateMPI(option%mycomm,surface_natural_vec_local_size,PETSC_DETERMINE,surface_natural_vec, &
476) ierr);CHKERRQ(ierr)
477) call VecGetArrayF90(surface_natural_vec,vec_ptr2,ierr);CHKERRQ(ierr)
478) call VecGetArrayF90(natural_vec,vec_ptr,ierr);CHKERRQ(ierr)
479) do ii = 1,surface_natural_vec_local_size
480) vec_ptr2(ii) = vec_ptr(ii)
481) enddo
482) call VecRestoreArrayF90(surface_natural_vec,vec_ptr2,ierr);CHKERRQ(ierr)
483) call VecRestoreArrayF90(natural_vec,vec_ptr,ierr);CHKERRQ(ierr)
484)
485) call VecGetArrayF90(surface_natural_vec,vec_ptr2,ierr);CHKERRQ(ierr)
486) call WriteTecplotDataSetNumPerLine(fid,surf_realization,vec_ptr2, &
487) TECPLOT_INTEGER, &
488) surface_natural_vec_local_size, &
489) FOUR_INTEGER)
490) call VecRestoreArrayF90(surface_natural_vec,vec_ptr2,ierr);CHKERRQ(ierr)
491) call VecDestroy(global_vec,ierr);CHKERRQ(ierr)
492) call VecDestroy(natural_vec,ierr);CHKERRQ(ierr)
493) call VecDestroy(surface_natural_vec,ierr);CHKERRQ(ierr)
494) call UGridDMDestroy(ugdm_element)
495)
496) end subroutine WriteTecplotUGridElements
497)
498) ! ************************************************************************** !
499)
500) subroutine WriteTecplotUGridVertices(fid,surf_realization)
501) !
502) ! This subroutine writes unstructure grid vertices for surface grid.
503) !
504) ! Author: Gautam Bisht, ORNL
505) ! Date: 05/30/12
506) !
507)
508) use Realization_Surface_class
509) use Grid_module
510) use Grid_Unstructured_Aux_module
511) use Option_module
512) use Patch_module
513) use Variables_module
514)
515) implicit none
516)
517) PetscInt :: fid
518) class(realization_surface_type) :: surf_realization
519)
520) type(grid_type), pointer :: grid
521) type(option_type), pointer :: option
522) type(patch_type), pointer :: patch
523) PetscReal, pointer :: vec_ptr(:)
524) Vec :: global_vertex_vec
525) PetscInt :: local_size
526) PetscErrorCode :: ierr
527)
528) patch => surf_realization%patch
529) grid => patch%grid
530) option => surf_realization%option
531)
532) call VecCreateMPI(option%mycomm,PETSC_DECIDE, &
533) grid%unstructured_grid%num_vertices_global, &
534) global_vertex_vec,ierr);CHKERRQ(ierr)
535) call VecGetLocalSize(global_vertex_vec,local_size,ierr);CHKERRQ(ierr)
536) call GetVertexCoordinates(grid, global_vertex_vec,X_COORDINATE,option)
537) call VecGetArrayF90(global_vertex_vec,vec_ptr,ierr);CHKERRQ(ierr)
538) call WriteTecplotDataSet(fid,surf_realization,vec_ptr,TECPLOT_REAL, &
539) local_size)
540) call VecRestoreArrayF90(global_vertex_vec,vec_ptr,ierr);CHKERRQ(ierr)
541)
542) call GetVertexCoordinates(grid,global_vertex_vec,Y_COORDINATE,option)
543) call VecGetArrayF90(global_vertex_vec,vec_ptr,ierr);CHKERRQ(ierr)
544) call WriteTecplotDataSet(fid,surf_realization,vec_ptr,TECPLOT_REAL, &
545) local_size)
546) call VecRestoreArrayF90(global_vertex_vec,vec_ptr,ierr);CHKERRQ(ierr)
547)
548) call GetVertexCoordinates(grid,global_vertex_vec, Z_COORDINATE,option)
549) call VecGetArrayF90(global_vertex_vec,vec_ptr,ierr);CHKERRQ(ierr)
550) call WriteTecplotDataSet(fid,surf_realization,vec_ptr,TECPLOT_REAL, &
551) local_size)
552) call VecRestoreArrayF90(global_vertex_vec,vec_ptr,ierr);CHKERRQ(ierr)
553)
554) call VecDestroy(global_vertex_vec, ierr);CHKERRQ(ierr)
555)
556) end subroutine WriteTecplotUGridVertices
557)
558) ! ************************************************************************** !
559)
560) subroutine OutputHydrograph(surf_realization)
561) !
562) ! This routine outputs hydrograph fluxes.
563) ! Author: Gautam Bisht, LBNL
564) !
565)
566) use Realization_Surface_class
567) use Patch_module
568) use Grid_module
569) use Option_module
570) use Coupler_module
571) use Connection_module
572) use Utility_module
573)
574) implicit none
575)
576) class(realization_surface_type) :: surf_realization
577) type(option_type), pointer :: option
578) type(patch_type), pointer :: patch
579) type(coupler_type), pointer :: boundary_condition
580) type(connection_set_type), pointer :: cur_connection_set
581) type(output_option_type), pointer :: output_option
582) PetscInt :: iconn
583) PetscInt :: sum_connection
584) PetscInt :: icol
585) PetscReal :: sum_flux, sum_flux_global
586)
587) character(len=MAXSTRINGLENGTH) :: filename
588) character(len=MAXWORDLENGTH) :: word, units
589) character(len=MAXSTRINGLENGTH) :: string
590) PetscErrorCode :: ierr
591)
592) PetscInt :: fid = 86
593)
594) patch => surf_realization%patch
595) option => surf_realization%option
596) output_option => surf_realization%output_option
597)
598) if (output_option%print_column_ids) then
599) icol = 1
600) else
601) icol = -1
602) endif
603)
604)
605) boundary_condition => patch%boundary_condition_list%first
606) sum_connection = 0
607) do
608) if (.not.associated(boundary_condition)) exit
609)
610) cur_connection_set => boundary_condition%connection_set
611)
612) filename = trim(boundary_condition%name) // '_hydrograph.dat'
613)
614) if (option%myrank == option%io_rank) then
615) if (hydrograph_first .or. .not.FileExists(filename)) then
616) open(unit=fid,file=filename,action="write",status="replace")
617)
618) ! write header
619) write(fid,'(a)',advance="no") ' "Time [' // &
620) trim(output_option%tunit) // ']"'
621)
622) if (option%iflowmode > 0) then
623) call OutputWriteToHeader(fid,'dt_flow',output_option%tunit,'',icol)
624) endif
625)
626) !string = 'Outflow '
627) !call OutputAppendToHeader(header,string,'[m^2/s]','',icol)
628) string = 'Outflow'
629) call OutputWriteToHeader(fid,string,'[m^3/s]','',icol)
630) write(fid,'(a)') ''
631) else
632) open(unit=fid,file=filename,action="write",status="old", &
633) position="append")
634) endif
635) endif
636)
637) 100 format(100es16.8)
638) 110 format(100es16.8)
639)
640) ! write time
641) if (option%myrank == option%io_rank) then
642) write(fid,100,advance="no") option%surf_flow_time/output_option%tconv
643) endif
644)
645) if (option%nflowdof > 0) then
646) if (option%myrank == option%io_rank) &
647) write(fid,100,advance="no") option%surf_flow_dt/output_option%tconv
648) endif
649)
650) sum_flux = 0.d0
651) do iconn = 1, cur_connection_set%num_connections
652) sum_connection = sum_connection + 1
653) !patch%boundary_velocities(1,sum_connection)
654) sum_flux = sum_flux + patch%boundary_flow_fluxes(RICHARDS_PRESSURE_DOF,sum_connection)
655) enddo
656)
657) call MPI_Reduce(sum_flux,sum_flux_global, &
658) ONE_INTEGER_MPI,MPI_DOUBLE_PRECISION,MPI_SUM, &
659) option%io_rank,option%mycomm,ierr)
660)
661) if (option%myrank == option%io_rank) then
662) ! change sign for positive in / negative out
663) write(fid,110,advance="no") -sum_flux_global
664) write(fid,'(a)') ''
665) close(fid)
666) endif
667)
668) boundary_condition => boundary_condition%next
669) enddo
670)
671) hydrograph_first = PETSC_FALSE
672)
673) end subroutine OutputHydrograph
674)
675) ! ************************************************************************** !
676)
677) subroutine OutputSurfaceHDF5UGridXDMF(surf_realization,realization, &
678) var_list_type)
679) !
680) ! This routine writes unstructured grid data in HDF5 XDMF format for
681) ! surface flow.
682) !
683) ! Author: Gautam Bisht, LBNL
684) ! Date: 10/29/2012
685) !
686)
687) use Realization_Surface_class
688) use Realization_Subsurface_class
689) use Discretization_module
690) use Option_module
691) use Grid_module
692) use Surface_Field_module
693) use Patch_module
694) use Reaction_Aux_module
695)
696) #if !defined(PETSC_HAVE_HDF5)
697) implicit none
698)
699) class(realization_surface_type) :: surf_realization
700) class(realization_subsurface_type) :: realization
701) PetscInt :: var_list_type
702)
703) call printMsg(surf_realization%option,'')
704) write(surf_realization%option%io_buffer, &
705) '("PFLOTRAN must be compiled with HDF5 to &
706) &write HDF5 formatted structured grids Darn.")')
707) call printErrMsg(surf_realization%option)
708) #else
709)
710) ! 64-bit stuff
711) #ifdef PETSC_USE_64BIT_INDICES
712) !#define HDF_NATIVE_INTEGER H5T_STD_I64LE
713) #define HDF_NATIVE_INTEGER H5T_NATIVE_INTEGER
714) #else
715) #define HDF_NATIVE_INTEGER H5T_NATIVE_INTEGER
716) #endif
717)
718) use hdf5
719) use HDF5_module
720) use HDF5_Aux_module
721)
722) implicit none
723)
724) class(realization_surface_type) :: surf_realization
725) class(realization_subsurface_type) :: realization
726) PetscInt :: var_list_type
727)
728) #if defined(SCORPIO_WRITE)
729) integer :: file_id
730) integer :: data_type
731) integer :: grp_id
732) integer :: file_space_id
733) integer :: memory_space_id
734) integer :: data_set_id
735) integer :: realization_set_id
736) integer :: prop_id
737) PetscMPIInt :: rank
738) integer :: rank_mpi,file_space_rank_mpi
739) integer :: dims(3)
740) integer :: start(3), length(3), stride(3)
741) #else
742) integer(HID_T) :: file_id
743) integer(HID_T) :: data_type
744) integer(HID_T) :: grp_id
745) integer(HID_T) :: file_space_id
746) integer(HID_T) :: realization_set_id
747) integer(HID_T) :: memory_space_id
748) integer(HID_T) :: data_set_id
749) integer(HID_T) :: prop_id
750) PetscMPIInt :: rank
751) PetscMPIInt :: rank_mpi,file_space_rank_mpi
752) integer(HSIZE_T) :: dims(3)
753) integer(HSIZE_T) :: start(3), length(3), stride(3)
754) #endif
755)
756) type(grid_type), pointer :: subsurf_grid
757) type(grid_type), pointer :: surf_grid
758) type(discretization_type), pointer :: surf_discretization
759) type(option_type), pointer :: option
760) type(surface_field_type), pointer :: surf_field
761) type(patch_type), pointer :: patch
762) type(output_option_type), pointer :: output_option
763) type(output_variable_type), pointer :: cur_variable
764)
765) Vec :: global_vec
766) Vec :: natural_vec
767) PetscReal, pointer :: v_ptr
768)
769) character(len=MAXSTRINGLENGTH) :: filename
770) character(len=MAXSTRINGLENGTH) :: xmf_filename, att_datasetname, group_name
771) character(len=MAXSTRINGLENGTH) :: string
772) character(len=MAXSTRINGLENGTH) :: string2
773) character(len=MAXSTRINGLENGTH) :: string3
774) character(len=MAXWORDLENGTH) :: word
775) character(len=2) :: free_mol_char, tot_mol_char, sec_mol_char
776) PetscReal, pointer :: array(:)
777) PetscInt :: i, istart
778) PetscInt :: nviz_flow, nviz_tran, nviz_dof
779) PetscInt :: current_component
780) PetscMPIInt, parameter :: ON=1, OFF=0
781) PetscFortranAddr :: app_ptr
782) PetscMPIInt :: hdf5_err
783) PetscBool :: first
784) PetscInt :: ivar, isubvar, var_type
785) PetscInt :: vert_count
786) PetscErrorCode :: ierr
787)
788) #ifdef SCORPIO_WRITE
789) option%io_buffer='OutputHDF5UGridXDMF not supported with SCORPIO_WRITE'
790) call printErrMsg(option)
791) #else
792)
793) surf_discretization => surf_realization%discretization
794) patch => surf_realization%patch
795) option => surf_realization%option
796) surf_field => surf_realization%surf_field
797) output_option => surf_realization%output_option
798)
799) ! xmf_filename = OutputFilename(output_option,option,'xmf','surf')
800) select case (var_list_type)
801) case (INSTANTANEOUS_VARS)
802) string2=''
803) write(string3,'(i4)') output_option%plot_number
804) xmf_filename = OutputFilename(output_option,option,'xmf','surf')
805) case (AVERAGED_VARS)
806) string2='-aveg'
807) write(string3,'(i4)') &
808) int(option%time/output_option%periodic_snap_output_time_incr)
809) xmf_filename = OutputFilename(output_option,option,'xmf','surf_aveg')
810) end select
811)
812) if (output_option%print_single_h5_file) then
813) first = surf_hdf5_first
814) filename = trim(option%global_prefix) // trim(string2) // &
815) trim(option%group_prefix) // '-surf.h5'
816) else
817) !string = OutputFilenameID(output_option,option)
818) !first = PETSC_TRUE
819) !filename = trim(option%global_prefix) // trim(option%group_prefix) // &
820) ! '-' // trim(string) // '-surf.h5'
821) string = OutputSurfaceHDF5FilenameID(output_option,option,var_list_type)
822) select case (var_list_type)
823) case (INSTANTANEOUS_VARS)
824) if (mod(output_option%plot_number,output_option%times_per_h5_file)==0) then
825) first = PETSC_TRUE
826) else
827) first = PETSC_FALSE
828) endif
829) case (AVERAGED_VARS)
830) if (mod((option%time-output_option%periodic_snap_output_time_incr)/ &
831) output_option%periodic_snap_output_time_incr, &
832) dble(output_option%times_per_h5_file))==0) then
833) first = PETSC_TRUE
834) else
835) first = PETSC_FALSE
836) endif
837) end select
838)
839) filename = trim(option%global_prefix) // trim(option%group_prefix) // &
840) trim(string2) // '-' // trim(string) // '-surf.h5'
841) endif
842)
843) subsurf_grid => realization%patch%grid
844) surf_grid => surf_realization%patch%grid
845)
846) ! initialize fortran interface
847) call h5open_f(hdf5_err)
848)
849) call h5pcreate_f(H5P_FILE_ACCESS_F,prop_id,hdf5_err)
850) #ifndef SERIAL_HDF5
851) call h5pset_fapl_mpio_f(prop_id,option%mycomm,MPI_INFO_NULL,hdf5_err)
852) #endif
853)
854) if (.not.first) then
855) call h5eset_auto_f(OFF,hdf5_err)
856) call h5fopen_f(filename,H5F_ACC_RDWR_F,file_id,hdf5_err,prop_id)
857) if (hdf5_err /= 0) first = PETSC_TRUE
858) call h5eset_auto_f(ON,hdf5_err)
859) endif
860) if (first) then
861) call h5fcreate_f(filename,H5F_ACC_TRUNC_F,file_id,hdf5_err, &
862) H5P_DEFAULT_F,prop_id)
863) endif
864) call h5pclose_f(prop_id,hdf5_err)
865)
866) if (first) then
867) option%io_buffer = '--> creating hdf5 output file: ' // trim(filename)
868) else
869) option%io_buffer = '--> appending to hdf5 output file: ' // trim(filename)
870) endif
871) call printMsg(option)
872)
873) if (first) then
874) ! create a group for the coordinates data set
875) string = "Domain"
876) call h5gcreate_f(file_id,string,grp_id,hdf5_err,OBJECT_NAMELEN_DEFAULT_F)
877) call WriteHDF5CoordinatesUGridXDMF(surf_realization,realization,option,grp_id)
878) call h5gclose_f(grp_id,hdf5_err)
879) endif
880)
881) if (option%myrank == option%io_rank) then
882) option%io_buffer = '--> write xmf output file: ' // trim(xmf_filename)
883) call printMsg(option)
884) open(unit=OUTPUT_UNIT,file=xmf_filename,action="write")
885) call OutputXMFHeader(OUTPUT_UNIT, &
886) option%time/output_option%tconv, &
887) surf_grid%nmax, &
888) surf_realization%output_option%surf_xmf_vert_len, &
889) subsurf_grid%unstructured_grid%num_vertices_global,filename)
890)
891) endif
892)
893) ! create a group for the data set
894) write(string,'(''Time'',es13.5,x,a1)') &
895) option%time/output_option%tconv,output_option%tunit
896) if (len_trim(output_option%plot_name) > 2) then
897) string = trim(string) // ' ' // output_option%plot_name
898) endif
899) string = trim(string3) // ' ' // trim(string)
900)
901) call h5eset_auto_f(OFF,hdf5_err)
902) call h5gopen_f(file_id,string,grp_id,hdf5_err)
903) group_name=string
904) if (hdf5_err /= 0) then
905) call h5gcreate_f(file_id,string,grp_id,hdf5_err,OBJECT_NAMELEN_DEFAULT_F)
906) endif
907) call h5eset_auto_f(ON,hdf5_err)
908)
909) ! write out data sets
910) call DiscretizationCreateVector(surf_discretization,ONEDOF,global_vec,GLOBAL, &
911) option)
912) call DiscretizationCreateVector(surf_discretization,ONEDOF,natural_vec,NATURAL, &
913) option)
914)
915) select case (var_list_type)
916) case (INSTANTANEOUS_VARS)
917) ! loop over snapshot variables and write to file
918) cur_variable => output_option%output_snap_variable_list%first
919) do
920) if (.not.associated(cur_variable)) exit
921) call OutputSurfaceGetVarFromArray(surf_realization,global_vec,cur_variable%ivar, &
922) cur_variable%isubvar)
923) call DiscretizationGlobalToNatural(surf_discretization,global_vec, &
924) natural_vec,ONEDOF)
925) string = cur_variable%name
926) if (len_trim(cur_variable%units) > 0) then
927) word = cur_variable%units
928) call HDF5MakeStringCompatible(word)
929) string = trim(string) // ' [' // trim(word) // ']'
930) endif
931) if (cur_variable%iformat == 0) then
932) call HDF5WriteDataSetFromVec(string,option, &
933) natural_vec,grp_id,H5T_NATIVE_DOUBLE)
934) else
935) call HDF5WriteDataSetFromVec(string,option, &
936) natural_vec,grp_id,H5T_NATIVE_INTEGER)
937) endif
938) att_datasetname = trim(filename) // ":/" // trim(group_name) // "/" // trim(string)
939) if (option%myrank == option%io_rank) then
940) call OutputXMFAttribute(OUTPUT_UNIT,surf_grid%nmax,string,att_datasetname)
941) endif
942) cur_variable => cur_variable%next
943) enddo
944) case(AVERAGED_VARS)
945) if (associated(output_option%aveg_output_variable_list%first)) then
946) cur_variable => output_option%aveg_output_variable_list%first
947) do ivar = 1,output_option%aveg_output_variable_list%nvars
948) string = 'Aveg. ' // cur_variable%name
949) if (len_trim(cur_variable%units) > 0) then
950) word = cur_variable%units
951) call HDF5MakeStringCompatible(word)
952) string = trim(string) // ' [' // trim(word) // ']'
953) endif
954)
955) call DiscretizationGlobalToNatural(surf_discretization,surf_field%avg_vars_vec(ivar), &
956) natural_vec,ONEDOF)
957) call HDF5WriteDataSetFromVec(string,option, &
958) natural_vec,grp_id,H5T_NATIVE_DOUBLE)
959) att_datasetname = trim(filename) // ":/" // trim(group_name) // "/" // trim(string)
960) if (option%myrank == option%io_rank) then
961) call OutputXMFAttribute(OUTPUT_UNIT,surf_grid%nmax,string,att_datasetname)
962) endif
963) cur_variable => cur_variable%next
964) enddo
965) endif
966) end select
967)
968) !Output flowrates
969) if (output_option%print_hdf5_mass_flowrate.or. &
970) output_option%print_hdf5_energy_flowrate.or. &
971) output_option%print_hdf5_aveg_mass_flowrate.or. &
972) output_option%print_hdf5_aveg_energy_flowrate) then
973) select case (var_list_type)
974) case (INSTANTANEOUS_VARS)
975) call OutputSurfaceGetFlowrates(surf_realization)
976) if (output_option%print_hdf5_mass_flowrate.or.&
977) output_option%print_hdf5_energy_flowrate) then
978) call WriteHDF5SurfaceFlowratesUGrid(surf_realization,grp_id,var_list_type)
979) endif
980) case (AVERAGED_VARS)
981) if (output_option%print_hdf5_aveg_mass_flowrate.or.&
982) output_option%print_hdf5_aveg_energy_flowrate) then
983) call WriteHDF5SurfaceFlowratesUGrid(surf_realization,grp_id,var_list_type)
984) endif
985) end select
986) endif
987)
988) call VecDestroy(global_vec,ierr);CHKERRQ(ierr)
989) call VecDestroy(natural_vec,ierr);CHKERRQ(ierr)
990) call h5gclose_f(grp_id,hdf5_err)
991)
992) call h5fclose_f(file_id,hdf5_err)
993) call h5close_f(hdf5_err)
994)
995) if (option%myrank == option%io_rank) then
996) call OutputXMFFooter(OUTPUT_UNIT)
997) close(OUTPUT_UNIT)
998) endif
999)
1000) surf_hdf5_first = PETSC_FALSE
1001)
1002) #endif
1003) !ifdef SCORPIO_WRITE
1004) #endif
1005) !ifdef PETSC_HAVE_HDF5
1006)
1007) end subroutine OutputSurfaceHDF5UGridXDMF
1008)
1009) #if defined(PETSC_HAVE_HDF5)
1010)
1011) ! ************************************************************************** !
1012)
1013) subroutine WriteHDF5CoordinatesUGridXDMF(surf_realization,realization, &
1014) option,file_id)
1015) !
1016) ! This routine writes unstructured coordinates to HDF5 file in XDMF format
1017) !
1018) ! Author: Gautam Bisht, LBNL
1019) ! Date: 10/29/2012
1020) !
1021)
1022) use hdf5
1023) use HDF5_module
1024) use Realization_Surface_class
1025) use Realization_Subsurface_class
1026) use Grid_module
1027) use Option_module
1028) use Grid_Unstructured_Aux_module
1029) use Variables_module
1030)
1031) implicit none
1032)
1033) class(realization_surface_type) :: surf_realization
1034) class(realization_subsurface_type) :: realization
1035) type(option_type), pointer :: option
1036)
1037) #if defined(SCORPIO_WRITE)
1038) integer :: file_id
1039) integer :: data_type
1040) integer :: grp_id
1041) integer :: file_space_id
1042) integer :: memory_space_id
1043) integer :: data_set_id
1044) integer :: realization_set_id
1045) integer :: prop_id
1046) integer :: dims(3)
1047) integer :: start(3), length(3), stride(3)
1048) integer :: rank_mpi,file_space_rank_mpi
1049) integer :: hdf5_flag
1050) integer, parameter :: ON=1, OFF=0
1051) #else
1052) integer(HID_T) :: file_id
1053) integer(HID_T) :: data_type
1054) integer(HID_T) :: grp_id
1055) integer(HID_T) :: file_space_id
1056) integer(HID_T) :: realization_set_id
1057) integer(HID_T) :: memory_space_id
1058) integer(HID_T) :: data_set_id
1059) integer(HID_T) :: prop_id
1060) integer(HSIZE_T) :: dims(3)
1061) integer(HSIZE_T) :: start(3), length(3), stride(3)
1062) PetscMPIInt :: rank_mpi,file_space_rank_mpi
1063) PetscMPIInt :: hdf5_flag
1064) PetscMPIInt, parameter :: ON=1, OFF=0
1065) #endif
1066)
1067) PetscInt :: istart
1068) PetscMPIInt :: hdf5_err
1069) type(grid_type), pointer :: surf_grid
1070) type(grid_type), pointer :: subsurf_grid
1071) character(len=MAXSTRINGLENGTH) :: string
1072)
1073) PetscInt :: local_size,vert_count,nverts
1074) PetscInt :: i,j
1075) PetscReal, pointer :: vec_x_ptr(:),vec_y_ptr(:),vec_z_ptr(:)
1076) PetscReal, pointer :: double_array(:)
1077) Vec :: global_x_vertex_vec,global_y_vertex_vec,global_z_vertex_vec
1078) Vec :: global_x_cell_vec,global_y_cell_vec,global_z_cell_vec
1079) Vec :: natural_x_cell_vec,natural_y_cell_vec,natural_z_cell_vec
1080) PetscErrorCode :: ierr
1081)
1082) PetscReal, pointer :: vec_ptr(:)
1083) Vec :: global_vec, natural_vec
1084) PetscInt, pointer :: int_array(:)
1085) type(ugdm_type),pointer :: ugdm_element,ugdm_cell
1086)
1087) PetscInt :: TRI_ID_XDMF = 4
1088) PetscInt :: QUA_ID_XDMF = 5
1089)
1090) surf_grid => surf_realization%patch%grid
1091) subsurf_grid => realization%patch%grid
1092)
1093) #if defined(SCORPIO_WRITE)
1094) write(*,*) 'SCORPIO_WRITE'
1095) option%io_buffer = 'WriteHDF5CoordinatesUGrid not supported for SCORPIO_WRITE'
1096) call printErrMsg(option)
1097) #else
1098)
1099) call VecCreateMPI(option%mycomm,PETSC_DECIDE, &
1100) subsurf_grid%unstructured_grid%num_vertices_global, &
1101) global_x_vertex_vec,ierr);CHKERRQ(ierr)
1102) call VecCreateMPI(option%mycomm,PETSC_DECIDE, &
1103) subsurf_grid%unstructured_grid%num_vertices_global, &
1104) global_y_vertex_vec,ierr);CHKERRQ(ierr)
1105) call VecCreateMPI(option%mycomm,PETSC_DECIDE, &
1106) subsurf_grid%unstructured_grid%num_vertices_global, &
1107) global_z_vertex_vec,ierr);CHKERRQ(ierr)
1108)
1109) call VecGetLocalSize(global_x_vertex_vec,local_size,ierr);CHKERRQ(ierr)
1110) call VecGetLocalSize(global_y_vertex_vec,local_size,ierr);CHKERRQ(ierr)
1111) call VecGetLocalSize(global_z_vertex_vec,local_size,ierr);CHKERRQ(ierr)
1112)
1113) call GetVertexCoordinates(subsurf_grid, global_x_vertex_vec,X_COORDINATE,option)
1114) call GetVertexCoordinates(subsurf_grid, global_y_vertex_vec,Y_COORDINATE,option)
1115) call GetVertexCoordinates(subsurf_grid, global_z_vertex_vec,Z_COORDINATE,option)
1116)
1117) call VecGetArrayF90(global_x_vertex_vec,vec_x_ptr,ierr);CHKERRQ(ierr)
1118) call VecGetArrayF90(global_y_vertex_vec,vec_y_ptr,ierr);CHKERRQ(ierr)
1119) call VecGetArrayF90(global_z_vertex_vec,vec_z_ptr,ierr);CHKERRQ(ierr)
1120)
1121) ! memory space which is a 1D vector
1122) rank_mpi = 1
1123) dims = 0
1124) dims(1) = local_size * 3
1125) call h5screate_simple_f(rank_mpi,dims,memory_space_id,hdf5_err,dims)
1126)
1127) ! file space which is a 2D block
1128) rank_mpi = 2
1129) dims = 0
1130) dims(2) = subsurf_grid%unstructured_grid%num_vertices_global
1131) dims(1) = 3
1132) call h5pcreate_f(H5P_DATASET_CREATE_F,prop_id,hdf5_err)
1133)
1134) string = "Vertices" // CHAR(0)
1135)
1136) call h5eset_auto_f(OFF,hdf5_err)
1137) call h5dopen_f(file_id,string,data_set_id,hdf5_err)
1138) hdf5_flag = hdf5_err
1139) call h5eset_auto_f(ON,hdf5_err)
1140) if (hdf5_flag < 0) then
1141) call h5screate_simple_f(rank_mpi,dims,file_space_id,hdf5_err,dims)
1142) call h5dcreate_f(file_id,string,H5T_NATIVE_DOUBLE,file_space_id, &
1143) data_set_id,hdf5_err,prop_id)
1144) else
1145) call h5dget_space_f(data_set_id,file_space_id,hdf5_err)
1146) endif
1147)
1148) call h5pclose_f(prop_id,hdf5_err)
1149)
1150) istart = 0
1151) call MPI_Exscan(local_size, istart, ONE_INTEGER_MPI, &
1152) MPIU_INTEGER, MPI_SUM, option%mycomm, ierr)
1153)
1154) start(2) = istart
1155) start(1) = 0
1156)
1157) length(2) = local_size
1158) length(1) = 3
1159)
1160) stride = 1
1161) call h5sselect_hyperslab_f(file_space_id,H5S_SELECT_SET_F,start,length, &
1162) hdf5_err,stride,stride)
1163)
1164) ! write the data
1165) call h5pcreate_f(H5P_DATASET_XFER_F,prop_id,hdf5_err)
1166) #ifndef SERIAL_HDF5
1167) call h5pset_dxpl_mpio_f(prop_id,H5FD_MPIO_INDEPENDENT_F, &
1168) hdf5_err)
1169) #endif
1170)
1171) allocate(double_array(local_size*3))
1172) do i=1,local_size
1173) double_array((i-1)*3+1) = vec_x_ptr(i)
1174) double_array((i-1)*3+2) = vec_y_ptr(i)
1175) double_array((i-1)*3+3) = vec_z_ptr(i)
1176) enddo
1177)
1178) call PetscLogEventBegin(logging%event_h5dwrite_f,ierr);CHKERRQ(ierr)
1179) call h5dwrite_f(data_set_id,H5T_NATIVE_DOUBLE,double_array,dims, &
1180) hdf5_err,memory_space_id,file_space_id,prop_id)
1181) call PetscLogEventEnd(logging%event_h5dwrite_f,ierr);CHKERRQ(ierr)
1182)
1183) deallocate(double_array)
1184) call h5pclose_f(prop_id,hdf5_err)
1185)
1186) call h5dclose_f(data_set_id,hdf5_err)
1187) call h5sclose_f(file_space_id,hdf5_err)
1188)
1189) call VecRestoreArrayF90(global_x_vertex_vec,vec_x_ptr,ierr);CHKERRQ(ierr)
1190) call VecRestoreArrayF90(global_y_vertex_vec,vec_y_ptr,ierr);CHKERRQ(ierr)
1191) call VecRestoreArrayF90(global_z_vertex_vec,vec_z_ptr,ierr);CHKERRQ(ierr)
1192)
1193)
1194) call VecDestroy(global_x_vertex_vec,ierr);CHKERRQ(ierr)
1195) call VecDestroy(global_y_vertex_vec,ierr);CHKERRQ(ierr)
1196) call VecDestroy(global_z_vertex_vec,ierr);CHKERRQ(ierr)
1197)
1198) !
1199) ! Write elements
1200) !
1201) call UGridCreateUGDM(surf_grid%unstructured_grid,ugdm_element,FOUR_INTEGER,option)
1202) call UGridDMCreateVector(surf_grid%unstructured_grid,ugdm_element,global_vec, &
1203) GLOBAL,option)
1204) call UGridDMCreateVector(surf_grid%unstructured_grid,ugdm_element,natural_vec, &
1205) NATURAL,option)
1206) call GetCellConnections(surf_grid,global_vec)
1207) call VecScatterBegin(ugdm_element%scatter_gton,global_vec,natural_vec, &
1208) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1209) call VecScatterEnd(ugdm_element%scatter_gton,global_vec,natural_vec, &
1210) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1211) call VecGetArrayF90(natural_vec,vec_ptr,ierr);CHKERRQ(ierr)
1212)
1213) local_size = surf_grid%unstructured_grid%nlmax
1214)
1215) vert_count=0
1216) do i=1,local_size*FOUR_INTEGER
1217) if (int(vec_ptr(i)) >0 ) vert_count=vert_count+1
1218) enddo
1219) vert_count=vert_count+surf_grid%nlmax
1220)
1221) ! memory space which is a 1D vector
1222) rank_mpi = 1
1223) dims(1) = vert_count
1224) call h5screate_simple_f(rank_mpi,dims,memory_space_id,hdf5_err,dims)
1225)
1226) call MPI_Allreduce(vert_count,dims(1),ONE_INTEGER_MPI,MPIU_INTEGER,MPI_SUM,option%mycomm,ierr)
1227) surf_realization%output_option%surf_xmf_vert_len=int(dims(1))
1228)
1229) ! file space which is a 2D block
1230) rank_mpi = 1
1231) call h5pcreate_f(H5P_DATASET_CREATE_F,prop_id,hdf5_err)
1232)
1233) string = "Cells" // CHAR(0)
1234)
1235) call h5eset_auto_f(OFF,hdf5_err)
1236) call h5dopen_f(file_id,string,data_set_id,hdf5_err)
1237) hdf5_flag = hdf5_err
1238) call h5eset_auto_f(ON,hdf5_err)
1239) if (hdf5_flag < 0) then
1240) call h5screate_simple_f(rank_mpi,dims,file_space_id,hdf5_err,dims)
1241) call h5dcreate_f(file_id,string,H5T_NATIVE_INTEGER,file_space_id, &
1242) data_set_id,hdf5_err,prop_id)
1243) else
1244) call h5dget_space_f(data_set_id,file_space_id,hdf5_err)
1245) endif
1246)
1247) call h5pclose_f(prop_id,hdf5_err)
1248)
1249) istart = 0
1250) call MPI_Exscan(vert_count, istart, ONE_INTEGER_MPI, &
1251) MPIU_INTEGER, MPI_SUM, option%mycomm, ierr)
1252)
1253) start(1) = istart
1254) length(1) = vert_count
1255) stride = 1
1256) call h5sselect_hyperslab_f(file_space_id,H5S_SELECT_SET_F,start,length, &
1257) hdf5_err,stride,stride)
1258)
1259) ! write the data
1260) call h5pcreate_f(H5P_DATASET_XFER_F,prop_id,hdf5_err)
1261) #ifndef SERIAL_HDF5
1262) call h5pset_dxpl_mpio_f(prop_id,H5FD_MPIO_INDEPENDENT_F, &
1263) hdf5_err)
1264) #endif
1265)
1266) allocate(int_array(vert_count))
1267)
1268) vert_count=0
1269) do i=1,local_size
1270) nverts=0
1271) do j=1,FOUR_INTEGER
1272) if (vec_ptr((i-1)*FOUR_INTEGER+j)>0) nverts=nverts+1
1273) enddo
1274) vert_count=vert_count+1
1275) select case (nverts)
1276) case (3) ! Triangle
1277) int_array(vert_count) = TRI_ID_XDMF
1278) case (4) ! Quadrilateral
1279) int_array(vert_count) = QUA_ID_XDMF
1280) end select
1281)
1282) do j=1,FOUR_INTEGER
1283) if (vec_ptr((i-1)*FOUR_INTEGER+j)>0) then
1284) vert_count=vert_count+1
1285) int_array(vert_count) = INT(vec_ptr((i-1)*FOUR_INTEGER+j))-1
1286) endif
1287) enddo
1288) enddo
1289)
1290) call PetscLogEventBegin(logging%event_h5dwrite_f,ierr);CHKERRQ(ierr)
1291) call h5dwrite_f(data_set_id,H5T_NATIVE_INTEGER,int_array,dims, &
1292) hdf5_err,memory_space_id,file_space_id,prop_id)
1293) call PetscLogEventEnd(logging%event_h5dwrite_f,ierr);CHKERRQ(ierr)
1294)
1295) deallocate(int_array)
1296) call h5pclose_f(prop_id,hdf5_err)
1297)
1298) call h5dclose_f(data_set_id,hdf5_err)
1299) call h5sclose_f(file_space_id,hdf5_err)
1300)
1301) call VecRestoreArrayF90(natural_vec,vec_ptr,ierr);CHKERRQ(ierr)
1302) call VecDestroy(global_vec,ierr);CHKERRQ(ierr)
1303) call VecDestroy(natural_vec,ierr);CHKERRQ(ierr)
1304) call UGridDMDestroy(ugdm_element)
1305)
1306) ! Cell center X/Y/Z
1307) call VecCreateMPI(option%mycomm,surf_grid%nlmax, &
1308) PETSC_DETERMINE, &
1309) global_x_cell_vec,ierr);CHKERRQ(ierr)
1310) call VecCreateMPI(option%mycomm,surf_grid%nlmax, &
1311) PETSC_DETERMINE, &
1312) global_y_cell_vec,ierr);CHKERRQ(ierr)
1313) call VecCreateMPI(option%mycomm,surf_grid%nlmax, &
1314) PETSC_DETERMINE, &
1315) global_z_cell_vec,ierr);CHKERRQ(ierr)
1316)
1317) call GetCellCoordinates(surf_grid, global_x_cell_vec,X_COORDINATE)
1318) call GetCellCoordinates(surf_grid, global_y_cell_vec,Y_COORDINATE)
1319) call GetCellCoordinates(surf_grid, global_z_cell_vec,Z_COORDINATE)
1320)
1321)
1322) call UGridCreateUGDM(surf_grid%unstructured_grid,ugdm_cell,ONE_INTEGER,option)
1323) call UGridDMCreateVector(surf_grid%unstructured_grid,ugdm_cell,natural_x_cell_vec, &
1324) NATURAL,option)
1325) call UGridDMCreateVector(surf_grid%unstructured_grid,ugdm_cell,natural_y_cell_vec, &
1326) NATURAL,option)
1327) call UGridDMCreateVector(surf_grid%unstructured_grid,ugdm_cell,natural_z_cell_vec, &
1328) NATURAL,option)
1329)
1330) call VecScatterBegin(ugdm_cell%scatter_gton,global_x_cell_vec,natural_x_cell_vec, &
1331) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1332) call VecScatterEnd(ugdm_cell%scatter_gton,global_x_cell_vec,natural_x_cell_vec, &
1333) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1334)
1335) call VecScatterBegin(ugdm_cell%scatter_gton,global_y_cell_vec,natural_y_cell_vec, &
1336) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1337) call VecScatterEnd(ugdm_cell%scatter_gton,global_y_cell_vec,natural_y_cell_vec, &
1338) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1339)
1340) call VecScatterBegin(ugdm_cell%scatter_gton,global_z_cell_vec,natural_z_cell_vec, &
1341) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1342) call VecScatterEnd(ugdm_cell%scatter_gton,global_z_cell_vec,natural_z_cell_vec, &
1343) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1344)
1345) call VecGetArrayF90(natural_x_cell_vec,vec_x_ptr,ierr);CHKERRQ(ierr)
1346) call VecGetArrayF90(natural_y_cell_vec,vec_y_ptr,ierr);CHKERRQ(ierr)
1347) call VecGetArrayF90(natural_z_cell_vec,vec_z_ptr,ierr);CHKERRQ(ierr)
1348) local_size = surf_grid%unstructured_grid%nlmax
1349)
1350) ! XC
1351) ! memory space which is a 1D vector
1352) rank_mpi = 1
1353) dims = 0
1354) dims(1) = local_size
1355) call h5screate_simple_f(rank_mpi,dims,memory_space_id,hdf5_err,dims)
1356)
1357) ! file space which is a 2D block
1358) rank_mpi = 1
1359) dims = 0
1360) dims(1) = surf_grid%nmax
1361) call h5pcreate_f(H5P_DATASET_CREATE_F,prop_id,hdf5_err)
1362)
1363) string = "XC" // CHAR(0)
1364)
1365) call h5eset_auto_f(OFF,hdf5_err)
1366) call h5dopen_f(file_id,string,data_set_id,hdf5_err)
1367) hdf5_flag = hdf5_err
1368) call h5eset_auto_f(ON,hdf5_err)
1369) if (hdf5_flag < 0) then
1370) call h5screate_simple_f(rank_mpi,dims,file_space_id,hdf5_err,dims)
1371) call h5dcreate_f(file_id,string,H5T_NATIVE_DOUBLE,file_space_id, &
1372) data_set_id,hdf5_err,prop_id)
1373) else
1374) call h5dget_space_f(data_set_id,file_space_id,hdf5_err)
1375) endif
1376)
1377) call h5pclose_f(prop_id,hdf5_err)
1378)
1379) istart = 0
1380) call MPI_Exscan(local_size, istart, ONE_INTEGER_MPI, &
1381) MPIU_INTEGER, MPI_SUM, option%mycomm, ierr)
1382) start(1) = istart
1383) length(1) = local_size
1384)
1385) stride = 1
1386) call h5sselect_hyperslab_f(file_space_id,H5S_SELECT_SET_F,start,length, &
1387) hdf5_err,stride,stride)
1388) ! write the data
1389) call h5pcreate_f(H5P_DATASET_XFER_F,prop_id,hdf5_err)
1390) #ifndef SERIAL_HDF5
1391) call h5pset_dxpl_mpio_f(prop_id,H5FD_MPIO_INDEPENDENT_F, &
1392) hdf5_err)
1393) #endif
1394)
1395) call PetscLogEventBegin(logging%event_h5dwrite_f,ierr);CHKERRQ(ierr)
1396) call h5dwrite_f(data_set_id,H5T_NATIVE_DOUBLE,vec_x_ptr,dims, &
1397) hdf5_err,memory_space_id,file_space_id,prop_id)
1398) call PetscLogEventEnd(logging%event_h5dwrite_f,ierr);CHKERRQ(ierr)
1399)
1400) call h5pclose_f(prop_id,hdf5_err)
1401)
1402) call h5dclose_f(data_set_id,hdf5_err)
1403) call h5sclose_f(file_space_id,hdf5_err)
1404)
1405) ! YC
1406) ! memory space which is a 1D vector
1407) rank_mpi = 1
1408) dims = 0
1409) dims(1) = local_size
1410) call h5screate_simple_f(rank_mpi,dims,memory_space_id,hdf5_err,dims)
1411)
1412) ! file space which is a 2D block
1413) rank_mpi = 1
1414) dims = 0
1415) dims(1) = surf_grid%nmax
1416) call h5pcreate_f(H5P_DATASET_CREATE_F,prop_id,hdf5_err)
1417)
1418) string = "YC" // CHAR(0)
1419)
1420) call h5eset_auto_f(OFF,hdf5_err)
1421) call h5dopen_f(file_id,string,data_set_id,hdf5_err)
1422) hdf5_flag = hdf5_err
1423) call h5eset_auto_f(ON,hdf5_err)
1424) if (hdf5_flag < 0) then
1425) call h5screate_simple_f(rank_mpi,dims,file_space_id,hdf5_err,dims)
1426) call h5dcreate_f(file_id,string,H5T_NATIVE_DOUBLE,file_space_id, &
1427) data_set_id,hdf5_err,prop_id)
1428) else
1429) call h5dget_space_f(data_set_id,file_space_id,hdf5_err)
1430) endif
1431)
1432) call h5pclose_f(prop_id,hdf5_err)
1433)
1434) istart = 0
1435) call MPI_Exscan(local_size, istart, ONE_INTEGER_MPI, &
1436) MPIU_INTEGER, MPI_SUM, option%mycomm, ierr)
1437) start(1) = istart
1438) length(1) = local_size
1439)
1440) stride = 1
1441) call h5sselect_hyperslab_f(file_space_id,H5S_SELECT_SET_F,start,length, &
1442) hdf5_err,stride,stride)
1443) ! write the data
1444) call h5pcreate_f(H5P_DATASET_XFER_F,prop_id,hdf5_err)
1445) #ifndef SERIAL_HDF5
1446) call h5pset_dxpl_mpio_f(prop_id,H5FD_MPIO_INDEPENDENT_F, &
1447) hdf5_err)
1448) #endif
1449)
1450) call PetscLogEventBegin(logging%event_h5dwrite_f,ierr);CHKERRQ(ierr)
1451) call h5dwrite_f(data_set_id,H5T_NATIVE_DOUBLE,vec_y_ptr,dims, &
1452) hdf5_err,memory_space_id,file_space_id,prop_id)
1453) call PetscLogEventEnd(logging%event_h5dwrite_f,ierr);CHKERRQ(ierr)
1454)
1455) call h5pclose_f(prop_id,hdf5_err)
1456)
1457) call h5dclose_f(data_set_id,hdf5_err)
1458) call h5sclose_f(file_space_id,hdf5_err)
1459)
1460) ! ZC
1461) ! memory space which is a 1D vector
1462) rank_mpi = 1
1463) dims = 0
1464) dims(1) = local_size
1465) call h5screate_simple_f(rank_mpi,dims,memory_space_id,hdf5_err,dims)
1466)
1467) ! file space which is a 2D block
1468) rank_mpi = 1
1469) dims = 0
1470) dims(1) = surf_grid%nmax
1471) call h5pcreate_f(H5P_DATASET_CREATE_F,prop_id,hdf5_err)
1472)
1473) string = "ZC" // CHAR(0)
1474)
1475) call h5eset_auto_f(OFF,hdf5_err)
1476) call h5dopen_f(file_id,string,data_set_id,hdf5_err)
1477) hdf5_flag = hdf5_err
1478) call h5eset_auto_f(ON,hdf5_err)
1479) if (hdf5_flag < 0) then
1480) call h5screate_simple_f(rank_mpi,dims,file_space_id,hdf5_err,dims)
1481) call h5dcreate_f(file_id,string,H5T_NATIVE_DOUBLE,file_space_id, &
1482) data_set_id,hdf5_err,prop_id)
1483) else
1484) call h5dget_space_f(data_set_id,file_space_id,hdf5_err)
1485) endif
1486)
1487) call h5pclose_f(prop_id,hdf5_err)
1488)
1489) istart = 0
1490) call MPI_Exscan(local_size, istart, ONE_INTEGER_MPI, &
1491) MPIU_INTEGER, MPI_SUM, option%mycomm, ierr)
1492) start(1) = istart
1493) length(1) = local_size
1494)
1495) stride = 1
1496) call h5sselect_hyperslab_f(file_space_id,H5S_SELECT_SET_F,start,length, &
1497) hdf5_err,stride,stride)
1498) ! write the data
1499) call h5pcreate_f(H5P_DATASET_XFER_F,prop_id,hdf5_err)
1500) #ifndef SERIAL_HDF5
1501) call h5pset_dxpl_mpio_f(prop_id,H5FD_MPIO_INDEPENDENT_F, &
1502) hdf5_err)
1503) #endif
1504)
1505) call PetscLogEventBegin(logging%event_h5dwrite_f,ierr);CHKERRQ(ierr)
1506) call h5dwrite_f(data_set_id,H5T_NATIVE_DOUBLE,vec_z_ptr,dims, &
1507) hdf5_err,memory_space_id,file_space_id,prop_id)
1508) call PetscLogEventEnd(logging%event_h5dwrite_f,ierr);CHKERRQ(ierr)
1509)
1510) call h5pclose_f(prop_id,hdf5_err)
1511)
1512) call h5dclose_f(data_set_id,hdf5_err)
1513) call h5sclose_f(file_space_id,hdf5_err)
1514)
1515)
1516) call VecRestoreArrayF90(natural_x_cell_vec,vec_x_ptr,ierr);CHKERRQ(ierr)
1517) call VecRestoreArrayF90(natural_y_cell_vec,vec_y_ptr,ierr);CHKERRQ(ierr)
1518) call VecRestoreArrayF90(natural_z_cell_vec,vec_z_ptr,ierr);CHKERRQ(ierr)
1519)
1520) call VecDestroy(global_x_cell_vec,ierr);CHKERRQ(ierr)
1521) call VecDestroy(global_y_cell_vec,ierr);CHKERRQ(ierr)
1522) call VecDestroy(global_z_cell_vec,ierr);CHKERRQ(ierr)
1523)
1524) call VecDestroy(natural_x_cell_vec,ierr);CHKERRQ(ierr)
1525) call VecDestroy(natural_y_cell_vec,ierr);CHKERRQ(ierr)
1526) call VecDestroy(natural_z_cell_vec,ierr);CHKERRQ(ierr)
1527)
1528) #endif
1529) ! ifdef SCORPIO_WRITE
1530)
1531) end subroutine WriteHDF5CoordinatesUGridXDMF
1532) #endif
1533)
1534) ! ************************************************************************** !
1535)
1536) subroutine OutputSurfaceGetVarFromArray(surf_realization,vec,ivar,isubvar,isubvar1)
1537) !
1538) ! ifdef PETSC_HAVE_HDF5
1539) ! This routine extracts variables indexed by ivar from a multivar array
1540) !
1541) ! Author: Gautam Bisht, LBNL
1542) ! Date: 01/30/13
1543) !
1544)
1545) use Realization_Surface_class, only : realization_surface_type, &
1546) RealizSurfGetVariable
1547) use Grid_module
1548) use Option_module
1549) use Field_module
1550)
1551) implicit none
1552)
1553) #include "petsc/finclude/petscvec.h"
1554) #include "petsc/finclude/petscvec.h90"
1555) #include "petsc/finclude/petsclog.h"
1556)
1557) class(realization_surface_type) :: surf_realization
1558) Vec :: vec
1559) PetscInt :: ivar
1560) PetscInt :: isubvar
1561) PetscInt, optional :: isubvar1
1562)
1563) PetscErrorCode :: ierr
1564)
1565) call PetscLogEventBegin(logging%event_output_get_var_from_array, &
1566) ierr);CHKERRQ(ierr)
1567)
1568) call RealizSurfGetVariable(surf_realization,vec,ivar,isubvar,isubvar1)
1569)
1570) call PetscLogEventEnd(logging%event_output_get_var_from_array, &
1571) ierr);CHKERRQ(ierr)
1572)
1573) end subroutine OutputSurfaceGetVarFromArray
1574)
1575) ! ************************************************************************** !
1576)
1577) subroutine OutputSurfaceAvegVars(surf_realization,realization)
1578) !
1579) ! This routine to temporally average variables and output them.
1580) !
1581) ! Author: Gautam Bisht, LBNL
1582) ! Date: 03/22/13
1583) !
1584)
1585) use Realization_Surface_class, only : realization_surface_type
1586) use Realization_Subsurface_class, only : realization_subsurface_type
1587) use Option_module, only : OptionCheckTouch, option_type, printMsg
1588) use Output_Aux_module
1589) use Output_Common_module, only : OutputGetVarFromArray
1590) use Surface_Field_module
1591)
1592) implicit none
1593)
1594) class(realization_surface_type) :: surf_realization
1595) class(realization_subsurface_type) :: realization
1596)
1597) type(option_type), pointer :: option
1598) type(output_option_type), pointer :: output_option
1599) type(output_variable_type), pointer :: cur_variable
1600) type(surface_field_type), pointer :: surf_field
1601)
1602) PetscReal :: dtime
1603) PetscBool :: aveg_plot_flag
1604) PetscInt :: ivar
1605) PetscReal,pointer :: aval_p(:),ival_p(:)
1606) PetscErrorCode :: ierr
1607) PetscLogDouble :: tstart, tend
1608)
1609) option => surf_realization%option
1610) output_option => surf_realization%output_option
1611) surf_field => surf_realization%surf_field
1612)
1613) ! Return back if called at the beginning of simulation
1614) if (option%time<1.d-10) return
1615)
1616) dtime = option%time-output_option%aveg_var_time
1617) output_option%aveg_var_dtime = output_option%aveg_var_dtime + dtime
1618) output_option%aveg_var_time = output_option%aveg_var_time + dtime
1619)
1620) if (abs(output_option%aveg_var_dtime - &
1621) output_option%periodic_snap_output_time_incr)<1.d0) then
1622) aveg_plot_flag=PETSC_TRUE
1623) else
1624) aveg_plot_flag=PETSC_FALSE
1625) endif
1626)
1627) if (.not.associated(output_option%aveg_output_variable_list%first)) then
1628) if (output_option%print_hdf5_aveg_mass_flowrate.or. &
1629) output_option%print_hdf5_aveg_energy_flowrate) then
1630) ! There is a possibility to output average-flowrates, thus
1631) ! call output subroutine depending on mesh type
1632) if (surf_realization%discretization%itype == UNSTRUCTURED_GRID) then
1633) call OutputSurfaceHDF5UGridXDMF(surf_realization,realization,AVERAGED_VARS)
1634) else
1635) ! call OutputHDF5(realization_base,AVERAGED_VARS)
1636) endif
1637) endif
1638) return
1639) endif
1640)
1641) ivar = 0
1642) cur_variable => output_option%aveg_output_variable_list%first
1643) do
1644) if (.not.associated(cur_variable)) exit
1645)
1646) ! Get the variable
1647) call OutputGetVarFromArray(surf_realization,surf_field%work, &
1648) cur_variable%ivar, &
1649) cur_variable%isubvar)
1650)
1651) ! Cumulatively add the variable*dtime
1652) ivar = ivar + 1
1653) call VecGetArrayF90(surf_field%work,ival_p,ierr);CHKERRQ(ierr)
1654) call VecGetArrayF90(surf_field%avg_vars_vec(ivar),aval_p, &
1655) ierr);CHKERRQ(ierr)
1656) aval_p = aval_p + ival_p*dtime
1657) call VecRestoreArrayF90(surf_field%work,ival_p,ierr);CHKERRQ(ierr)
1658) call VecRestoreArrayF90(surf_field%avg_vars_vec(ivar),aval_p, &
1659) ierr);CHKERRQ(ierr)
1660)
1661) ! Check if it is time to output the temporally average variable
1662) if (aveg_plot_flag) then
1663)
1664) ! Divide vector values by 'time'
1665) call VecGetArrayF90(surf_field%avg_vars_vec(ivar),aval_p, &
1666) ierr);CHKERRQ(ierr)
1667) aval_p = aval_p/output_option%periodic_snap_output_time_incr
1668) call VecRestoreArrayF90(surf_field%avg_vars_vec(ivar),aval_p, &
1669) ierr);CHKERRQ(ierr)
1670)
1671) endif
1672)
1673) cur_variable => cur_variable%next
1674) enddo
1675)
1676) if (aveg_plot_flag) then
1677)
1678) if (surf_realization%output_option%print_hdf5) then
1679) call PetscTime(tstart,ierr);CHKERRQ(ierr)
1680) call PetscLogEventBegin(logging%event_output_hdf5,ierr);CHKERRQ(ierr)
1681) if (surf_realization%discretization%itype == UNSTRUCTURED_GRID) then
1682) call OutputSurfaceHDF5UGridXDMF(surf_realization,realization,AVERAGED_VARS)
1683) else
1684) ! call OutputHDF5(surf_realization,AVERAGED_VARS)
1685) endif
1686) call PetscLogEventEnd(logging%event_output_hdf5,ierr);CHKERRQ(ierr)
1687) call PetscTime(tend,ierr);CHKERRQ(ierr)
1688) write(option%io_buffer,'(f10.2," Seconds to write HDF5 file.")') tend-tstart
1689) call printMsg(option)
1690) endif
1691)
1692) ! Reset the vectors to zero
1693) do ivar=1,output_option%aveg_output_variable_list%nvars
1694) call VecSet(surf_field%avg_vars_vec(ivar),0.d0,ierr);CHKERRQ(ierr)
1695) enddo
1696)
1697) output_option%aveg_var_dtime=0.d0
1698)
1699) endif
1700)
1701) end subroutine OutputSurfaceAvegVars
1702)
1703) ! ************************************************************************** !
1704)
1705) subroutine OutputSurfaceVariableRead(input,option,output_variable_list)
1706) !
1707) ! This routine reads variable from input file.
1708) !
1709) ! Author: Gautam Bisht, LBNL
1710) ! Date: 03/22/13
1711) !
1712)
1713) use Option_module
1714) use Input_Aux_module
1715) use String_module
1716) use Variables_module
1717)
1718) implicit none
1719)
1720) type(option_type) :: option
1721) type(input_type), pointer :: input
1722) type(output_variable_list_type), pointer :: output_variable_list
1723)
1724) character(len=MAXWORDLENGTH) :: word
1725) character(len=MAXWORDLENGTH) :: name, units
1726) type(output_variable_type), pointer :: output_variable
1727)
1728) do
1729) call InputReadPflotranString(input,option)
1730) if (InputError(input)) exit
1731) if (InputCheckExit(input,option)) exit
1732)
1733) call InputReadWord(input,option,word,PETSC_TRUE)
1734) call InputErrorMsg(input,option,'keyword','VARIABLES')
1735) call StringToUpper(word)
1736)
1737) select case(trim(word))
1738) case ('SURFACE_LIQUID_HEAD')
1739) name = 'H'
1740) units = 'm'
1741) call OutputVariableAddToList(output_variable_list,name,OUTPUT_GENERIC, &
1742) units,SURFACE_LIQUID_HEAD)
1743)
1744) case ('TEMPERATURE')
1745) name = 'Temperature'
1746) units = 'C'
1747) call OutputVariableAddToList(output_variable_list,name,OUTPUT_GENERIC, &
1748) units,TEMPERATURE)
1749) case ('PROCESS_ID')
1750) units = ''
1751) name = 'Process ID'
1752) output_variable => OutputVariableCreate(name,OUTPUT_DISCRETE, &
1753) units,PROCESS_ID)
1754) output_variable%plot_only = PETSC_TRUE ! toggle output off for observation
1755) output_variable%iformat = 1 ! integer
1756) call OutputVariableAddToList(output_variable_list,output_variable)
1757) case default
1758) call InputKeywordUnrecognized(word,'SURFACE,VARIABLES',option)
1759) end select
1760) enddo
1761)
1762) end subroutine OutputSurfaceVariableRead
1763)
1764) ! ************************************************************************** !
1765)
1766) function OutputSurfaceHDF5FilenameID(output_option,option,var_list_type)
1767) !
1768) ! This subroutine creates an ID for HDF5 filename for:
1769) ! - Instantaneous, or
1770) ! - Temporally averaged variables.
1771) !
1772) ! Author: Gautam Bisht, LBNL
1773) ! Date: 03/22/13
1774) !
1775)
1776) use Option_module
1777)
1778) implicit none
1779)
1780) type(option_type) :: option
1781) type(output_option_type) :: output_option
1782) PetscInt :: var_list_type
1783)
1784) character(len=MAXWORDLENGTH) :: OutputSurfaceHDF5FilenameID
1785) PetscInt :: file_number
1786)
1787) select case(var_list_type)
1788) case (INSTANTANEOUS_VARS)
1789) file_number = floor(real(output_option%plot_number)/ &
1790) output_option%times_per_h5_file)
1791) case (AVERAGED_VARS)
1792) file_number = floor((option%time - &
1793) output_option%periodic_snap_output_time_incr)/ &
1794) output_option%periodic_snap_output_time_incr/ &
1795) output_option%times_per_h5_file)
1796) end select
1797)
1798) if (file_number < 10) then
1799) write(OutputSurfaceHDF5FilenameID,'("00",i1)') file_number
1800) else if (output_option%plot_number < 100) then
1801) write(OutputSurfaceHDF5FilenameID,'("0",i2)') file_number
1802) else if (output_option%plot_number < 1000) then
1803) write(OutputSurfaceHDF5FilenameID,'(i3)') file_number
1804) else if (output_option%plot_number < 10000) then
1805) write(OutputSurfaceHDF5FilenameID,'(i4)') file_number
1806) endif
1807)
1808) OutputSurfaceHDF5FilenameID = adjustl(OutputSurfaceHDF5FilenameID)
1809)
1810) end function OutputSurfaceHDF5FilenameID
1811)
1812) #if defined(PETSC_HAVE_HDF5)
1813)
1814) ! ************************************************************************** !
1815)
1816) subroutine OutputSurfaceGetFlowrates(surf_realization)
1817) !
1818) ! This returns mass/energy flowrate at all faces of a control volume for
1819) ! surface realizaton.
1820) !
1821) ! Author: Gautam Bisht, LBNL
1822) ! Date: 03/25/2013
1823) !
1824)
1825) use hdf5
1826) use HDF5_module
1827) use Realization_Surface_class, only : realization_surface_type
1828) use Patch_module
1829) use Grid_module
1830) use Option_module
1831) use Grid_Unstructured_Aux_module
1832) use Grid_Unstructured_Cell_module
1833) use Variables_module
1834) use Connection_module
1835) use Coupler_module
1836) use HDF5_Aux_module
1837) use Output_Aux_module
1838) use Surface_Field_module
1839)
1840) implicit none
1841)
1842) #include "petsc/finclude/petscvec.h"
1843) #include "petsc/finclude/petscvec.h90"
1844) #include "petsc/finclude/petsclog.h"
1845) #include "petsc/finclude/petscsys.h"
1846)
1847) class(realization_surface_type) :: surf_realization
1848) type(option_type), pointer :: option
1849)
1850) type(patch_type), pointer :: patch
1851) type(grid_type), pointer :: grid
1852) type(grid_unstructured_type),pointer :: ugrid
1853) type(connection_set_list_type), pointer :: connection_set_list
1854) type(connection_set_type), pointer :: cur_connection_set
1855) type(coupler_type), pointer :: boundary_condition
1856) type(ugdm_type),pointer :: ugdm
1857) type(output_option_type), pointer :: output_option
1858) type(surface_field_type), pointer :: surf_field
1859)
1860) PetscInt :: local_id
1861) PetscInt :: ghosted_id
1862) PetscInt :: idual
1863) PetscInt :: iconn
1864) PetscInt :: face_id
1865) PetscInt :: local_id_up,local_id_dn
1866) PetscInt :: ghosted_id_up,ghosted_id_dn
1867) PetscInt :: iface_up,iface_dn
1868) PetscInt :: dof
1869) PetscInt :: sum_connection
1870) PetscInt :: offset
1871) PetscInt :: cell_type
1872) PetscInt :: local_size
1873) PetscInt :: i
1874) PetscInt :: iface
1875) PetscInt :: ndof
1876)
1877) PetscReal, pointer :: flowrates(:,:,:)
1878) PetscReal, pointer :: vec_ptr(:)
1879) PetscReal, pointer :: vec_ptr2(:)
1880) PetscReal, pointer :: vec_ptr3(:)
1881) PetscReal, pointer :: double_array(:)
1882) PetscReal :: dtime
1883)
1884) Vec :: natural_flowrates_vec
1885)
1886) PetscMPIInt :: hdf5_err
1887) PetscErrorCode :: ierr
1888)
1889) character(len=MAXSTRINGLENGTH) :: string
1890) character(len=MAXWORDLENGTH) :: unit_string
1891)
1892) patch => surf_realization%patch
1893) grid => patch%grid
1894) ugrid => grid%unstructured_grid
1895) output_option =>surf_realization%output_option
1896) option => surf_realization%option
1897) surf_field => surf_realization%surf_field
1898)
1899) ! Create UGDM for
1900) call UGridCreateUGDM(grid%unstructured_grid,ugdm, &
1901) (option%nflowdof*MAX_FACE_PER_CELL_SURF + 1),option)
1902)
1903) ! Create a flowrate vector in natural order
1904) call UGridDMCreateVector(grid%unstructured_grid,ugdm,natural_flowrates_vec, &
1905) NATURAL,option)
1906)
1907)
1908) allocate(flowrates(option%nflowdof,MAX_FACE_PER_CELL_SURF,ugrid%nlmax))
1909) flowrates = 0.d0
1910) call VecGetArrayF90(surf_field%flowrate_inst,vec_ptr,ierr);CHKERRQ(ierr)
1911) vec_ptr = 0.d0
1912)
1913) offset = 1+option%nflowdof*MAX_FACE_PER_CELL_SURF
1914) ! Save the number of faces of all cell
1915) do local_id = 1,grid%nlmax
1916) ghosted_id = grid%nL2G(local_id)
1917) cell_type = ugrid%cell_type(ghosted_id)
1918) vec_ptr((local_id-1)*offset+1) = UCellGetNFaces(cell_type,option)
1919) enddo
1920)
1921) ! Interior Flowrates Terms -----------------------------------
1922) connection_set_list => grid%internal_connection_set_list
1923) cur_connection_set => connection_set_list%first
1924) sum_connection = 0
1925) do
1926) if (.not.associated(cur_connection_set)) exit
1927) do iconn = 1, cur_connection_set%num_connections
1928) sum_connection = sum_connection + 1
1929) face_id = cur_connection_set%face_id(iconn)
1930) ghosted_id_up = cur_connection_set%id_up(iconn)
1931) ghosted_id_dn = cur_connection_set%id_dn(iconn)
1932) local_id_up = grid%nG2L(ghosted_id_up)
1933) local_id_dn = grid%nG2L(ghosted_id_dn)
1934) do iface_up = 1,MAX_FACE_PER_CELL_SURF
1935) if (face_id==ugrid%cell_to_face_ghosted(iface_up,local_id_up)) exit
1936) enddo
1937) iface_dn=-1
1938) if (local_id_dn>0) then
1939) do iface_dn = 1,MAX_FACE_PER_CELL_SURF
1940) if (face_id==ugrid%cell_to_face_ghosted(iface_dn,local_id_dn)) exit
1941) enddo
1942) endif
1943)
1944) do dof=1,option%nflowdof
1945) ! Save flowrate for iface_up of local_id_up cell using flowrate up-->dn
1946) flowrates(dof,iface_up,local_id_up) = patch%internal_flow_fluxes(dof,sum_connection)
1947) vec_ptr((local_id_up-1)*offset + (dof-1)*MAX_FACE_PER_CELL_SURF + iface_up + 1) = &
1948) patch%internal_flow_fluxes(dof,sum_connection)
1949)
1950) if (iface_dn>0) then
1951) ! Save flowrate for iface_dn of local_id_dn cell using -ve flowrate up-->dn
1952) flowrates(dof,iface_dn,local_id_dn) = -patch%internal_flow_fluxes(dof,sum_connection)
1953) vec_ptr((local_id_dn-1)*offset + (dof-1)*MAX_FACE_PER_CELL_SURF + iface_dn + 1) = &
1954) -patch%internal_flow_fluxes(dof,sum_connection)
1955) endif
1956) enddo
1957)
1958) enddo
1959) cur_connection_set => cur_connection_set%next
1960) enddo
1961)
1962) ! Boundary Flowrates Terms -----------------------------------
1963) boundary_condition => patch%boundary_condition_list%first
1964) sum_connection = 0
1965) do
1966) if (.not.associated(boundary_condition)) exit
1967)
1968) cur_connection_set => boundary_condition%connection_set
1969) sum_connection = 0
1970)
1971) do iconn = 1, cur_connection_set%num_connections
1972) sum_connection = sum_connection + 1
1973) face_id = cur_connection_set%face_id(iconn)
1974) ghosted_id_dn = cur_connection_set%id_dn(iconn)
1975) local_id_dn = grid%nG2L(ghosted_id_dn)
1976) do iface_dn = 1,MAX_FACE_PER_CELL_SURF
1977) if (face_id==ugrid%cell_to_face_ghosted(iface_dn,local_id_dn)) exit
1978) enddo
1979)
1980) do dof=1,option%nflowdof
1981) ! Save flowrate for iface_dn of local_id_dn cell using -ve flowrate up-->dn
1982) flowrates(dof,iface_dn,local_id_dn) = -patch%boundary_flow_fluxes(dof,sum_connection)
1983) vec_ptr((local_id_dn-1)*offset + (dof-1)*MAX_FACE_PER_CELL_SURF + iface_dn + 1) = &
1984) -patch%boundary_flow_fluxes(dof,sum_connection)
1985) enddo
1986) enddo
1987) boundary_condition => boundary_condition%next
1988) enddo
1989)
1990) deallocate(flowrates)
1991)
1992) call VecRestoreArrayF90(surf_field%flowrate_inst,vec_ptr,ierr);CHKERRQ(ierr)
1993)
1994) ! Scatter flowrate from Global --> Natural order
1995) call VecScatterBegin(ugdm%scatter_gton,surf_field%flowrate_inst,natural_flowrates_vec, &
1996) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1997) call VecScatterEnd(ugdm%scatter_gton,surf_field%flowrate_inst,natural_flowrates_vec, &
1998) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1999)
2000) call VecGetArrayF90(natural_flowrates_vec,vec_ptr,ierr);CHKERRQ(ierr)
2001) call VecGetArrayF90(surf_field%flowrate_inst,vec_ptr2,ierr);CHKERRQ(ierr)
2002)
2003) ! Copy the vectors
2004) vec_ptr2 = vec_ptr
2005)
2006) if (output_option%print_hdf5_aveg_mass_flowrate.or. &
2007) output_option%print_hdf5_aveg_energy_flowrate) then
2008)
2009) dtime = option%time-output_option%aveg_var_time
2010) call VecGetArrayF90(surf_field%flowrate_aveg,vec_ptr3,ierr);CHKERRQ(ierr)
2011) vec_ptr3 = vec_ptr3 + vec_ptr2/dtime
2012) call VecRestoreArrayF90(surf_field%flowrate_aveg,vec_ptr3, &
2013) ierr);CHKERRQ(ierr)
2014)
2015) endif
2016)
2017) call VecRestoreArrayF90(natural_flowrates_vec,vec_ptr,ierr);CHKERRQ(ierr)
2018) call VecRestoreArrayF90(surf_field%flowrate_inst,vec_ptr2, &
2019) ierr);CHKERRQ(ierr)
2020)
2021) call VecDestroy(natural_flowrates_vec,ierr);CHKERRQ(ierr)
2022) call UGridDMDestroy(ugdm)
2023)
2024) end subroutine OutputSurfaceGetFlowrates
2025)
2026) ! ************************************************************************** !
2027)
2028) subroutine WriteHDF5SurfaceFlowratesUGrid(surf_realization,file_id,var_list_type)
2029) !
2030) ! This routine writes (mass/energy) flowrate for unstructured grid for
2031) ! surface realization.
2032) !
2033) ! Author: Gautam Bisht, LBNL
2034) ! Date: 03/25/2013
2035) !
2036)
2037) use hdf5
2038) use HDF5_module
2039) use Realization_Surface_class, only : realization_surface_type
2040) use Patch_module
2041) use Grid_module
2042) use Option_module
2043) use Grid_Unstructured_Aux_module
2044) use Grid_Unstructured_Cell_module
2045) use Variables_module
2046) use Connection_module
2047) use Coupler_module
2048) use HDF5_Aux_module
2049) use Output_Aux_module
2050) use Surface_Field_module
2051)
2052) implicit none
2053)
2054) #include "petsc/finclude/petscvec.h"
2055) #include "petsc/finclude/petscvec.h90"
2056) #include "petsc/finclude/petsclog.h"
2057) #include "petsc/finclude/petscsys.h"
2058)
2059) class(realization_surface_type) :: surf_realization
2060) type(option_type), pointer :: option
2061) PetscInt :: var_list_type
2062)
2063) #if defined(SCORPIO_WRITE)
2064) integer :: file_id
2065) integer :: data_type
2066) integer :: grp_id
2067) integer :: file_space_id
2068) integer :: memory_space_id
2069) integer :: data_set_id
2070) integer :: realization_set_id
2071) integer :: prop_id
2072) integer :: dims(3)
2073) integer :: start(3), length(3), stride(3)
2074) integer :: rank_mpi,file_space_rank_mpi
2075) integer :: hdf5_flag
2076) integer, parameter :: ON=1, OFF=0
2077) #else
2078) integer(HID_T) :: file_id
2079) integer(HID_T) :: data_type
2080) integer(HID_T) :: grp_id
2081) integer(HID_T) :: file_space_id
2082) integer(HID_T) :: realization_set_id
2083) integer(HID_T) :: memory_space_id
2084) integer(HID_T) :: data_set_id
2085) integer(HID_T) :: prop_id
2086) integer(HSIZE_T) :: dims(3)
2087) integer(HSIZE_T) :: start(3), length(3), stride(3)
2088) PetscMPIInt :: rank_mpi,file_space_rank_mpi
2089) PetscMPIInt :: hdf5_flag
2090) PetscMPIInt, parameter :: ON=1, OFF=0
2091) #endif
2092)
2093) type(patch_type), pointer :: patch
2094) type(grid_type), pointer :: grid
2095) type(grid_unstructured_type),pointer :: ugrid
2096) type(connection_set_list_type), pointer :: connection_set_list
2097) type(connection_set_type), pointer :: cur_connection_set
2098) type(coupler_type), pointer :: boundary_condition
2099) type(ugdm_type),pointer :: ugdm
2100) type(output_option_type), pointer :: output_option
2101) type(surface_field_type), pointer :: surf_field
2102)
2103) PetscInt :: local_id
2104) PetscInt :: ghosted_id
2105) PetscInt :: idual
2106) PetscInt :: iconn
2107) PetscInt :: face_id
2108) PetscInt :: local_id_up,local_id_dn
2109) PetscInt :: ghosted_id_up,ghosted_id_dn
2110) PetscInt :: istart
2111) PetscInt :: iface_up,iface_dn
2112) PetscInt :: dof
2113) PetscInt :: sum_connection
2114) PetscInt :: offset
2115) PetscInt :: cell_type
2116) PetscInt :: local_size
2117) PetscInt :: i
2118) PetscInt :: iface
2119) PetscInt :: ndof
2120)
2121) PetscReal, pointer :: flowrates(:,:,:)
2122) PetscReal, pointer :: vec_ptr1(:)
2123) PetscReal, pointer :: vec_ptr2(:)
2124) PetscReal, pointer :: double_array(:)
2125) PetscReal :: dtime
2126)
2127) PetscBool :: mass_flowrate
2128) PetscBool :: energy_flowrate
2129)
2130) Vec :: global_flowrates_vec
2131) Vec :: natural_flowrates_vec
2132)
2133) PetscMPIInt :: hdf5_err
2134) PetscErrorCode :: ierr
2135)
2136) character(len=MAXSTRINGLENGTH) :: string
2137) character(len=MAXWORDLENGTH) :: unit_string
2138)
2139) patch => surf_realization%patch
2140) grid => patch%grid
2141) ugrid => grid%unstructured_grid
2142) option => surf_realization%option
2143) output_option =>surf_realization%output_option
2144) surf_field => surf_realization%surf_field
2145)
2146) #if defined(SCORPIO_WRITE)
2147) write(*,*) 'SCORPIO_WRITE'
2148) option%io_buffer = 'WriteHDF5FlowratesUGrid not supported for SCORPIO_WRITE'
2149) call printErrMsg(option)
2150) #else
2151) select case(option%iflowmode)
2152) case (RICHARDS_MODE)
2153) ndof=1
2154) case (TH_MODE)
2155) ndof=1
2156) if (output_option%print_hdf5_mass_flowrate.and.output_option%print_hdf5_energy_flowrate) ndof = 2
2157) case default
2158) option%io_buffer='FLOWRATE output not supported in this mode'
2159) call printErrMsg(option)
2160) end select
2161)
2162) call VecGetLocalSize(surf_field%flowrate_inst,local_size,ierr);CHKERRQ(ierr)
2163) local_size = local_size/(MAX_FACE_PER_CELL_SURF+1)/option%nflowdof
2164)
2165) allocate(double_array(local_size*(MAX_FACE_PER_CELL_SURF+1)))
2166) double_array = 0.d0
2167)
2168) offset = option%nflowdof*MAX_FACE_PER_CELL_SURF+1
2169)
2170) select case (var_list_type)
2171) case (INSTANTANEOUS_VARS)
2172) call VecGetArrayF90(surf_field%flowrate_inst,vec_ptr1, &
2173) ierr);CHKERRQ(ierr)
2174) mass_flowrate = output_option%print_hdf5_mass_flowrate
2175) energy_flowrate = output_option%print_hdf5_energy_flowrate
2176) case (AVERAGED_VARS)
2177) call VecGetArrayF90(surf_field%flowrate_inst,vec_ptr1, &
2178) ierr);CHKERRQ(ierr)
2179) call VecGetArrayF90(surf_field%flowrate_aveg,vec_ptr2, &
2180) ierr);CHKERRQ(ierr)
2181) mass_flowrate = output_option%print_hdf5_aveg_mass_flowrate
2182) energy_flowrate = output_option%print_hdf5_aveg_energy_flowrate
2183) end select
2184)
2185)
2186) do dof = 1,option%nflowdof
2187)
2188) if (dof==1 .and. (.not.mass_flowrate)) exit
2189) if (dof==2 .and. (.not.energy_flowrate)) exit
2190)
2191) select case(option%iflowmode)
2192) case(RICHARDS_MODE)
2193) string = "Mass_Flowrate [kg_s]" // CHAR(0)
2194) case(TH_MODE)
2195) if (dof==1) then
2196) string = "Mass_Flowrate [kg_s]" // CHAR(0)
2197) else
2198) string = "Energy_Flowrate [J_s]" // CHAR(0)
2199) endif
2200) case default
2201) option%io_buffer='FLOWRATE output not implemented in this mode.'
2202) call printErrMsg(option)
2203) end select
2204)
2205) if (var_list_type==AVERAGED_VARS) string = 'Aveg_' // trim(string) // CHAR(0)
2206)
2207) ! memory space which is a 1D vector
2208) rank_mpi = 1
2209) dims = 0
2210) dims(1) = local_size*(MAX_FACE_PER_CELL_SURF+1)
2211) call h5screate_simple_f(rank_mpi,dims,memory_space_id,hdf5_err,dims)
2212)
2213) ! file space which is a 2D block
2214) rank_mpi = 2
2215) dims = 0
2216) dims(2) = ugrid%nmax
2217) dims(1) = MAX_FACE_PER_CELL_SURF + 1
2218) call h5pcreate_f(H5P_DATASET_CREATE_F,prop_id,hdf5_err)
2219)
2220) call h5eset_auto_f(OFF,hdf5_err)
2221) call h5dopen_f(file_id,trim(string),data_set_id,hdf5_err)
2222) hdf5_flag = hdf5_err
2223) call h5eset_auto_f(ON,hdf5_err)
2224) if (hdf5_flag < 0) then
2225) call h5screate_simple_f(rank_mpi,dims,file_space_id,hdf5_err,dims)
2226) call h5dcreate_f(file_id,trim(string),H5T_NATIVE_DOUBLE,file_space_id, &
2227) data_set_id,hdf5_err,prop_id)
2228) else
2229) call h5dget_space_f(data_set_id,file_space_id,hdf5_err)
2230) endif
2231)
2232) call h5pclose_f(prop_id,hdf5_err)
2233)
2234) istart = 0
2235) call MPI_Exscan(local_size, istart, ONE_INTEGER_MPI, &
2236) MPIU_INTEGER, MPI_SUM, option%mycomm, ierr)
2237)
2238) start(2) = istart
2239) start(1) = 0
2240)
2241) length(2) = local_size
2242) length(1) = MAX_FACE_PER_CELL_SURF + 1
2243)
2244) stride = 1
2245) call h5sselect_hyperslab_f(file_space_id,H5S_SELECT_SET_F,start,length, &
2246) hdf5_err,stride,stride)
2247) ! write the data
2248) call h5pcreate_f(H5P_DATASET_XFER_F,prop_id,hdf5_err)
2249) #ifndef SERIAL_HDF5
2250) call h5pset_dxpl_mpio_f(prop_id,H5FD_MPIO_INDEPENDENT_F, &
2251) hdf5_err)
2252) #endif
2253)
2254) select case (var_list_type)
2255) case (INSTANTANEOUS_VARS)
2256)
2257) do i=1,local_size
2258) ! Num. of faces for each cell (Note: Use vec_ptr1 not vec_ptr2)
2259) double_array((i-1)*(MAX_FACE_PER_CELL_SURF+1)+1) = &
2260) vec_ptr1((i-1)*offset+(dof-1)*MAX_FACE_PER_CELL_SURF+1)
2261) ! Flowrate values for each face
2262) do iface = 1,MAX_FACE_PER_CELL_SURF
2263) double_array((i-1)*(MAX_FACE_PER_CELL_SURF+1)+iface+1) = &
2264) vec_ptr1((i-1)*offset + (dof-1)*MAX_FACE_PER_CELL_SURF + iface + 1)
2265) enddo
2266) enddo
2267)
2268) case (AVERAGED_VARS)
2269)
2270) do i=1,local_size
2271) ! Num. of faces for each cell (Note: Use vec_ptr1 not vec_ptr2)
2272) double_array((i-1)*(MAX_FACE_PER_CELL_SURF+1)+1) = &
2273) vec_ptr1((i-1)*offset+(dof-1)*MAX_FACE_PER_CELL_SURF+1)
2274) ! Divide the flowrate values by integration 'time'
2275) do iface = 1,MAX_FACE_PER_CELL_SURF
2276) double_array((i-1)*(MAX_FACE_PER_CELL_SURF+1)+iface+1) = &
2277) vec_ptr2((i-1)*offset + &
2278) (dof-1)*MAX_FACE_PER_CELL_SURF + iface + 1)/ &
2279) output_option%periodic_snap_output_time_incr
2280) enddo
2281) enddo
2282) end select
2283)
2284) !call PetscLogEventBegin(logging%event_h5dwrite_f,ierr)
2285) call h5dwrite_f(data_set_id,H5T_NATIVE_DOUBLE,double_array,dims, &
2286) hdf5_err,memory_space_id,file_space_id,prop_id)
2287) !call PetscLogEventEnd(logging%event_h5dwrite_f,ierr)
2288)
2289) call h5pclose_f(prop_id,hdf5_err)
2290)
2291) call h5dclose_f(data_set_id,hdf5_err)
2292) call h5sclose_f(file_space_id,hdf5_err)
2293)
2294) enddo
2295)
2296) #endif
2297) ! #ifdef SCORPIO_WRITE
2298)
2299) end subroutine WriteHDF5SurfaceFlowratesUGrid
2300) #endif
2301)
2302) end module Output_Surface_module