output_common.F90 coverage: 56.52 %func 51.04 %block
1) module Output_Common_module
2)
3) use Logging_module
4) use Output_Aux_module
5)
6) !note: only realization_base_type can be used throughout this module.
7) use Realization_Base_class, only : realization_base_type
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)
19) PetscInt, save, public :: max_local_size_saved = -1
20)
21) !geh: would prefer that this be local to Output_Tecplot_module, but needed
22) ! in Output_Surface_module
23) PetscInt, parameter, public :: TECPLOT_INTEGER = 0
24) PetscInt, parameter, public :: TECPLOT_REAL = 1
25)
26) public :: OutputCommonInit, &
27) OutputGetVarFromArray, &
28) OutputGetVarFromArrayAtCoord, &
29) OutputGetCellCenteredVelocities, &
30) ConvertArrayToNatural, &
31) GetCellCoordinates, &
32) GetVertexCoordinates, &
33) OutputFilenameID, &
34) OutputFilename, &
35) GetCellConnections, &
36) OutputXMFHeader, &
37) OutputXMFAttribute, &
38) OutputXMFFooter, &
39) OutputGetFaceVelUGrid, &
40) OutputGetFaceFlowrateUGrid, &
41) ExplicitGetCellCoordinates, &
42) OutputGetExplicitFlowrates, &
43) GetCellConnectionsExplicit, &
44) OutputXMFHeaderExplicit, &
45) OutputXMFAttributeExplicit, &
46) OutputGetExplicitIDsFlowrates, &
47) OutputGetExplicitAuxVars, &
48) OutputGetExplicitCellInfo
49)
50) contains
51)
52) ! ************************************************************************** !
53)
54) subroutine OutputCommonInit()
55) !
56) ! Initializes module variables for common formats
57) !
58) ! Author: Glenn Hammond
59) ! Date: 01/16/13
60) !
61)
62) use Option_module
63)
64) implicit none
65)
66) ! set size to -1 in order to re-initialize parallel communication blocks
67) max_local_size_saved = -1
68)
69) end subroutine OutputCommonInit
70)
71) ! ************************************************************************** !
72)
73) function OutputFilenameID(output_option,option)
74) !
75) ! Creates an ID for filename
76) !
77) ! Author: Glenn Hammond
78) ! Date: 01/13/12
79) !
80)
81) use Option_module
82)
83) implicit none
84)
85) type(option_type) :: option
86) type(output_option_type) :: output_option
87)
88) character(len=MAXWORDLENGTH) :: OutputFilenameID
89)
90) if (output_option%plot_number < 10) then
91) write(OutputFilenameID,'("00",i1)') output_option%plot_number
92) else if (output_option%plot_number < 100) then
93) write(OutputFilenameID,'("0",i2)') output_option%plot_number
94) else if (output_option%plot_number < 1000) then
95) write(OutputFilenameID,'(i3)') output_option%plot_number
96) else if (output_option%plot_number < 10000) then
97) write(OutputFilenameID,'(i4)') output_option%plot_number
98) endif
99)
100) OutputFilenameID = adjustl(OutputFilenameID)
101)
102) end function OutputFilenameID
103)
104) ! ************************************************************************** !
105)
106) function OutputFilename(output_option,option,suffix,optional_string)
107) !
108) ! Creates a filename for a Tecplot file
109) !
110) ! Author: Glenn Hammond
111) ! Date: 01/13/12
112) !
113)
114) use Option_module
115)
116) implicit none
117)
118) type(option_type) :: option
119) type(output_option_type) :: output_option
120) character(len=*) :: suffix
121) character(len=*) :: optional_string
122)
123) character(len=MAXSTRINGLENGTH) :: OutputFilename
124)
125) character(len=MAXWORDLENGTH) :: final_suffix
126) character(len=MAXSTRINGLENGTH) :: final_optional_string
127)
128)
129) if (len_trim(optional_string) > 0) then
130) final_optional_string = '-' // optional_string
131) else
132) final_optional_string = ''
133) endif
134) final_suffix = '.' // suffix
135)
136) ! open file
137) if (len_trim(output_option%plot_name) > 2) then
138) OutputFilename = trim(output_option%plot_name) // &
139) trim(final_optional_string) // &
140) final_suffix
141) else
142) OutputFilename = trim(option%global_prefix) // &
143) trim(option%group_prefix) // &
144) trim(final_optional_string) // &
145) '-' // &
146) trim(OutputFilenameID(output_option,option)) // &
147) final_suffix
148) endif
149)
150) end function OutputFilename
151)
152) ! ************************************************************************** !
153)
154) subroutine OutputGetVarFromArray(realization_base,vec,ivar,isubvar,isubvar1)
155) !
156) ! Extracts variables indexed by ivar from a multivar array
157) !
158) ! Author: Glenn Hammond
159) ! Date: 10/25/07
160) !
161)
162) use Realization_Base_class, only : RealizationGetVariable
163) use Grid_module
164) use Option_module
165) use Field_module
166)
167) implicit none
168)
169) #include "petsc/finclude/petscvec.h"
170) #include "petsc/finclude/petscvec.h90"
171) #include "petsc/finclude/petsclog.h"
172)
173) class(realization_base_type) :: realization_base
174) Vec :: vec
175) PetscInt :: ivar
176) PetscInt :: isubvar
177) PetscInt, optional :: isubvar1
178)
179) PetscErrorCode :: ierr
180)
181) call PetscLogEventBegin(logging%event_output_get_var_from_array, &
182) ierr);CHKERRQ(ierr)
183)
184) call RealizationGetVariable(realization_base,vec,ivar,isubvar,isubvar1)
185)
186) call PetscLogEventEnd(logging%event_output_get_var_from_array, &
187) ierr);CHKERRQ(ierr)
188)
189) end subroutine OutputGetVarFromArray
190)
191) ! ************************************************************************** !
192)
193) subroutine ConvertArrayToNatural(indices,array,local_size,global_size,option)
194) !
195) ! Converts an array to natural ordering
196) !
197) ! Author: Glenn Hammond
198) ! Date: 10/25/07
199) !
200)
201) use Option_module
202)
203) implicit none
204)
205) #include "petsc/finclude/petscvec.h"
206) #include "petsc/finclude/petscvec.h90"
207)
208) PetscInt :: local_size, global_size
209) PetscInt :: indices(:)
210) PetscReal, pointer :: array(:)
211) type(option_type) :: option
212)
213) Vec :: natural_vec
214) PetscInt, allocatable :: indices_zero_based(:)
215) PetscReal, pointer :: vec_ptr(:)
216) PetscErrorCode :: ierr
217)
218) call VecCreate(option%mycomm,natural_vec,ierr);CHKERRQ(ierr)
219) call VecSetSizes(natural_vec,PETSC_DECIDE,global_size,ierr);CHKERRQ(ierr)
220) call VecSetType(natural_vec,VECMPI,ierr);CHKERRQ(ierr)
221)
222) allocate(indices_zero_based(local_size))
223) indices_zero_based(1:local_size) = indices(1:local_size)-1
224)
225) call VecSetValues(natural_vec,local_size,indices_zero_based, &
226) array,INSERT_VALUES,ierr);CHKERRQ(ierr)
227)
228) call VecAssemblyBegin(natural_vec,ierr);CHKERRQ(ierr)
229) call VecAssemblyEnd(natural_vec,ierr);CHKERRQ(ierr)
230)
231) call VecGetLocalSize(natural_vec,local_size,ierr);CHKERRQ(ierr)
232) deallocate(array)
233) allocate(array(local_size))
234)
235) call VecGetArrayF90(natural_vec,vec_ptr,ierr);CHKERRQ(ierr)
236) array(1:local_size) = vec_ptr(1:local_size)
237) call VecRestoreArrayF90(natural_vec,vec_ptr,ierr);CHKERRQ(ierr)
238)
239) call VecDestroy(natural_vec,ierr);CHKERRQ(ierr)
240)
241) end subroutine ConvertArrayToNatural
242)
243) ! ************************************************************************** !
244)
245) function OutputGetVarFromArrayAtCoord(realization_base,ivar,isubvar,x,y,z, &
246) num_cells,ghosted_ids,isubvar1)
247) !
248) ! Extracts variables indexed by ivar from a multivar array
249) !
250) ! Author: Glenn Hammond
251) ! Date: 02/11/08
252) !
253)
254) use Realization_Base_class, only : RealizGetVariableValueAtCell
255) use Grid_module
256) use Option_module
257)
258) implicit none
259)
260) PetscReal :: OutputGetVarFromArrayAtCoord
261) class(realization_base_type) :: realization_base
262) PetscInt :: ivar
263) PetscInt :: isubvar
264) PetscInt, optional :: isubvar1
265) PetscReal :: x,y,z
266) PetscInt :: num_cells
267) PetscInt :: ghosted_ids(num_cells)
268)
269) type(grid_type), pointer :: grid
270) PetscInt :: icell, ghosted_id
271) PetscReal :: dx, dy, dz
272) PetscReal :: value, sum_value
273) PetscReal :: weight, sum_weight, sum_root
274)
275) sum_value = 0.d0
276) sum_weight = 0.d0
277)
278) grid => realization_base%patch%grid
279)
280) do icell=1, num_cells
281) ghosted_id = ghosted_ids(icell)
282) dx = x-grid%x(ghosted_id)
283) dy = y-grid%y(ghosted_id)
284) dz = z-grid%z(ghosted_id)
285) sum_root = sqrt(dx*dx+dy*dy+dz*dz)
286) value = 0.d0
287) value = RealizGetVariableValueAtCell(realization_base,ivar,isubvar,ghosted_id, &
288) isubvar1)
289) if (sum_root < 1.d-40) then ! bail because it is right on this coordinate
290) sum_weight = 1.d0
291) sum_value = value
292) exit
293) endif
294) weight = 1.d0/sum_root
295) sum_weight = sum_weight + weight
296) sum_value = sum_value + weight * value
297) enddo
298)
299) OutputGetVarFromArrayAtCoord = sum_value/sum_weight
300)
301) end function OutputGetVarFromArrayAtCoord
302)
303) ! ************************************************************************** !
304)
305) subroutine OutputGetCellCenteredVelocities(realization_base,vec_x,vec_y, &
306) vec_z, iphase)
307) !
308) ! Computes the cell-centered velocity component
309) ! as an averages of cell face velocities
310) !
311) ! Author: Glenn Hammond
312) ! Date: 10/25/07; refactored 01/31/14
313) !
314) use Logging_module
315) use Patch_module
316) use Grid_module
317)
318) implicit none
319)
320) #include "petsc/finclude/petscvec.h"
321) #include "petsc/finclude/petscvec.h90"
322)
323) class(realization_base_type) :: realization_base
324) Vec :: vec_x,vec_y,vec_z
325) PetscInt :: direction
326) PetscInt :: iphase
327)
328) PetscErrorCode :: ierr
329)
330) PetscReal, pointer :: vec_x_ptr(:),vec_y_ptr(:),vec_z_ptr(:)
331) PetscReal, allocatable :: velocities(:,:)
332)
333) call PetscLogEventBegin(logging%event_output_get_cell_vel, &
334) ierr);CHKERRQ(ierr)
335)
336) allocate(velocities(3,realization_base%patch%grid%nlmax))
337) call PatchGetCellCenteredVelocities(realization_base%patch,iphase,velocities)
338)
339) call VecGetArrayF90(vec_x,vec_x_ptr,ierr);CHKERRQ(ierr)
340) call VecGetArrayF90(vec_y,vec_y_ptr,ierr);CHKERRQ(ierr)
341) call VecGetArrayF90(vec_z,vec_z_ptr,ierr);CHKERRQ(ierr)
342)
343) vec_x_ptr(:) = velocities(X_DIRECTION,:)*realization_base%output_option%tconv
344) vec_y_ptr(:) = velocities(Y_DIRECTION,:)*realization_base%output_option%tconv
345) vec_z_ptr(:) = velocities(Z_DIRECTION,:)*realization_base%output_option%tconv
346)
347) call VecRestoreArrayF90(vec_x,vec_x_ptr,ierr);CHKERRQ(ierr)
348) call VecRestoreArrayF90(vec_y,vec_y_ptr,ierr);CHKERRQ(ierr)
349) call VecRestoreArrayF90(vec_z,vec_z_ptr,ierr);CHKERRQ(ierr)
350)
351) deallocate(velocities)
352)
353) call PetscLogEventEnd(logging%event_output_get_cell_vel,ierr);CHKERRQ(ierr)
354)
355) end subroutine OutputGetCellCenteredVelocities
356)
357) ! ************************************************************************** !
358)
359) subroutine GetCellCoordinates(grid,vec,direction)
360) !
361) ! Extracts coordinates of cells into a PetscVec
362) !
363) ! Author: Glenn Hammond
364) ! Date: 10/25/07
365) !
366)
367) use Grid_module
368) use Variables_module
369)
370) implicit none
371)
372) #include "petsc/finclude/petscvec.h"
373) #include "petsc/finclude/petscvec.h90"
374)
375) type(grid_type) :: grid
376) Vec :: vec
377) PetscInt :: direction
378) PetscErrorCode :: ierr
379)
380) PetscInt :: i
381) PetscReal, pointer :: vec_ptr(:)
382)
383) call VecGetArrayF90(vec,vec_ptr,ierr);CHKERRQ(ierr)
384)
385) if (direction == X_COORDINATE) then
386) do i = 1,grid%nlmax
387) vec_ptr(i) = grid%x(grid%nL2G(i))
388) enddo
389) else if (direction == Y_COORDINATE) then
390) do i = 1,grid%nlmax
391) vec_ptr(i) = grid%y(grid%nL2G(i))
392) enddo
393) else if (direction == Z_COORDINATE) then
394) do i = 1,grid%nlmax
395) vec_ptr(i) = grid%z(grid%nL2G(i))
396) enddo
397) endif
398)
399) call VecRestoreArrayF90(vec,vec_ptr,ierr);CHKERRQ(ierr)
400)
401) end subroutine GetCellCoordinates
402)
403) ! ************************************************************************** !
404)
405) subroutine GetVertexCoordinates(grid,vec,direction,option)
406) !
407) ! Extracts vertex coordinates of cells into a PetscVec
408) !
409) ! Author: Gautam Bisht
410) ! Date: 11/01/2011
411) !
412)
413) use Grid_module
414) use Option_module
415) use Variables_module, only : X_COORDINATE, Y_COORDINATE, Z_COORDINATE
416)
417) implicit none
418)
419) #include "petsc/finclude/petscvec.h"
420) #include "petsc/finclude/petscvec.h90"
421)
422) type(grid_type) :: grid
423) Vec :: vec
424) PetscInt :: direction
425) type(option_type) :: option
426)
427) PetscInt :: ivertex
428) PetscReal, pointer :: vec_ptr(:)
429) PetscInt, allocatable :: indices(:)
430) PetscReal, allocatable :: values(:)
431) PetscErrorCode :: ierr
432)
433) if (option%mycommsize == 1) then
434) call VecGetArrayF90(vec,vec_ptr,ierr);CHKERRQ(ierr)
435) select case(direction)
436) case(X_COORDINATE)
437) do ivertex = 1,grid%unstructured_grid%num_vertices_local
438) vec_ptr(ivertex) = grid%unstructured_grid%vertices(ivertex)%x
439) enddo
440) case(Y_COORDINATE)
441) do ivertex = 1,grid%unstructured_grid%num_vertices_local
442) vec_ptr(ivertex) = grid%unstructured_grid%vertices(ivertex)%y
443) enddo
444) case(Z_COORDINATE)
445) do ivertex = 1,grid%unstructured_grid%num_vertices_local
446) vec_ptr(ivertex) = grid%unstructured_grid%vertices(ivertex)%z
447) enddo
448) end select
449) call VecRestoreArrayF90(vec,vec_ptr,ierr);CHKERRQ(ierr)
450) else
451) ! initialize to UNINITIALIZED_INTEGER to catch bugs
452) call VecSet(vec,UNINITIALIZED_DOUBLE,ierr);CHKERRQ(ierr)
453) allocate(values(grid%unstructured_grid%num_vertices_local))
454) allocate(indices(grid%unstructured_grid%num_vertices_local))
455) select case(direction)
456) case(X_COORDINATE)
457) do ivertex = 1,grid%unstructured_grid%num_vertices_local
458) values(ivertex) = grid%unstructured_grid%vertices(ivertex)%x
459) enddo
460) case(Y_COORDINATE)
461) do ivertex = 1,grid%unstructured_grid%num_vertices_local
462) values(ivertex) = grid%unstructured_grid%vertices(ivertex)%y
463) enddo
464) case(Z_COORDINATE)
465) do ivertex = 1,grid%unstructured_grid%num_vertices_local
466) values(ivertex) = grid%unstructured_grid%vertices(ivertex)%z
467) enddo
468) end select
469) indices(:) = grid%unstructured_grid%vertex_ids_natural(:)-1
470) call VecSetValues(vec,grid%unstructured_grid%num_vertices_local, &
471) indices,values,INSERT_VALUES,ierr);CHKERRQ(ierr)
472) call VecAssemblyBegin(vec,ierr);CHKERRQ(ierr)
473) deallocate(values)
474) deallocate(indices)
475) call VecAssemblyEnd(vec,ierr);CHKERRQ(ierr)
476) endif
477)
478) end subroutine GetVertexCoordinates
479)
480) ! ************************************************************************** !
481)
482) subroutine ExplicitGetCellCoordinates(grid,vec,direction,option)
483) !
484) ! Extracts cell coordinates for explicit grid
485) ! into a PetscVec
486) !
487) ! Author: Satish Karra, LANL
488) ! Date: 12/11/12
489) !
490)
491) use Grid_module
492) use Option_module
493) use Variables_module, only : X_COORDINATE, Y_COORDINATE, Z_COORDINATE
494)
495) implicit none
496)
497) #include "petsc/finclude/petscvec.h"
498) #include "petsc/finclude/petscvec.h90"
499)
500) type(grid_type) :: grid
501) Vec :: vec
502) PetscInt :: direction
503) type(option_type) :: option
504)
505) PetscInt :: ivertex
506) PetscReal, pointer :: vec_ptr(:)
507) PetscInt, allocatable :: indices(:)
508) PetscReal, allocatable :: values(:)
509) PetscErrorCode :: ierr
510)
511) if (option%mycommsize == 1) then
512) call VecGetArrayF90(vec,vec_ptr,ierr);CHKERRQ(ierr)
513) select case(direction)
514) case(X_COORDINATE)
515) do ivertex = 1,grid%ngmax
516) vec_ptr(ivertex) = grid%unstructured_grid%explicit_grid%vertex_coordinates(ivertex)%x
517) enddo
518) case(Y_COORDINATE)
519) do ivertex = 1,grid%ngmax
520) vec_ptr(ivertex) = grid%unstructured_grid%explicit_grid%vertex_coordinates(ivertex)%y
521) enddo
522) case(Z_COORDINATE)
523) do ivertex = 1,grid%ngmax
524) vec_ptr(ivertex) = grid%unstructured_grid%explicit_grid%vertex_coordinates(ivertex)%z
525) enddo
526) end select
527) call VecRestoreArrayF90(vec,vec_ptr,ierr);CHKERRQ(ierr)
528) else
529) ! initialize to UNINITIALIZED_INTEGER to catch bugs
530) call VecSet(vec,UNINITIALIZED_DOUBLE,ierr);CHKERRQ(ierr)
531) allocate(values(grid%nlmax))
532) allocate(indices(grid%nlmax))
533) select case(direction)
534) case(X_COORDINATE)
535) do ivertex = 1,grid%nlmax
536) values(ivertex) = grid%unstructured_grid%explicit_grid%vertex_coordinates(ivertex)%x
537) enddo
538) case(Y_COORDINATE)
539) do ivertex = 1,grid%nlmax
540) values(ivertex) = grid%unstructured_grid%explicit_grid%vertex_coordinates(ivertex)%y
541) enddo
542) case(Z_COORDINATE)
543) do ivertex = 1,grid%nlmax
544) values(ivertex) = grid%unstructured_grid%explicit_grid%vertex_coordinates(ivertex)%z
545) enddo
546) end select
547) indices(:) = grid%unstructured_grid%cell_ids_natural(:)-1
548) call VecSetValues(vec,grid%nlmax,indices,values,INSERT_VALUES, &
549) ierr);CHKERRQ(ierr)
550) call VecAssemblyBegin(vec,ierr);CHKERRQ(ierr)
551) deallocate(values)
552) deallocate(indices)
553) call VecAssemblyEnd(vec,ierr);CHKERRQ(ierr)
554) endif
555)
556)
557) end subroutine ExplicitGetCellCoordinates
558)
559) ! ************************************************************************** !
560)
561) subroutine GetCellConnections(grid, vec)
562) !
563) ! This routine returns a vector containing vertex ids in natural order of
564) ! local cells for unstructured grid.
565) !
566) ! Author: Gautam Bisht, ORNL
567) ! Date: 05/31/12
568) !
569)
570) use Grid_module
571) use Grid_Unstructured_Aux_module
572) use Grid_Unstructured_Cell_module
573)
574) implicit none
575)
576) #include "petsc/finclude/petscvec.h"
577) #include "petsc/finclude/petscvec.h90"
578)
579) type(grid_type) :: grid
580) type(grid_unstructured_type),pointer :: ugrid
581) Vec :: vec
582) PetscInt :: local_id
583) PetscInt :: ghosted_id
584) PetscInt :: offset
585) PetscInt :: ivertex
586) PetscReal, pointer :: vec_ptr(:)
587) PetscErrorCode :: ierr
588)
589) ugrid => grid%unstructured_grid
590)
591) call VecGetArrayF90( vec, vec_ptr, ierr);CHKERRQ(ierr)
592)
593) ! initialize
594) vec_ptr = UNINITIALIZED_DOUBLE
595) do local_id=1, ugrid%nlmax
596) ghosted_id = local_id
597) select case(ugrid%cell_type(ghosted_id))
598) case(HEX_TYPE)
599) offset = (local_id-1)*8
600) do ivertex = 1, 8
601) vec_ptr(offset + ivertex) = &
602) ugrid%vertex_ids_natural(ugrid%cell_vertices(ivertex,local_id))
603) enddo
604) case(WEDGE_TYPE)
605) offset = (local_id-1)*8
606) do ivertex = 1, 6
607) vec_ptr(offset + ivertex) = &
608) ugrid%vertex_ids_natural(ugrid%cell_vertices(ivertex,local_id))
609) enddo
610) vec_ptr(offset + 7) = 0
611) vec_ptr(offset + 8) = 0
612) case (PYR_TYPE)
613) offset = (local_id-1)*8
614) do ivertex = 1, 5
615) vec_ptr(offset + ivertex) = &
616) ugrid%vertex_ids_natural(ugrid%cell_vertices(ivertex,local_id))
617) enddo
618) do ivertex = 6, 8
619) vec_ptr(offset + ivertex) = 0
620) enddo
621) case (TET_TYPE)
622) offset = (local_id-1)*8
623) do ivertex = 1, 4
624) vec_ptr(offset + ivertex) = &
625) ugrid%vertex_ids_natural(ugrid%cell_vertices(ivertex,local_id))
626) enddo
627) do ivertex = 5, 8
628) vec_ptr(offset + ivertex) = 0
629) enddo
630) case (QUAD_TYPE)
631) offset = (local_id-1)*4
632) do ivertex = 1, 4
633) vec_ptr(offset + ivertex) = &
634) ugrid%vertex_ids_natural(ugrid%cell_vertices(ivertex,local_id))
635) enddo
636) case (TRI_TYPE)
637) offset = (local_id-1)*4
638) do ivertex = 1, 3
639) vec_ptr(offset + ivertex) = &
640) ugrid%vertex_ids_natural(ugrid%cell_vertices(ivertex,local_id))
641) enddo
642) ivertex = 4
643) vec_ptr(offset + 4) = 0
644) end select
645) enddo
646)
647) call VecRestoreArrayF90( vec, vec_ptr, ierr);CHKERRQ(ierr)
648)
649) end subroutine GetCellConnections
650)
651) ! ************************************************************************** !
652)
653) subroutine GetCellConnectionsExplicit(grid, vec)
654) !
655) ! returns a vector containing vertex ids in natural order of
656) ! local cells for unstructured grid of explicit type
657) !
658) ! Author: Satish Karra
659) ! Date: 07/16/13
660) !
661)
662) use Grid_module
663) use Grid_Unstructured_Aux_module
664) use Grid_Unstructured_Cell_module
665)
666) implicit none
667)
668) #include "petsc/finclude/petscvec.h"
669) #include "petsc/finclude/petscvec.h90"
670)
671) type(grid_type) :: grid
672) type(grid_unstructured_type),pointer :: ugrid
673) type(unstructured_explicit_type), pointer :: explicit_grid
674) Vec :: vec
675) PetscInt :: offset
676) PetscInt :: ivertex, iconn
677) PetscReal, pointer :: vec_ptr(:)
678) PetscErrorCode :: ierr
679)
680) ugrid => grid%unstructured_grid
681) explicit_grid => ugrid%explicit_grid
682)
683) call VecGetArrayF90( vec, vec_ptr, ierr);CHKERRQ(ierr)
684)
685) ! initialize
686) vec_ptr = UNINITIALIZED_DOUBLE
687) do iconn = 1, explicit_grid%num_elems
688) select case(explicit_grid%cell_connectivity(0,iconn))
689) case(8)
690) offset = (iconn-1)*8
691) do ivertex = 1, 8
692) vec_ptr(offset + ivertex) = &
693) explicit_grid%cell_connectivity(ivertex,iconn)
694) enddo
695) case(6)
696) offset = (iconn-1)*8
697) do ivertex = 1, 6
698) vec_ptr(offset + ivertex) = &
699) explicit_grid%cell_connectivity(ivertex,iconn)
700) enddo
701) vec_ptr(offset + 7) = 0
702) vec_ptr(offset + 8) = 0
703) case (5)
704) offset = (iconn-1)*8
705) do ivertex = 1, 5
706) vec_ptr(offset + ivertex) = &
707) explicit_grid%cell_connectivity(ivertex,iconn)
708) enddo
709) do ivertex = 6, 8
710) vec_ptr(offset + ivertex) = 0
711) enddo
712) case (4)
713) if (grid%unstructured_grid%grid_type /= TWO_DIM_GRID) then
714) offset = (iconn-1)*8
715) do ivertex = 1, 4
716) vec_ptr(offset + ivertex) = &
717) explicit_grid%cell_connectivity(ivertex,iconn)
718) enddo
719) do ivertex = 5, 8
720) vec_ptr(offset + ivertex) = 0
721) enddo
722) else
723) offset = (iconn-1)*8
724) do ivertex = 1, 4
725) vec_ptr(offset + ivertex) = &
726) explicit_grid%cell_connectivity(ivertex,iconn)
727) enddo
728) do ivertex = 5, 8
729) vec_ptr(offset + ivertex) = 0
730) enddo
731) endif
732) case (3)
733) offset = (iconn-1)*8
734) do ivertex = 1, 3
735) vec_ptr(offset + ivertex) = &
736) explicit_grid%cell_connectivity(ivertex,iconn)
737) enddo
738) do ivertex = 4, 8
739) vec_ptr(offset + ivertex) = 0
740) enddo
741) end select
742) enddo
743)
744) call VecRestoreArrayF90( vec, vec_ptr, ierr);CHKERRQ(ierr)
745)
746) end subroutine GetCellConnectionsExplicit
747)
748) ! ************************************************************************** !
749)
750) subroutine OutputXMFHeader(fid,time,nmax,xmf_vert_len,ngvert,filename)
751) !
752) ! This subroutine writes header to a .xmf file
753) !
754) ! Author: Gautam Bisht, LBNL
755) ! Date: 10/29/12
756) !
757)
758) implicit none
759)
760) PetscInt :: fid, vert_count
761) PetscReal :: time
762) PetscInt :: nmax,xmf_vert_len,ngvert
763) character(len=MAXSTRINGLENGTH) :: filename
764)
765) character(len=MAXSTRINGLENGTH) :: string, string2
766) character(len=MAXWORDLENGTH) :: word
767)
768) string="<?xml version=""1.0"" ?>"
769) write(fid,'(a)') trim(string)
770)
771) string="<!DOCTYPE Xdmf SYSTEM ""Xdmf.dtd"" []>"
772) write(fid,'(a)') trim(string)
773)
774) string="<Xdmf>"
775) write(fid,'(a)') trim(string)
776)
777) string=" <Domain>"
778) write(fid,'(a)') trim(string)
779)
780) string=" <Grid Name=""Mesh"">"
781) write(fid,'(a)') trim(string)
782)
783) write(string2,'(es13.5)') time
784) string=" <Time Value = """ // trim(adjustl(string2)) // """ />"
785) write(fid,'(a)') trim(string)
786)
787) write(string2,*) nmax
788) string=" <Topology Type=""Mixed"" NumberOfElements=""" // &
789) trim(adjustl(string2)) // """ >"
790) write(fid,'(a)') trim(string)
791)
792) write(string2,*) xmf_vert_len
793) string=" <DataItem Format=""HDF"" DataType=""Int"" Dimensions=""" // &
794) trim(adjustl(string2)) // """>"
795) write(fid,'(a)') trim(string)
796)
797) string=" "//trim(filename) //":/Domain/Cells"
798) write(fid,'(a)') trim(string)
799)
800) string=" </DataItem>"
801) write(fid,'(a)') trim(string)
802)
803) string=" </Topology>"
804) write(fid,'(a)') trim(string)
805)
806) string=" <Geometry GeometryType=""XYZ"">"
807) write(fid,'(a)') trim(string)
808)
809) write(string2,*) ngvert
810) string=" <DataItem Format=""HDF"" Dimensions=""" // trim(adjustl(string2)) // " 3"">"
811) write(fid,'(a)') trim(string)
812)
813) string=" "//trim(filename) //":/Domain/Vertices"
814) write(fid,'(a)') trim(string)
815)
816) string=" </DataItem>"
817) write(fid,'(a)') trim(string)
818)
819) string=" </Geometry>"
820) write(fid,'(a)') trim(string)
821)
822) string=" <Attribute Name=""XC"" AttributeType=""Scalar"" Center=""Cell"">"
823) write(fid,'(a)') trim(string)
824)
825) write(string2,*) nmax
826) string=" <DataItem Dimensions=""" // trim(adjustl(string2)) // " 1"" Format=""HDF""> "
827) write(fid,'(a)') trim(string)
828)
829) string=" "//trim(filename) //":/Domain/XC"
830) write(fid,'(a)') trim(string)
831)
832) string=" </DataItem> "
833) write(fid,'(a)') trim(string)
834)
835) string=" </Attribute>"
836) write(fid,'(a)') trim(string)
837)
838) string=" <Attribute Name=""YC"" AttributeType=""Scalar"" Center=""Cell"">"
839) write(fid,'(a)') trim(string)
840)
841) write(string2,*) nmax
842) string=" <DataItem Dimensions=""" // trim(adjustl(string2)) // " 1"" Format=""HDF""> "
843) write(fid,'(a)') trim(string)
844)
845) string=" "//trim(filename) //":/Domain/YC"
846) write(fid,'(a)') trim(string)
847)
848) string=" </DataItem> "
849) write(fid,'(a)') trim(string)
850)
851) string=" </Attribute>"
852) write(fid,'(a)') trim(string)
853)
854) string=" <Attribute Name=""ZC"" AttributeType=""Scalar"" Center=""Cell"">"
855) write(fid,'(a)') trim(string)
856)
857) write(string2,*) nmax
858) string=" <DataItem Dimensions=""" // trim(adjustl(string2)) // " 1"" Format=""HDF""> "
859) write(fid,'(a)') trim(string)
860)
861) string=" "//trim(filename) //":/Domain/ZC"
862) write(fid,'(a)') trim(string)
863)
864) string=" </DataItem> "
865) write(fid,'(a)') trim(string)
866)
867) string=" </Attribute>"
868) write(fid,'(a)') trim(string)
869)
870) end subroutine OutputXMFHeader
871)
872) ! ************************************************************************** !
873)
874) subroutine OutputXMFHeaderExplicit(fid,time,nmax,xmf_vert_len,ngvert,filename)
875) !
876) ! Header for xdmf output with explicit grid
877) !
878) ! Author: Satish Karra
879) ! Date: 07/17/13
880) !
881)
882) implicit none
883)
884) PetscInt :: fid, vert_count
885) PetscReal :: time
886) PetscInt :: nmax,xmf_vert_len,ngvert
887) character(len=MAXSTRINGLENGTH) :: filename
888)
889) character(len=MAXSTRINGLENGTH) :: string, string2
890) character(len=MAXWORDLENGTH) :: word
891)
892) string="<?xml version=""1.0"" ?>"
893) write(fid,'(a)') trim(string)
894)
895) string="<!DOCTYPE Xdmf SYSTEM ""Xdmf.dtd"" []>"
896) write(fid,'(a)') trim(string)
897)
898) string="<Xdmf>"
899) write(fid,'(a)') trim(string)
900)
901) string=" <Domain>"
902) write(fid,'(a)') trim(string)
903)
904) string=" <Grid Name=""Mesh"">"
905) write(fid,'(a)') trim(string)
906)
907) write(string2,'(es13.5)') time
908) string=" <Time Value = """ // trim(adjustl(string2)) // """ />"
909) write(fid,'(a)') trim(string)
910)
911) write(string2,*) nmax
912) string=" <Topology Type=""Mixed"" NumberOfElements=""" // &
913) trim(adjustl(string2)) // """ >"
914) write(fid,'(a)') trim(string)
915)
916) write(string2,*) xmf_vert_len
917) string=" <DataItem Format=""HDF"" DataType=""Int"" Dimensions=""" // &
918) trim(adjustl(string2)) // """>"
919) write(fid,'(a)') trim(string)
920)
921) string=" "//trim(filename) //":/Domain/Cells"
922) write(fid,'(a)') trim(string)
923)
924) string=" </DataItem>"
925) write(fid,'(a)') trim(string)
926)
927) string=" </Topology>"
928) write(fid,'(a)') trim(string)
929)
930) string=" <Geometry GeometryType=""XYZ"">"
931) write(fid,'(a)') trim(string)
932)
933) write(string2,*) ngvert
934) string=" <DataItem Format=""HDF"" Dimensions=""" // trim(adjustl(string2)) // " 3"">"
935) write(fid,'(a)') trim(string)
936)
937) string=" "//trim(filename) //":/Domain/Vertices"
938) write(fid,'(a)') trim(string)
939)
940) string=" </DataItem>"
941) write(fid,'(a)') trim(string)
942)
943) string=" </Geometry>"
944) write(fid,'(a)') trim(string)
945)
946) end subroutine OutputXMFHeaderExplicit
947)
948) ! ************************************************************************** !
949)
950) subroutine OutputXMFFooter(fid)
951) !
952) ! This subroutine writes footer to a .xmf file
953) !
954) ! Author: Gautam Bisht, LBNL
955) ! Date: 10/29/12
956) !
957)
958) implicit none
959)
960) PetscInt :: fid
961)
962) character(len=MAXSTRINGLENGTH) :: string
963)
964) string=" </Grid>"
965) write(fid,'(a)') trim(string)
966)
967) string=" </Domain>"
968) write(fid,'(a)') trim(string)
969)
970) string="</Xdmf>"
971) write(fid,'(a)') trim(string)
972)
973) end subroutine OutputXMFFooter
974)
975) ! ************************************************************************** !
976)
977) subroutine OutputXMFAttribute(fid,nmax,attname,att_datasetname)
978) !
979) ! This subroutine writes an attribute to a .xmf file
980) !
981) ! Author: Gautam Bisht, LBNL
982) ! Date: 10/29/12
983) !
984)
985) implicit none
986)
987) PetscInt :: fid,nmax
988)
989) character(len=MAXSTRINGLENGTH) :: attname, att_datasetname
990) character(len=MAXSTRINGLENGTH) :: string,string2
991) string=" <Attribute Name=""" // trim(attname) // &
992) """ AttributeType=""Scalar"" Center=""Cell"">"
993) write(fid,'(a)') trim(string)
994)
995) ! write(string2,*) grid%nmax
996) write(string2,*) nmax
997) string=" <DataItem Dimensions=""" // trim(adjustl(string2)) // " 1"" Format=""HDF""> "
998) write(fid,'(a)') trim(string)
999)
1000) string=" " // trim(att_datasetname)
1001) write(fid,'(a)') trim(string)
1002)
1003) string=" </DataItem> "
1004) write(fid,'(a)') trim(string)
1005)
1006) string=" </Attribute>"
1007) write(fid,'(a)') trim(string)
1008)
1009) end subroutine OutputXMFAttribute
1010)
1011) ! ************************************************************************** !
1012)
1013) subroutine OutputXMFAttributeExplicit(fid,nmax,attname,att_datasetname)
1014) !
1015) ! Header for xdmf attribute with explicit grid
1016) !
1017) ! Author: Satish Karra
1018) ! Date: 07/17/13
1019) !
1020)
1021) implicit none
1022)
1023) PetscInt :: fid,nmax
1024)
1025) character(len=MAXSTRINGLENGTH) :: attname, att_datasetname
1026) character(len=MAXSTRINGLENGTH) :: string,string2
1027) string=" <Attribute Name=""" // trim(attname) // &
1028) """ AttributeType=""Scalar"" Center=""Node"">"
1029) write(fid,'(a)') trim(string)
1030)
1031) ! write(string2,*) grid%nmax
1032) write(string2,*) nmax
1033) string=" <DataItem Dimensions=""" // trim(adjustl(string2)) // " 1"" Format=""HDF""> "
1034) write(fid,'(a)') trim(string)
1035)
1036) string=" " // trim(att_datasetname)
1037) write(fid,'(a)') trim(string)
1038)
1039) string=" </DataItem> "
1040) write(fid,'(a)') trim(string)
1041)
1042) string=" </Attribute>"
1043) write(fid,'(a)') trim(string)
1044)
1045) end subroutine OutputXMFAttributeExplicit
1046)
1047) ! ************************************************************************** !
1048)
1049) subroutine OutputGetFaceVelUGrid(realization_base)
1050) !
1051) ! This subroutine saves:
1052) ! - Face elocities at x/y/z directions, or
1053) !
1054) ! Author: Gautam Bisht, LBNL
1055) ! Date: 06/15/2016
1056) !
1057)
1058) use HDF5_module
1059) use Realization_Base_class, only : realization_base_type
1060) use Patch_module
1061) use Grid_module
1062) use Option_module
1063) use Grid_Unstructured_Aux_module
1064) use Grid_Unstructured_Cell_module
1065) use Variables_module
1066) use Connection_module
1067) use Coupler_module
1068) use HDF5_Aux_module
1069) use Output_Aux_module
1070) use Field_module
1071)
1072) implicit none
1073)
1074) #include "petsc/finclude/petscvec.h"
1075) #include "petsc/finclude/petscvec.h90"
1076) #include "petsc/finclude/petsclog.h"
1077) #include "petsc/finclude/petscsys.h"
1078)
1079) class(realization_base_type) :: realization_base
1080)
1081) type(option_type), pointer :: option
1082) type(patch_type), pointer :: patch
1083) type(grid_type), pointer :: grid
1084) type(grid_unstructured_type),pointer :: ugrid
1085) type(connection_set_list_type), pointer :: connection_set_list
1086) type(connection_set_type), pointer :: cur_connection_set
1087) type(coupler_type), pointer :: boundary_condition
1088) type(ugdm_type),pointer :: ugdm
1089) type(output_option_type), pointer :: output_option
1090) type(field_type), pointer :: field
1091)
1092) PetscInt :: local_id
1093) PetscInt :: ghosted_id
1094) PetscInt :: idual
1095) PetscInt :: iconn
1096) PetscInt :: face_id
1097) PetscInt :: local_id_up,local_id_dn
1098) PetscInt :: ghosted_id_up,ghosted_id_dn
1099) PetscInt :: iface_up,iface_dn
1100) PetscInt :: dof
1101) PetscInt :: sum_connection
1102) PetscInt :: offset
1103) PetscInt :: cell_type
1104) PetscInt :: local_size
1105) PetscInt :: i
1106) PetscInt :: iface
1107) PetscInt :: ndof
1108) PetscInt :: idx
1109)
1110) PetscReal, pointer :: flowrates(:,:,:)
1111) PetscReal, pointer :: vx(:,:,:)
1112) PetscReal, pointer :: vy(:,:,:)
1113) PetscReal, pointer :: vz(:,:,:)
1114) PetscReal, pointer :: vec_ptr(:)
1115) PetscReal, pointer :: vec_ptr2(:)
1116) PetscReal, pointer :: vec_ptr3(:)
1117) PetscReal, pointer :: vx_ptr(:)
1118) PetscReal, pointer :: vy_ptr(:)
1119) PetscReal, pointer :: vz_ptr(:)
1120) PetscReal, pointer :: double_array(:)
1121) PetscReal :: vel_vector(3)
1122) PetscReal :: dtime
1123)
1124) Vec :: natural_flowrates_vec
1125) Vec :: natural_vx_vec
1126) Vec :: natural_vy_vec
1127) Vec :: natural_vz_vec
1128)
1129) PetscMPIInt :: hdf5_err
1130) PetscErrorCode :: ierr
1131)
1132) character(len=MAXSTRINGLENGTH) :: string
1133) character(len=MAXWORDLENGTH) :: unit_string
1134)
1135) patch => realization_base%patch
1136) grid => patch%grid
1137) ugrid => grid%unstructured_grid
1138) output_option =>realization_base%output_option
1139) option => realization_base%option
1140) field => realization_base%field
1141)
1142) ! Create UGDM for
1143) call UGridCreateUGDM(grid%unstructured_grid,ugdm, &
1144) (option%nflowspec*MAX_FACE_PER_CELL + 1),option)
1145)
1146) ! Create vectors in natural order for velocity in x/y/z direction
1147) call UGridDMCreateVector(grid%unstructured_grid,ugdm,natural_vx_vec, &
1148) NATURAL,option)
1149) call UGridDMCreateVector(grid%unstructured_grid,ugdm,natural_vy_vec, &
1150) NATURAL,option)
1151) call UGridDMCreateVector(grid%unstructured_grid,ugdm,natural_vz_vec, &
1152) NATURAL,option)
1153)
1154) allocate(vx(option%nflowspec,MAX_FACE_PER_CELL,ugrid%nlmax))
1155) allocate(vy(option%nflowspec,MAX_FACE_PER_CELL,ugrid%nlmax))
1156) allocate(vz(option%nflowspec,MAX_FACE_PER_CELL,ugrid%nlmax))
1157)
1158) vx = 0.d0
1159) vy = 0.d0
1160) vz = 0.d0
1161)
1162) call VecGetArrayF90(field%vx_face_inst,vx_ptr,ierr);CHKERRQ(ierr)
1163) call VecGetArrayF90(field%vy_face_inst,vy_ptr,ierr);CHKERRQ(ierr)
1164) call VecGetArrayF90(field%vz_face_inst,vz_ptr,ierr);CHKERRQ(ierr)
1165)
1166) vx_ptr = 0.d0
1167) vy_ptr = 0.d0
1168) vz_ptr = 0.d0
1169)
1170) offset = 1 + option%nflowspec*MAX_FACE_PER_CELL
1171)
1172) ! Save the number of faces of all cell
1173) do local_id = 1,grid%nlmax
1174) ghosted_id = grid%nL2G(local_id)
1175) cell_type = ugrid%cell_type(ghosted_id)
1176)
1177) vx_ptr((local_id-1)*offset+1) = UCellGetNFaces(cell_type,option)
1178) vy_ptr((local_id-1)*offset+1) = UCellGetNFaces(cell_type,option)
1179) vz_ptr((local_id-1)*offset+1) = UCellGetNFaces(cell_type,option)
1180) enddo
1181)
1182) ! Interior Flowrates Terms -----------------------------------
1183) connection_set_list => grid%internal_connection_set_list
1184) cur_connection_set => connection_set_list%first
1185) sum_connection = 0
1186) do
1187) if (.not.associated(cur_connection_set)) exit
1188)
1189) do iconn = 1, cur_connection_set%num_connections
1190) sum_connection = sum_connection + 1
1191) face_id = cur_connection_set%face_id(iconn)
1192) ghosted_id_up = cur_connection_set%id_up(iconn)
1193) ghosted_id_dn = cur_connection_set%id_dn(iconn)
1194) local_id_up = grid%nG2L(ghosted_id_up)
1195) local_id_dn = grid%nG2L(ghosted_id_dn)
1196) do iface_up = 1,MAX_FACE_PER_CELL
1197) if (face_id==ugrid%cell_to_face_ghosted(iface_up,local_id_up)) exit
1198) enddo
1199)
1200) iface_dn=-1
1201) if (local_id_dn>0) then
1202) do iface_dn = 1,MAX_FACE_PER_CELL
1203) if (face_id==ugrid%cell_to_face_ghosted(iface_dn,local_id_dn)) exit
1204) enddo
1205) endif
1206)
1207) do dof=1,option%nflowspec
1208)
1209) ! Save velocity for iface_up of local_id_up cell using flowrate up-->dn
1210) vel_vector = cur_connection_set%dist(1:3,iconn)* &
1211) patch%internal_velocities(dof,sum_connection)
1212)
1213) vx(dof,iface_up,local_id_up) = vel_vector(1)
1214) vy(dof,iface_up,local_id_up) = vel_vector(2)
1215) vz(dof,iface_up,local_id_up) = vel_vector(3)
1216)
1217) idx = (local_id_up-1)*offset + (dof-1)*MAX_FACE_PER_CELL + iface_up + 1
1218)
1219) vx_ptr(idx) = vel_vector(1)
1220) vy_ptr(idx) = vel_vector(2)
1221) vz_ptr(idx) = vel_vector(3)
1222)
1223) if (iface_dn>0) then
1224)
1225) ! Save velocity for iface_dn of local_id_dn cell using -ve flowrate up-->dn
1226)
1227) idx = (local_id_dn-1)*offset + (dof-1)*MAX_FACE_PER_CELL + iface_dn + 1
1228)
1229) vx(dof,iface_dn,local_id_dn) = -vel_vector(1)
1230) vy(dof,iface_dn,local_id_dn) = -vel_vector(2)
1231) vz(dof,iface_dn,local_id_dn) = -vel_vector(3)
1232)
1233) vx_ptr(idx) = -vel_vector(1)
1234) vx_ptr(idx) = -vel_vector(2)
1235) vx_ptr(idx) = -vel_vector(3)
1236)
1237) endif
1238)
1239) enddo ! dof-loop
1240)
1241) enddo ! iconn-loop
1242)
1243) cur_connection_set => cur_connection_set%next
1244)
1245) enddo
1246)
1247) ! Boundary Flowrates Terms -----------------------------------
1248) boundary_condition => patch%boundary_condition_list%first
1249) sum_connection = 0
1250) do
1251) if (.not.associated(boundary_condition)) exit
1252)
1253) cur_connection_set => boundary_condition%connection_set
1254)
1255) do iconn = 1, cur_connection_set%num_connections
1256)
1257) sum_connection = sum_connection + 1
1258) face_id = cur_connection_set%face_id(iconn)
1259) ghosted_id_dn = cur_connection_set%id_dn(iconn)
1260) local_id_dn = grid%nG2L(ghosted_id_dn)
1261)
1262) do iface_dn = 1,MAX_FACE_PER_CELL
1263) if (face_id==ugrid%cell_to_face_ghosted(iface_dn,local_id_dn)) exit
1264) enddo
1265)
1266) do dof=1,option%nflowspec
1267)
1268) ! Save velocity for iface_dn of local_id_dn cell using -ve flowrate up-->dn
1269) vel_vector = cur_connection_set%dist(1:3,iconn)* &
1270) patch%boundary_velocities(dof,sum_connection)
1271)
1272) idx = (local_id_dn-1)*offset + (dof-1)*MAX_FACE_PER_CELL + iface_dn + 1
1273)
1274) vx(dof,iface_dn,local_id_dn) = -vel_vector(1)
1275) vy(dof,iface_dn,local_id_dn) = -vel_vector(2)
1276) vz(dof,iface_dn,local_id_dn) = -vel_vector(3)
1277)
1278) vx_ptr(idx) = -vel_vector(1)
1279) vx_ptr(idx) = -vel_vector(2)
1280) vx_ptr(idx) = -vel_vector(3)
1281)
1282) enddo ! dof-loop
1283)
1284) enddo ! iconn-loop
1285)
1286) boundary_condition => boundary_condition%next
1287)
1288) enddo
1289)
1290) deallocate(vx)
1291) deallocate(vy)
1292) deallocate(vz)
1293)
1294) call VecRestoreArrayF90(field%vx_face_inst,vx_ptr,ierr);CHKERRQ(ierr)
1295) call VecRestoreArrayF90(field%vy_face_inst,vy_ptr,ierr);CHKERRQ(ierr)
1296) call VecRestoreArrayF90(field%vz_face_inst,vz_ptr,ierr);CHKERRQ(ierr)
1297)
1298) ! Scatter flowrate from Global --> Natural order
1299) call VecScatterBegin(ugdm%scatter_gton,field%vx_face_inst,natural_vx_vec, &
1300) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1301) call VecScatterEnd(ugdm%scatter_gton,field%vx_face_inst,natural_vx_vec, &
1302) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1303)
1304) call VecScatterBegin(ugdm%scatter_gton,field%vy_face_inst,natural_vy_vec, &
1305) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1306) call VecScatterEnd(ugdm%scatter_gton,field%vy_face_inst,natural_vy_vec, &
1307) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1308)
1309) call VecScatterBegin(ugdm%scatter_gton,field%vz_face_inst,natural_vz_vec, &
1310) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1311) call VecScatterEnd(ugdm%scatter_gton,field%vz_face_inst,natural_vz_vec, &
1312) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1313)
1314) ! X-direction
1315) call VecGetArrayF90(natural_vx_vec,vec_ptr,ierr);CHKERRQ(ierr)
1316) call VecGetArrayF90(field%vx_face_inst,vec_ptr2,ierr);CHKERRQ(ierr)
1317)
1318) ! Copy the vectors
1319) vec_ptr2 = vec_ptr
1320) call VecRestoreArrayF90(natural_vx_vec,vec_ptr,ierr);CHKERRQ(ierr)
1321) call VecRestoreArrayF90(field%vx_face_inst,vec_ptr2,ierr);CHKERRQ(ierr)
1322)
1323) ! Y-direction
1324) call VecGetArrayF90(natural_vy_vec,vec_ptr,ierr);CHKERRQ(ierr)
1325) call VecGetArrayF90(field%vy_face_inst,vec_ptr2,ierr);CHKERRQ(ierr)
1326)
1327) ! Copy the vectors
1328) vec_ptr2 = vec_ptr
1329) call VecRestoreArrayF90(natural_vy_vec,vec_ptr,ierr);CHKERRQ(ierr)
1330) call VecRestoreArrayF90(field%vy_face_inst,vec_ptr2,ierr);CHKERRQ(ierr)
1331)
1332) ! Z-direction
1333) call VecGetArrayF90(natural_vz_vec,vec_ptr,ierr);CHKERRQ(ierr)
1334) call VecGetArrayF90(field%vz_face_inst,vec_ptr2,ierr);CHKERRQ(ierr)
1335)
1336) ! Copy the vectors
1337) vec_ptr2 = vec_ptr
1338) call VecRestoreArrayF90(natural_vz_vec,vec_ptr,ierr);CHKERRQ(ierr)
1339) call VecRestoreArrayF90(field%vz_face_inst,vec_ptr2,ierr);CHKERRQ(ierr)
1340)
1341) call VecDestroy(natural_vx_vec,ierr);CHKERRQ(ierr)
1342) call VecDestroy(natural_vy_vec,ierr);CHKERRQ(ierr)
1343) call VecDestroy(natural_vz_vec,ierr);CHKERRQ(ierr)
1344)
1345) call UGridDMDestroy(ugdm)
1346)
1347) end subroutine OutputGetFaceVelUGrid
1348)
1349) ! ************************************************************************** !
1350)
1351) subroutine OutputGetFaceFlowrateUGrid(realization_base)
1352) !
1353) ! This subroutine saves:
1354) ! - Mass/energy flowrate at all faces of a control volume
1355) !
1356) ! Author: Gautam Bisht, LBNL
1357) ! Date: 06/15/2016
1358) !
1359)
1360) use HDF5_module
1361) use Realization_Base_class, only : realization_base_type
1362) use Patch_module
1363) use Grid_module
1364) use Option_module
1365) use Grid_Unstructured_Aux_module
1366) use Grid_Unstructured_Cell_module
1367) use Variables_module
1368) use Connection_module
1369) use Coupler_module
1370) use HDF5_Aux_module
1371) use Output_Aux_module
1372) use Field_module
1373)
1374) implicit none
1375)
1376) #include "petsc/finclude/petscvec.h"
1377) #include "petsc/finclude/petscvec.h90"
1378) #include "petsc/finclude/petsclog.h"
1379) #include "petsc/finclude/petscsys.h"
1380)
1381) class(realization_base_type) :: realization_base
1382)
1383) type(option_type), pointer :: option
1384) type(patch_type), pointer :: patch
1385) type(grid_type), pointer :: grid
1386) type(grid_unstructured_type),pointer :: ugrid
1387) type(connection_set_list_type), pointer :: connection_set_list
1388) type(connection_set_type), pointer :: cur_connection_set
1389) type(coupler_type), pointer :: boundary_condition
1390) type(ugdm_type),pointer :: ugdm
1391) type(output_option_type), pointer :: output_option
1392) type(field_type), pointer :: field
1393)
1394) PetscInt :: local_id
1395) PetscInt :: ghosted_id
1396) PetscInt :: idual
1397) PetscInt :: iconn
1398) PetscInt :: face_id
1399) PetscInt :: local_id_up,local_id_dn
1400) PetscInt :: ghosted_id_up,ghosted_id_dn
1401) PetscInt :: iface_up,iface_dn
1402) PetscInt :: dof
1403) PetscInt :: sum_connection
1404) PetscInt :: offset
1405) PetscInt :: cell_type
1406) PetscInt :: local_size
1407) PetscInt :: i
1408) PetscInt :: iface
1409) PetscInt :: ndof
1410) PetscInt :: idx
1411)
1412) PetscReal, pointer :: flowrates(:,:,:)
1413) PetscReal, pointer :: vx(:,:,:)
1414) PetscReal, pointer :: vy(:,:,:)
1415) PetscReal, pointer :: vz(:,:,:)
1416) PetscReal, pointer :: vec_ptr(:)
1417) PetscReal, pointer :: vec_ptr2(:)
1418) PetscReal, pointer :: vec_ptr3(:)
1419) PetscReal, pointer :: vx_ptr(:)
1420) PetscReal, pointer :: vy_ptr(:)
1421) PetscReal, pointer :: vz_ptr(:)
1422) PetscReal, pointer :: double_array(:)
1423) PetscReal :: vel_vector(3)
1424) PetscReal :: dtime
1425)
1426) Vec :: natural_flowrates_vec
1427) Vec :: natural_vx_vec
1428) Vec :: natural_vy_vec
1429) Vec :: natural_vz_vec
1430)
1431) PetscMPIInt :: hdf5_err
1432) PetscErrorCode :: ierr
1433)
1434) character(len=MAXSTRINGLENGTH) :: string
1435) character(len=MAXWORDLENGTH) :: unit_string
1436)
1437) patch => realization_base%patch
1438) grid => patch%grid
1439) ugrid => grid%unstructured_grid
1440) output_option =>realization_base%output_option
1441) option => realization_base%option
1442) field => realization_base%field
1443)
1444) ! Create UGDM for
1445) call UGridCreateUGDM(grid%unstructured_grid,ugdm, &
1446) (option%nflowdof*MAX_FACE_PER_CELL + 1),option)
1447)
1448) ! Create a flowrate vector in natural order
1449) call UGridDMCreateVector(grid%unstructured_grid,ugdm,natural_flowrates_vec, &
1450) NATURAL,option)
1451)
1452) allocate(flowrates(option%nflowdof,MAX_FACE_PER_CELL,ugrid%nlmax))
1453) flowrates = 0.d0
1454)
1455) call VecGetArrayF90(field%flowrate_inst,vec_ptr,ierr);CHKERRQ(ierr)
1456) vec_ptr = 0.d0
1457)
1458) offset = 1 + option%nflowdof*MAX_FACE_PER_CELL
1459) ! Save the number of faces of all cell
1460) do local_id = 1,grid%nlmax
1461) ghosted_id = grid%nL2G(local_id)
1462) cell_type = ugrid%cell_type(ghosted_id)
1463) vec_ptr((local_id-1)*offset+1) = UCellGetNFaces(cell_type,option)
1464) enddo
1465)
1466) ! Interior Flowrates Terms -----------------------------------
1467) connection_set_list => grid%internal_connection_set_list
1468) cur_connection_set => connection_set_list%first
1469) sum_connection = 0
1470) do
1471) if (.not.associated(cur_connection_set)) exit
1472)
1473) do iconn = 1, cur_connection_set%num_connections
1474) sum_connection = sum_connection + 1
1475) face_id = cur_connection_set%face_id(iconn)
1476) ghosted_id_up = cur_connection_set%id_up(iconn)
1477) ghosted_id_dn = cur_connection_set%id_dn(iconn)
1478) local_id_up = grid%nG2L(ghosted_id_up)
1479) local_id_dn = grid%nG2L(ghosted_id_dn)
1480)
1481) do iface_up = 1,MAX_FACE_PER_CELL
1482) if (face_id==ugrid%cell_to_face_ghosted(iface_up,local_id_up)) exit
1483) enddo
1484)
1485) iface_dn=-1
1486) if (local_id_dn>0) then
1487) do iface_dn = 1,MAX_FACE_PER_CELL
1488) if (face_id==ugrid%cell_to_face_ghosted(iface_dn,local_id_dn)) exit
1489) enddo
1490) endif
1491)
1492) do dof=1,option%nflowdof
1493)
1494) ! Save flowrate for iface_up of local_id_up cell using flowrate up-->dn
1495) flowrates(dof,iface_up,local_id_up) = &
1496) patch%internal_flow_fluxes(dof,sum_connection)
1497)
1498) idx = (local_id_up-1)*offset + (dof-1)*MAX_FACE_PER_CELL + iface_up + 1
1499) vec_ptr(idx) = patch%internal_flow_fluxes(dof,sum_connection)
1500)
1501) if (iface_dn>0) then
1502) ! Save flowrate for iface_dn of local_id_dn cell using -ve flowrate up-->dn
1503) flowrates(dof,iface_dn,local_id_dn) = -patch%internal_flow_fluxes(dof,sum_connection)
1504)
1505) idx = (local_id_dn-1)*offset + (dof-1)*MAX_FACE_PER_CELL + iface_dn + 1
1506) vec_ptr(idx) = -patch%internal_flow_fluxes(dof,sum_connection)
1507) endif
1508)
1509) enddo ! dof-loop
1510)
1511) enddo ! iconn-loop
1512)
1513) cur_connection_set => cur_connection_set%next
1514) enddo
1515)
1516) ! Boundary Flowrates Terms -----------------------------------
1517) boundary_condition => patch%boundary_condition_list%first
1518) sum_connection = 0
1519) do
1520) if (.not.associated(boundary_condition)) exit
1521)
1522) cur_connection_set => boundary_condition%connection_set
1523)
1524) do iconn = 1, cur_connection_set%num_connections
1525) sum_connection = sum_connection + 1
1526) face_id = cur_connection_set%face_id(iconn)
1527) ghosted_id_dn = cur_connection_set%id_dn(iconn)
1528) local_id_dn = grid%nG2L(ghosted_id_dn)
1529) do iface_dn = 1,MAX_FACE_PER_CELL
1530) if (face_id==ugrid%cell_to_face_ghosted(iface_dn,local_id_dn)) exit
1531) enddo
1532)
1533) do dof=1,option%nflowdof
1534)
1535) ! Save flowrate for iface_dn of local_id_dn cell using -ve flowrate up-->dn
1536) idx = (local_id_dn-1)*offset + (dof-1)*MAX_FACE_PER_CELL + iface_dn + 1
1537) flowrates(dof,iface_dn,local_id_dn) = &
1538) -patch%boundary_flow_fluxes(dof,sum_connection)
1539) vec_ptr(idx) = &
1540) -patch%boundary_flow_fluxes(dof,sum_connection)
1541)
1542) enddo ! dof-loop
1543)
1544) enddo ! iconn-loop
1545)
1546) boundary_condition => boundary_condition%next
1547) enddo
1548)
1549) deallocate(flowrates)
1550) call VecRestoreArrayF90(field%flowrate_inst,vec_ptr,ierr);CHKERRQ(ierr)
1551)
1552) ! Scatter flowrate from Global --> Natural order
1553) call VecScatterBegin(ugdm%scatter_gton,field%flowrate_inst,natural_flowrates_vec, &
1554) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1555) call VecScatterEnd(ugdm%scatter_gton,field%flowrate_inst,natural_flowrates_vec, &
1556) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1557)
1558) call VecGetArrayF90(natural_flowrates_vec,vec_ptr,ierr);CHKERRQ(ierr)
1559) call VecGetArrayF90(field%flowrate_inst,vec_ptr2,ierr);CHKERRQ(ierr)
1560)
1561) ! Copy the vectors
1562) vec_ptr2 = vec_ptr
1563)
1564) if (output_option%print_hdf5_aveg_mass_flowrate.or. &
1565) output_option%print_hdf5_aveg_energy_flowrate) then
1566)
1567) dtime = option%time-output_option%aveg_var_time
1568) call VecGetArrayF90(field%flowrate_aveg,vec_ptr3,ierr);CHKERRQ(ierr)
1569) vec_ptr3 = vec_ptr3 + vec_ptr2/dtime
1570) call VecRestoreArrayF90(field%flowrate_aveg,vec_ptr3,ierr);CHKERRQ(ierr)
1571) endif
1572)
1573) call VecRestoreArrayF90(natural_flowrates_vec,vec_ptr,ierr);CHKERRQ(ierr)
1574) call VecRestoreArrayF90(field%flowrate_inst,vec_ptr2,ierr);CHKERRQ(ierr)
1575)
1576) call VecDestroy(natural_flowrates_vec,ierr);CHKERRQ(ierr)
1577)
1578) call UGridDMDestroy(ugdm)
1579)
1580) end subroutine OutputGetFaceFlowrateUGrid
1581)
1582) ! ************************************************************************** !
1583)
1584) subroutine OutputGetExplicitIDsFlowrates(realization_base,count,vec_proc, &
1585) ids_up,ids_dn)
1586) !
1587) ! Calculates the ids of the nodes of a
1588) ! connection for flow rats output
1589) !
1590) ! Author: Satish Karra, LANL
1591) ! Date: 04/24/13
1592) !
1593)
1594) use Realization_Base_class, only : realization_base_type
1595) use Patch_module
1596) use Grid_module
1597) use Option_module
1598) use Grid_Unstructured_Aux_module
1599) use Field_module
1600) use Connection_module
1601)
1602) implicit none
1603)
1604) #include "petsc/finclude/petscvec.h"
1605) #include "petsc/finclude/petscvec.h90"
1606) #include "petsc/finclude/petsclog.h"
1607) #include "petsc/finclude/petscsys.h"
1608)
1609) class(realization_base_type) :: realization_base
1610) type(option_type), pointer :: option
1611) type(patch_type), pointer :: patch
1612) type(grid_type), pointer :: grid
1613) type(grid_unstructured_type),pointer :: ugrid
1614) type(field_type), pointer :: field
1615) type(ugdm_type), pointer :: ugdm
1616) type(connection_set_list_type), pointer :: connection_set_list
1617) type(connection_set_type), pointer :: cur_connection_set
1618)
1619)
1620) PetscReal, pointer :: vec_ptr(:)
1621) PetscReal, pointer :: vec_ptr2(:)
1622) PetscReal, pointer :: vec_proc_ptr(:)
1623) PetscInt, pointer :: ids_up(:),ids_dn(:)
1624) PetscInt :: offset
1625) PetscInt :: istart,iend
1626) PetscInt :: iconn
1627) PetscErrorCode :: ierr
1628) PetscReal :: val
1629) PetscInt :: ghosted_id_up, ghosted_id_dn
1630) PetscInt :: local_id_up, local_id_dn
1631) PetscReal :: proc_up, proc_dn, conn_proc
1632) PetscInt :: sum_connection, count
1633) Vec :: global_vec
1634) Vec :: local_vec
1635) Vec :: vec_proc
1636) PetscInt :: idof
1637)
1638) patch => realization_base%patch
1639) grid => patch%grid
1640) ugrid => grid%unstructured_grid
1641) option => realization_base%option
1642) field => realization_base%field
1643)
1644) call VecCreateMPI(option%mycomm, &
1645) size(grid%unstructured_grid%explicit_grid%connections,2), &
1646) PETSC_DETERMINE,vec_proc,ierr);CHKERRQ(ierr)
1647) call VecSet(vec_proc,0.d0,ierr);CHKERRQ(ierr)
1648)
1649) call UGridCreateUGDM(grid%unstructured_grid,ugdm,ONE_INTEGER,option)
1650) call UGridDMCreateVector(grid%unstructured_grid,ugdm,global_vec, &
1651) GLOBAL,option)
1652) call UGridDMCreateVector(grid%unstructured_grid,ugdm,local_vec, &
1653) LOCAL,option)
1654) call VecGetArrayF90(global_vec,vec_ptr,ierr);CHKERRQ(ierr)
1655) vec_ptr = option%myrank
1656) call VecRestoreArrayF90(global_vec,vec_ptr,ierr);CHKERRQ(ierr)
1657)
1658) call VecScatterBegin(ugdm%scatter_gtol,global_vec,local_vec, &
1659) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1660) call VecScatterEnd(ugdm%scatter_gtol,global_vec,local_vec, &
1661) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
1662)
1663) call VecGetArrayF90(local_vec,vec_ptr2,ierr);CHKERRQ(ierr)
1664) call VecGetArrayF90(vec_proc,vec_proc_ptr,ierr);CHKERRQ(ierr)
1665)
1666) connection_set_list => grid%internal_connection_set_list
1667) cur_connection_set => connection_set_list%first
1668) sum_connection = 0
1669) count = 0
1670) do
1671) if (.not.associated(cur_connection_set)) exit
1672) do iconn = 1, cur_connection_set%num_connections
1673) sum_connection = sum_connection + 1
1674) ghosted_id_up = cur_connection_set%id_up(iconn)
1675) ghosted_id_dn = cur_connection_set%id_dn(iconn)
1676) local_id_up = grid%nG2L(ghosted_id_up)
1677) local_id_dn = grid%nG2L(ghosted_id_dn)
1678) proc_up = vec_ptr2(ghosted_id_up)
1679) proc_dn = vec_ptr2(ghosted_id_dn)
1680) proc_up = min(option%myrank,int(proc_up))
1681) proc_dn = min(option%myrank,int(proc_dn))
1682) conn_proc = min(proc_up,proc_dn)
1683) vec_proc_ptr(sum_connection) = conn_proc
1684) enddo
1685) cur_connection_set => cur_connection_set%next
1686) enddo
1687)
1688) call VecRestoreArrayF90(local_vec,vec_ptr2,ierr);CHKERRQ(ierr)
1689) call VecRestoreArrayF90(vec_proc,vec_proc_ptr,ierr);CHKERRQ(ierr)
1690)
1691) call VecAssemblyBegin(vec_proc,ierr);CHKERRQ(ierr)
1692) call VecAssemblyEnd(vec_proc,ierr);CHKERRQ(ierr)
1693)
1694) ! Count the number of connections on a local process
1695) call VecGetArrayF90(vec_proc,vec_proc_ptr,ierr);CHKERRQ(ierr)
1696) cur_connection_set => connection_set_list%first
1697) sum_connection = 0
1698) count = 0
1699) do
1700) if (.not.associated(cur_connection_set)) exit
1701) do iconn = 1, cur_connection_set%num_connections
1702) sum_connection = sum_connection + 1
1703) if (option%myrank == int(vec_proc_ptr(sum_connection))) &
1704) count = count + 1
1705) enddo
1706) cur_connection_set => cur_connection_set%next
1707) enddo
1708) call VecRestoreArrayF90(vec_proc,vec_proc_ptr,ierr);CHKERRQ(ierr)
1709)
1710)
1711) ! Count the number of connections on a local process
1712) allocate(ids_up(count))
1713) allocate(ids_dn(count))
1714) call VecGetArrayF90(vec_proc,vec_proc_ptr,ierr);CHKERRQ(ierr)
1715) connection_set_list => grid%internal_connection_set_list
1716) cur_connection_set => connection_set_list%first
1717) sum_connection = 0
1718) count = 0
1719) do
1720) if (.not.associated(cur_connection_set)) exit
1721) do iconn = 1, cur_connection_set%num_connections
1722) sum_connection = sum_connection + 1
1723) ghosted_id_up = cur_connection_set%id_up(iconn)
1724) ghosted_id_dn = cur_connection_set%id_dn(iconn)
1725) local_id_up = grid%nG2L(ghosted_id_up)
1726) local_id_dn = grid%nG2L(ghosted_id_dn)
1727) if (option%myrank == int(vec_proc_ptr(sum_connection))) then
1728) count = count + 1
1729) ids_up(count) = grid%nG2A(ghosted_id_up)
1730) ids_dn(count) = grid%nG2A(ghosted_id_dn)
1731) endif
1732) enddo
1733) cur_connection_set => cur_connection_set%next
1734) enddo
1735) call VecRestoreArrayF90(vec_proc,vec_proc_ptr,ierr);CHKERRQ(ierr)
1736)
1737) end subroutine OutputGetExplicitIDsFlowrates
1738)
1739) ! ************************************************************************** !
1740)
1741) subroutine OutputGetExplicitFlowrates(realization_base,count,vec_proc, &
1742) flowrates,darcy,area)
1743) !
1744) ! Forms a vector of magnitude of flowrates
1745) ! which will be printed out to file for particle tracking.
1746) !
1747) ! Author: Satish Karra, LANL
1748) ! Date: 04/24/13
1749) !
1750)
1751) use Realization_Base_class, only : realization_base_type
1752) use Patch_module
1753) use Grid_module
1754) use Option_module
1755) use Grid_Unstructured_Aux_module
1756) use Field_module
1757) use Connection_module
1758)
1759) implicit none
1760)
1761) #include "petsc/finclude/petscvec.h"
1762) #include "petsc/finclude/petscvec.h90"
1763) #include "petsc/finclude/petsclog.h"
1764) #include "petsc/finclude/petscsys.h"
1765)
1766) class(realization_base_type) :: realization_base
1767) type(option_type), pointer :: option
1768) type(patch_type), pointer :: patch
1769) type(grid_type), pointer :: grid
1770) type(grid_unstructured_type),pointer :: ugrid
1771) type(field_type), pointer :: field
1772) type(connection_set_list_type), pointer :: connection_set_list
1773) type(connection_set_type), pointer :: cur_connection_set
1774)
1775)
1776) PetscReal, pointer :: vec_proc_ptr(:)
1777) PetscReal, pointer :: flowrates(:,:)
1778) PetscReal, pointer :: darcy(:), area(:)
1779) PetscInt :: offset
1780) PetscInt :: iconn
1781) PetscErrorCode :: ierr
1782) PetscInt :: ghosted_id_up, ghosted_id_dn
1783) PetscInt :: local_id_up, local_id_dn
1784) PetscReal :: proc_up, proc_dn, conn_proc
1785) PetscInt :: sum_connection, count
1786) Vec :: vec_proc
1787) PetscInt :: idof
1788)
1789) patch => realization_base%patch
1790) grid => patch%grid
1791) ugrid => grid%unstructured_grid
1792) option => realization_base%option
1793) field => realization_base%field
1794)
1795) ! Count the number of connections on a local process
1796) allocate(flowrates(count,option%nflowdof))
1797) allocate(darcy(count))
1798) allocate(area(count))
1799) call VecGetArrayF90(vec_proc,vec_proc_ptr,ierr);CHKERRQ(ierr)
1800) connection_set_list => grid%internal_connection_set_list
1801) cur_connection_set => connection_set_list%first
1802) sum_connection = 0
1803) count = 0
1804) do
1805) if (.not.associated(cur_connection_set)) exit
1806) do iconn = 1, cur_connection_set%num_connections
1807) sum_connection = sum_connection + 1
1808) ghosted_id_up = cur_connection_set%id_up(iconn)
1809) ghosted_id_dn = cur_connection_set%id_dn(iconn)
1810) local_id_up = grid%nG2L(ghosted_id_up)
1811) local_id_dn = grid%nG2L(ghosted_id_dn)
1812) if (option%myrank == int(vec_proc_ptr(sum_connection))) then
1813) count = count + 1
1814) do idof = 1,option%nflowdof
1815) flowrates(count,option%nflowdof) = &
1816) patch%internal_flow_fluxes(idof,sum_connection)
1817) enddo
1818) darcy(count) = patch%internal_velocities(1,sum_connection)
1819) area(count) = cur_connection_set%area(iconn)
1820) endif
1821) enddo
1822) cur_connection_set => cur_connection_set%next
1823) enddo
1824) call VecRestoreArrayF90(vec_proc,vec_proc_ptr,ierr);CHKERRQ(ierr)
1825)
1826) end subroutine OutputGetExplicitFlowrates
1827)
1828) ! ************************************************************************** !
1829)
1830) subroutine OutputGetExplicitAuxVars(realization_base,count,vec_proc,density)
1831) !
1832) ! Calculates density at the face
1833) ! between a connection
1834) !
1835) ! Author: Satish Karra, LANL
1836) ! Date: 07/17/13
1837) !
1838)
1839) use Realization_Base_class, only : realization_base_type
1840) use Patch_module
1841) use Grid_module
1842) use Option_module
1843) use Grid_Unstructured_Aux_module
1844) use Field_module
1845) use Connection_module
1846) use Global_Aux_module
1847) use Richards_Aux_module
1848) use Material_Aux_class
1849)
1850) implicit none
1851)
1852) #include "petsc/finclude/petscvec.h"
1853) #include "petsc/finclude/petscvec.h90"
1854) #include "petsc/finclude/petsclog.h"
1855) #include "petsc/finclude/petscsys.h"
1856)
1857) class(realization_base_type) :: realization_base
1858) type(option_type), pointer :: option
1859) type(patch_type), pointer :: patch
1860) type(grid_type), pointer :: grid
1861) type(grid_unstructured_type),pointer :: ugrid
1862) type(field_type), pointer :: field
1863) type(connection_set_list_type), pointer :: connection_set_list
1864) type(connection_set_type), pointer :: cur_connection_set
1865) type(global_auxvar_type), pointer :: global_auxvar(:)
1866) type(material_parameter_type), pointer :: material_parameter
1867)
1868)
1869) PetscReal, pointer :: vec_proc_ptr(:)
1870) PetscReal, pointer :: flowrates(:,:)
1871) PetscReal, pointer :: darcy(:)
1872) PetscReal, pointer :: density(:)
1873) PetscInt :: offset
1874) PetscInt :: iconn
1875) PetscErrorCode :: ierr
1876) PetscInt :: ghosted_id_up, ghosted_id_dn
1877) PetscInt :: local_id_up, local_id_dn
1878) PetscReal :: proc_up, proc_dn, conn_proc
1879) PetscInt :: sum_connection, count
1880) Vec :: vec_proc
1881) PetscInt :: idof
1882) PetscInt :: icap_up, icap_dn
1883) PetscReal :: sir_up, sir_dn
1884) PetscReal, parameter :: eps = 1.D-8
1885) PetscReal :: upweight
1886)
1887)
1888) patch => realization_base%patch
1889) option => realization_base%option
1890) field => realization_base%field
1891) grid => patch%grid
1892) global_auxvar => patch%aux%Global%auxvars
1893) material_parameter => patch%aux%Material%material_parameter
1894)
1895) allocate(density(count))
1896) call VecGetArrayF90(vec_proc,vec_proc_ptr,ierr);CHKERRQ(ierr)
1897) connection_set_list => grid%internal_connection_set_list
1898) cur_connection_set => connection_set_list%first
1899) sum_connection = 0
1900) count = 0
1901) do
1902) if (.not.associated(cur_connection_set)) exit
1903) do iconn = 1, cur_connection_set%num_connections
1904) sum_connection = sum_connection + 1
1905) ghosted_id_up = cur_connection_set%id_up(iconn)
1906) ghosted_id_dn = cur_connection_set%id_dn(iconn)
1907) local_id_up = grid%nG2L(ghosted_id_up)
1908) local_id_dn = grid%nG2L(ghosted_id_dn)
1909) icap_up = patch%sat_func_id(ghosted_id_up)
1910) icap_dn = patch%sat_func_id(ghosted_id_dn)
1911) if (option%myrank == int(vec_proc_ptr(sum_connection))) then
1912) count = count + 1
1913) sir_up = material_parameter%soil_residual_saturation(1,icap_up)
1914) sir_dn = material_parameter%soil_residual_saturation(1,icap_dn)
1915)
1916) if (global_auxvar(ghosted_id_up)%sat(1) > sir_up .or. &
1917) global_auxvar(ghosted_id_dn)%sat(1) > sir_dn) then
1918) if (global_auxvar(ghosted_id_up)%sat(1) <eps) then
1919) upweight = 0.d0
1920) else if (global_auxvar(ghosted_id_dn)%sat(1) <eps) then
1921) upweight = 1.d0
1922) endif
1923)
1924) density(count) = upweight*global_auxvar(ghosted_id_up)%den(1)+ &
1925) (1.D0 - upweight)*global_auxvar(ghosted_id_dn)%den(1)
1926) endif
1927) endif
1928)
1929) enddo
1930) cur_connection_set => cur_connection_set%next
1931) enddo
1932)
1933) call VecRestoreArrayF90(vec_proc,vec_proc_ptr,ierr);CHKERRQ(ierr)
1934)
1935)
1936) end subroutine OutputGetExplicitAuxVars
1937)
1938) ! ************************************************************************** !
1939)
1940) subroutine OutputGetExplicitCellInfo(realization_base,num_cells,ids,sat,por, &
1941) density,pressure)
1942) !
1943) ! Calculates porosity, saturation, density
1944) ! and pressure in a cell (explicit)
1945) !
1946) ! Author: Satish Karra, LANL
1947) ! Date: 08/21/13
1948) !
1949)
1950) use Realization_Base_class, only : realization_base_type
1951) use Patch_module
1952) use Grid_module
1953) use Option_module
1954) use Grid_Unstructured_Aux_module
1955) use Field_module
1956) use Connection_module
1957) use Global_Aux_module
1958)
1959) implicit none
1960)
1961) #include "petsc/finclude/petscvec.h"
1962) #include "petsc/finclude/petscvec.h90"
1963) #include "petsc/finclude/petsclog.h"
1964) #include "petsc/finclude/petscsys.h"
1965)
1966) class(realization_base_type) :: realization_base
1967) type(option_type), pointer :: option
1968) type(patch_type), pointer :: patch
1969) type(grid_type), pointer :: grid
1970) type(grid_unstructured_type),pointer :: ugrid
1971) type(field_type), pointer :: field
1972) type(global_auxvar_type), pointer :: global_auxvar(:)
1973)
1974)
1975) PetscErrorCode :: ierr
1976) PetscInt :: num_cells
1977) PetscReal, pointer :: sat(:)
1978) PetscReal, pointer :: por(:)
1979) PetscReal, pointer :: density(:)
1980) PetscReal, pointer :: pressure(:)
1981) PetscInt, pointer :: ids(:)
1982) PetscInt :: local_id, ghosted_id
1983)
1984) patch => realization_base%patch
1985) option => realization_base%option
1986) field => realization_base%field
1987) grid => patch%grid
1988) global_auxvar => patch%aux%Global%auxvars
1989)
1990) num_cells = grid%nlmax
1991) allocate(sat(num_cells))
1992) allocate(por(num_cells))
1993) allocate(ids(num_cells))
1994) allocate(density(num_cells))
1995) allocate(pressure(num_cells))
1996)
1997) do local_id = 1, num_cells
1998) ghosted_id = grid%nL2G(local_id)
1999) ids(local_id) = grid%nG2A(ghosted_id)
2000) sat(local_id) = global_auxvar(ghosted_id)%sat(1)
2001) por(local_id) = patch%aux%Material%auxvars(ghosted_id)%porosity
2002) density(local_id) = global_auxvar(ghosted_id)%den(1)
2003) pressure(local_id) = global_auxvar(ghosted_id)%pres(1)
2004) enddo
2005)
2006) end subroutine OutputGetExplicitCellInfo
2007)
2008) end module Output_Common_module