output_observation.F90 coverage: 59.09 %func 23.68 %block
1) module Output_Observation_module
2)
3) use Logging_module
4) use Output_Aux_module
5) use Output_Common_module
6)
7) use PFLOTRAN_Constants_module
8)
9) implicit none
10)
11) private
12)
13) #include "petsc/finclude/petscsys.h"
14)
15) ! flags signifying the first time a routine is called during a given
16) ! simulation
17) PetscBool :: observation_first
18) PetscBool :: check_for_obs_points
19) PetscBool :: calculate_velocities ! true if any obs. pt. prints velocity
20) PetscBool :: secondary_observation_first
21) PetscBool :: secondary_check_for_obs_points
22) PetscBool :: mass_balance_first
23) PetscBool :: integral_flux_first
24)
25) public :: OutputObservation, &
26) OutputObservationInit, &
27) OutputMassBalance, &
28) OutputIntegralFlux
29)
30) contains
31)
32) ! ************************************************************************** !
33)
34) subroutine OutputObservationInit(num_steps)
35) !
36) ! Initializes module variables for observation
37) !
38) ! Author: Glenn Hammond
39) ! Date: 01/16/13
40) !
41)
42) use Option_module
43)
44) implicit none
45)
46) PetscInt :: num_steps
47)
48) check_for_obs_points = PETSC_TRUE
49) calculate_velocities = PETSC_FALSE
50) secondary_check_for_obs_points = PETSC_TRUE
51) if (num_steps == 0) then
52) observation_first = PETSC_TRUE
53) secondary_observation_first = PETSC_TRUE
54) mass_balance_first = PETSC_TRUE
55) integral_flux_first = PETSC_TRUE
56) else
57) observation_first = PETSC_FALSE
58) secondary_observation_first = PETSC_FALSE
59) mass_balance_first = PETSC_FALSE
60) integral_flux_first = PETSC_FALSE
61) endif
62)
63) end subroutine OutputObservationInit
64)
65) ! ************************************************************************** !
66)
67) subroutine OutputObservation(realization_base)
68) !
69) ! Main driver for all observation output subroutines
70) !
71) ! Author: Glenn Hammond
72) ! Date: 02/11/08
73) !
74)
75) use Realization_Base_class, only : realization_base_type
76) use Option_module
77)
78) implicit none
79)
80) class(realization_base_type) :: realization_base
81)
82) ! if (realization_base%output_option%print_hdf5) then
83) ! call OutputObservationHDF5(realization)
84) ! call OutputObservationTecplot(realization)
85) ! endif
86)
87) ! if (realization_base%output_option%print_tecplot .or. &
88) ! realization_base%output_option%print_hdf5) then
89) if (realization_base%output_option%print_observation) then
90) call OutputObservationTecplotColumnTXT(realization_base)
91) call OutputIntegralFlux(realization_base)
92) if (realization_base%option%use_mc) then
93) call OutputObservationTecplotSecTXT(realization_base)
94) endif
95) endif
96)
97) end subroutine OutputObservation
98)
99) ! ************************************************************************** !
100)
101) subroutine OutputObservationTecplotColumnTXT(realization_base)
102) !
103) ! Print to observation data to text file
104) !
105) ! Author: Glenn Hammond
106) ! Date: 02/11/08
107) !
108)
109) use Realization_Base_class, only : realization_base_type
110) use Discretization_module
111) use Grid_module
112) use Option_module
113) use Field_module
114) use Patch_module
115) use Observation_module
116) use Utility_module
117)
118) implicit none
119)
120) class(realization_base_type) :: realization_base
121)
122) PetscInt :: fid, icell
123) character(len=MAXSTRINGLENGTH) :: filename
124) character(len=MAXSTRINGLENGTH) :: string, string2
125) type(grid_type), pointer :: grid
126) type(option_type), pointer :: option
127) type(field_type), pointer :: field
128) type(patch_type), pointer :: patch
129) type(output_option_type), pointer :: output_option
130) type(observation_type), pointer :: observation
131) PetscBool, save :: open_file = PETSC_FALSE
132) PetscReal, allocatable :: velocities(:,:,:)
133) PetscInt :: local_id
134) PetscInt :: icolumn
135) PetscErrorCode :: ierr
136)
137) call PetscLogEventBegin(logging%event_output_observation,ierr);CHKERRQ(ierr)
138)
139) patch => realization_base%patch
140) grid => patch%grid
141) option => realization_base%option
142) field => realization_base%field
143) output_option => realization_base%output_option
144)
145) if (check_for_obs_points) then
146) open_file = PETSC_FALSE
147) observation => patch%observation_list%first
148) do
149) if (.not.associated(observation)) exit
150) if (observation%print_velocities) calculate_velocities = PETSC_TRUE
151) if (observation%itype == OBSERVATION_SCALAR .or. &
152) (observation%itype == OBSERVATION_FLUX .and. &
153) option%myrank == option%io_rank)) then
154) open_file = PETSC_TRUE
155) exit
156) endif
157) observation => observation%next
158) enddo
159) check_for_obs_points = PETSC_FALSE
160) endif
161)
162) if (calculate_velocities) then
163) allocate(velocities(3,realization_base%patch%grid%nlmax,option%nphase))
164) call PatchGetCellCenteredVelocities(realization_base%patch, &
165) ONE_INTEGER,velocities(:,:,1))
166) if (option%nphase > 1) then
167) call PatchGetCellCenteredVelocities(realization_base%patch, &
168) TWO_INTEGER,velocities(:,:,2))
169) endif
170) endif
171)
172) if (open_file) then
173) write(string,'(i6)') option%myrank
174) filename = trim(option%global_prefix) // trim(option%group_prefix) // &
175) '-obs-' // trim(adjustl(string)) // '.tec'
176)
177) ! open file
178) fid = 86
179) if (observation_first .or. .not.FileExists(filename)) then
180) open(unit=fid,file=filename,action="write",status="replace")
181) ! write header
182) ! write title
183) write(fid,'(a)',advance="no") ' "Time [' // trim(output_option%tunit) // &
184) ']"'
185) observation => patch%observation_list%first
186)
187) ! must initialize icolumn here so that icolumn does not restart with
188) ! each observation point
189) if (output_option%print_column_ids) then
190) icolumn = 1
191) else
192) icolumn = -1
193) endif
194)
195) do
196) if (.not.associated(observation)) exit
197)
198) select case(observation%itype)
199) case(OBSERVATION_SCALAR)
200) if (associated(observation%region%coordinates) .and. &
201) .not.observation%at_cell_center) then
202) ! option%io_buffer = 'Writing of data at coordinates not ' // &
203) ! 'functioning properly for minerals. Perhaps due to ' // &
204) ! 'non-ghosting of vol frac....>? - geh'
205) ! call printErrMsg(option)
206) call WriteObservationHeaderForCoord(fid,realization_base, &
207) observation%region, &
208) observation%print_velocities, &
209) icolumn)
210) else
211) do icell=1,observation%region%num_cells
212) call WriteObservationHeaderForCell(fid,realization_base, &
213) observation%region,icell, &
214) observation%print_velocities, &
215) icolumn)
216) enddo
217) endif
218) case(OBSERVATION_FLUX)
219) if (option%myrank == option%io_rank) then
220) call WriteObservationHeaderForBC(fid,realization_base, &
221) observation%linkage_name)
222) endif
223) end select
224) observation => observation%next
225) enddo
226) write(fid,'(a)',advance="yes") ""
227) else
228) open(unit=fid,file=filename,action="write",status="old", &
229) position="append")
230) endif
231)
232) observation => patch%observation_list%first
233) write(fid,'(1es14.6)',advance="no") option%time/output_option%tconv
234) do
235) if (.not.associated(observation)) exit
236) select case(observation%itype)
237) case(OBSERVATION_SCALAR)
238) if (associated(observation%region%coordinates) .and. &
239) .not.observation%at_cell_center) then
240) call WriteObservationDataForCoord(fid,realization_base, &
241) observation%region)
242) if (observation%print_velocities) then
243) call WriteVelocityAtCoord(fid,realization_base, &
244) observation%region)
245) endif
246) else
247) do icell=1,observation%region%num_cells
248) local_id = observation%region%cell_ids(icell)
249) call WriteObservationDataForCell(fid,realization_base,local_id)
250) if (observation%print_velocities) then
251) call WriteVelocityAtCell2(fid,realization_base,local_id, &
252) velocities)
253) endif
254) enddo
255) endif
256) case(OBSERVATION_FLUX)
257) call WriteObservationDataForBC(fid,realization_base, &
258) patch, &
259) observation%connection_set)
260) end select
261) observation => observation%next
262) enddo
263) write(fid,'(a)',advance="yes") ""
264) close(fid)
265)
266) endif
267)
268) observation_first = PETSC_FALSE
269) if (allocated(velocities)) deallocate(velocities)
270)
271) call PetscLogEventEnd(logging%event_output_observation,ierr);CHKERRQ(ierr)
272)
273) end subroutine OutputObservationTecplotColumnTXT
274)
275) ! ************************************************************************** !
276)
277) subroutine WriteObservationHeaderForCell(fid,realization_base,region,icell, &
278) print_velocities, &
279) icolumn)
280) !
281) ! Print a header for data at a cell
282) !
283) ! Author: Glenn Hammond
284) ! Date: 02/11/08
285) !
286)
287) use Realization_Base_class, only : realization_base_type
288) use Grid_module
289) use Option_module
290) use Output_Aux_module
291) use Patch_module
292) use Region_module
293) use Utility_module, only : BestFloat
294)
295) implicit none
296)
297) PetscInt :: fid
298) class(realization_base_type) :: realization_base
299) type(region_type) :: region
300) PetscInt :: icell
301) PetscBool :: print_velocities
302) PetscInt :: icolumn
303)
304) PetscInt :: local_id
305) character(len=MAXSTRINGLENGTH) :: cell_string
306) character(len=MAXWORDLENGTH) :: x_string, y_string, z_string
307) type(grid_type), pointer :: grid
308)
309) grid => realization_base%patch%grid
310)
311) local_id = region%cell_ids(icell)
312) write(cell_string,*) grid%nG2A(grid%nL2G(region%cell_ids(icell)))
313) cell_string = trim(region%name) // ' (' // trim(adjustl(cell_string)) // ')'
314)
315) ! add coordinate of cell center
316) x_string = BestFloat(grid%x(grid%nL2G(local_id)),1.d4,1.d-2)
317) y_string = BestFloat(grid%y(grid%nL2G(local_id)),1.d4,1.d-2)
318) z_string = BestFloat(grid%z(grid%nL2G(local_id)),1.d4,1.d-2)
319) cell_string = trim(cell_string) // ' (' // trim(adjustl(x_string)) // &
320) ' ' // trim(adjustl(y_string)) // &
321) ' ' // trim(adjustl(z_string)) // ')'
322)
323) call WriteObservationHeader(fid,realization_base,cell_string, &
324) print_velocities,icolumn)
325)
326) end subroutine WriteObservationHeaderForCell
327)
328) ! ************************************************************************** !
329)
330) subroutine WriteObservationHeaderForCoord(fid,realization_base,region, &
331) print_velocities, &
332) icolumn)
333) !
334) ! Print a header for data at a coordinate
335) !
336) ! Author: Glenn Hammond
337) ! Date: 04/11/08
338) !
339)
340) use Realization_Base_class, only : realization_base_type
341) use Option_module
342) use Patch_module
343) use Region_module
344) use Utility_module, only : BestFloat
345)
346) implicit none
347)
348) PetscInt :: fid
349) class(realization_base_type) :: realization_base
350) type(region_type) :: region
351) PetscBool :: print_velocities
352) PetscInt :: icolumn
353)
354) character(len=MAXSTRINGLENGTH) :: cell_string
355) character(len=MAXWORDLENGTH) :: x_string, y_string, z_string
356)
357) cell_string = trim(region%name)
358)
359) x_string = BestFloat(region%coordinates(ONE_INTEGER)%x,1.d4,1.d-2)
360) y_string = BestFloat(region%coordinates(ONE_INTEGER)%y,1.d4,1.d-2)
361) z_string = BestFloat(region%coordinates(ONE_INTEGER)%z,1.d4,1.d-2)
362) cell_string = trim(cell_string) // ' (' // trim(adjustl(x_string)) // ' ' // &
363) trim(adjustl(y_string)) // ' ' // &
364) trim(adjustl(z_string)) // ')'
365)
366) call WriteObservationHeader(fid,realization_base,cell_string, &
367) print_velocities,icolumn)
368)
369) end subroutine WriteObservationHeaderForCoord
370)
371) ! ************************************************************************** !
372)
373) subroutine WriteObservationHeader(fid,realization_base,cell_string, &
374) print_velocities,icolumn)
375) !
376) ! Print a header for data
377) !
378) ! Author: Glenn Hammond
379) ! Date: 10/27/11
380) !
381)
382) use Realization_Base_class, only : realization_base_type
383) use Option_module
384) use Reaction_Aux_module
385)
386) implicit none
387)
388) PetscInt :: fid
389) class(realization_base_type) :: realization_base
390) type(reaction_type), pointer :: reaction
391) PetscBool :: print_velocities
392) character(len=MAXSTRINGLENGTH) :: cell_string
393) PetscInt :: icolumn
394)
395) PetscInt :: variable_count
396) character(len=MAXSTRINGLENGTH) :: string
397) type(option_type), pointer :: option
398) type(output_option_type), pointer :: output_option
399)
400) option => realization_base%option
401) output_option => realization_base%output_option
402)
403) call OutputWriteVariableListToHeader(fid, &
404) output_option%output_obs_variable_list, &
405) cell_string,icolumn,PETSC_FALSE, &
406) variable_count)
407)
408) if (print_velocities) then
409) write(string,'(''m/'',a,'' '')') trim(realization_base%output_option%tunit)
410) call OutputWriteToHeader(fid,'qlx',string,cell_string,icolumn)
411) call OutputWriteToHeader(fid,'qly',string,cell_string,icolumn)
412) call OutputWriteToHeader(fid,'qlz',string,cell_string,icolumn)
413)
414) if (option%nphase > 1) then
415) call OutputWriteToHeader(fid,'qgx',string,cell_string,icolumn)
416) call OutputWriteToHeader(fid,'qgy',string,cell_string,icolumn)
417) call OutputWriteToHeader(fid,'qgz',string,cell_string,icolumn)
418) endif
419) endif
420)
421) end subroutine WriteObservationHeader
422)
423) ! ************************************************************************** !
424)
425) subroutine OutputObservationTecplotSecTXT(realization_base)
426) !
427) ! Print to secondary continuum observation
428) ! data to text file
429) !
430) ! Author: Satish Karra, LANL
431) ! Date: 04/08/13
432) !
433)
434) use Realization_Base_class, only : realization_base_type
435) use Discretization_module
436) use Grid_module
437) use Option_module
438) use Field_module
439) use Patch_module
440) use Observation_module
441) use Utility_module
442)
443) implicit none
444)
445) class(realization_base_type) :: realization_base
446)
447) PetscInt :: fid, icell
448) character(len=MAXSTRINGLENGTH) :: filename
449) character(len=MAXSTRINGLENGTH) :: string, string2
450) type(grid_type), pointer :: grid
451) type(option_type), pointer :: option
452) type(field_type), pointer :: field
453) type(patch_type), pointer :: patch
454) type(output_option_type), pointer :: output_option
455) type(observation_type), pointer :: observation
456) PetscBool, save :: open_file = PETSC_FALSE
457) PetscInt :: local_id
458) PetscInt :: icolumn
459) PetscErrorCode :: ierr
460)
461) call PetscLogEventBegin(logging%event_output_observation,ierr);CHKERRQ(ierr)
462)
463) patch => realization_base%patch
464) grid => patch%grid
465) option => realization_base%option
466) field => realization_base%field
467) output_option => realization_base%output_option
468)
469) if (secondary_check_for_obs_points) then
470) open_file = PETSC_FALSE
471) observation => patch%observation_list%first
472) do
473) if (.not.associated(observation)) exit
474) if (observation%itype == OBSERVATION_SCALAR .or. &
475) (observation%itype == OBSERVATION_FLUX .and. &
476) option%myrank == option%io_rank)) then
477) open_file = PETSC_TRUE
478) exit
479) endif
480) observation => observation%next
481) enddo
482) secondary_check_for_obs_points = PETSC_FALSE
483) endif
484)
485)
486) if (open_file) then
487) write(string,'(i6)') option%myrank
488) filename = trim(option%global_prefix) // trim(option%group_prefix) // &
489) '-obs-sec-' // trim(adjustl(string)) // '.tec'
490)
491) ! open file
492) fid = 86
493) if (secondary_observation_first .or. .not.FileExists(filename)) then
494) open(unit=fid,file=filename,action="write",status="replace")
495) ! write header
496) ! write title
497) write(fid,'(a)',advance="no") ' "Time [' // trim(output_option%tunit) // &
498) ']"'
499) observation => patch%observation_list%first
500)
501) ! must initialize icolumn here so that icolumn does not restart with
502) ! each observation point
503) if (output_option%print_column_ids) then
504) icolumn = 1
505) else
506) icolumn = -1
507) endif
508)
509) do
510) if (.not.associated(observation)) exit
511)
512) select case(observation%itype)
513) case(OBSERVATION_SCALAR)
514) if (associated(observation%region%coordinates) .and. &
515) .not.observation%at_cell_center) then
516) option%io_buffer = 'Writing of data at coordinates not ' // &
517) 'functioning properly for minerals. Perhaps due to ' // &
518) 'non-ghosting of vol frac....>? - geh'
519) call printErrMsg(option)
520) call WriteObservationHeaderForCoordSec(fid,realization_base, &
521) observation%region, &
522) observation% &
523) print_secondary_data, &
524) icolumn)
525) else
526) do icell=1,observation%region%num_cells
527) call WriteObservationHeaderForCellSec(fid,realization_base, &
528) observation%region,icell, &
529) observation% &
530) print_secondary_data, &
531) icolumn)
532) enddo
533) endif
534) end select
535) observation => observation%next
536) enddo
537) write(fid,'(a)',advance="yes") ""
538) else
539) open(unit=fid,file=filename,action="write",status="old", &
540) position="append")
541) endif
542)
543) observation => patch%observation_list%first
544) write(fid,'(1es14.6)',advance="no") option%time/output_option%tconv
545) do
546) if (.not.associated(observation)) exit
547) select case(observation%itype)
548) case(OBSERVATION_SCALAR)
549) do icell=1,observation%region%num_cells
550) local_id = observation%region%cell_ids(icell)
551) if (observation%print_secondary_data(1)) then
552) call WriteObservationSecondaryDataAtCell(fid, &
553) realization_base, &
554) local_id, &
555) PRINT_SEC_TEMP)
556) endif
557) if (observation%print_secondary_data(2)) then
558) call WriteObservationSecondaryDataAtCell(fid, &
559) realization_base, &
560) local_id, &
561) PRINT_SEC_CONC)
562) endif
563) if (observation%print_secondary_data(3)) then
564) call WriteObservationSecondaryDataAtCell(fid, &
565) realization_base, &
566) local_id, &
567) PRINT_SEC_MIN_VOLFRAC)
568) endif
569) if (observation%print_secondary_data(4)) then
570) call WriteObservationSecondaryDataAtCell(fid, &
571) realization_base, &
572) local_id, &
573) PRINT_SEC_MIN_RATE)
574) endif
575) if (observation%print_secondary_data(5)) then
576) call WriteObservationSecondaryDataAtCell(fid, &
577) realization_base, &
578) local_id, &
579) PRINT_SEC_MIN_SI)
580) endif
581) enddo
582) end select
583) observation => observation%next
584) enddo
585) write(fid,'(a)',advance="yes") ""
586) close(fid)
587)
588) endif
589)
590) secondary_observation_first = PETSC_FALSE
591)
592) call PetscLogEventEnd(logging%event_output_observation,ierr);CHKERRQ(ierr)
593)
594) end subroutine OutputObservationTecplotSecTXT
595)
596) ! ************************************************************************** !
597)
598) subroutine WriteObservationHeaderForCellSec(fid,realization_base,region,icell, &
599) print_secondary_data, &
600) icolumn)
601) !
602) ! Print a header for data at a cell for
603) ! secondary continuum
604) !
605) ! Author: Satish Karra, LANL
606) ! Date: 04/08/13
607) !
608)
609) use Realization_Base_class, only : realization_base_type
610) use Grid_module
611) use Option_module
612) use Output_Aux_module
613) use Patch_module
614) use Region_module
615) use Utility_module, only : BestFloat
616)
617) implicit none
618)
619) PetscInt :: fid
620) class(realization_base_type) :: realization_base
621) type(region_type) :: region
622) PetscInt :: icell
623) PetscBool :: print_secondary_data(5)
624) PetscInt :: icolumn
625)
626) PetscInt :: local_id
627) character(len=MAXSTRINGLENGTH) :: cell_string
628) character(len=MAXWORDLENGTH) :: x_string, y_string, z_string
629) type(grid_type), pointer :: grid
630)
631) grid => realization_base%patch%grid
632)
633) local_id = region%cell_ids(icell)
634) write(cell_string,*) grid%nG2A(grid%nL2G(region%cell_ids(icell)))
635) cell_string = trim(region%name) // ' (' // trim(adjustl(cell_string)) // ')'
636)
637) ! add coordinate of cell center
638) x_string = BestFloat(grid%x(grid%nL2G(local_id)),1.d4,1.d-2)
639) y_string = BestFloat(grid%y(grid%nL2G(local_id)),1.d4,1.d-2)
640) z_string = BestFloat(grid%z(grid%nL2G(local_id)),1.d4,1.d-2)
641) cell_string = trim(cell_string) // ' (' // trim(adjustl(x_string)) // &
642) ' ' // trim(adjustl(y_string)) // &
643) ' ' // trim(adjustl(z_string)) // ')'
644)
645) call WriteObservationHeaderSec(fid,realization_base,cell_string, &
646) print_secondary_data,icolumn)
647)
648) end subroutine WriteObservationHeaderForCellSec
649)
650) ! ************************************************************************** !
651)
652) subroutine WriteObservationHeaderForCoordSec(fid,realization_base,region, &
653) print_secondary_data, &
654) icolumn)
655) !
656) ! Print a header for data at a coordinate
657) ! for secondary continuum
658) !
659) ! Author: Satish Karra, LANL
660) ! Date: 04/08/13
661) !
662)
663) use Realization_Base_class, only : realization_base_type
664) use Option_module
665) use Patch_module
666) use Region_module
667) use Utility_module, only : BestFloat
668)
669) implicit none
670)
671) PetscInt :: fid
672) class(realization_base_type) :: realization_base
673) type(region_type) :: region
674) PetscBool :: print_secondary_data(5)
675) PetscInt :: icolumn
676)
677) character(len=MAXSTRINGLENGTH) :: cell_string
678) character(len=MAXWORDLENGTH) :: x_string, y_string, z_string
679)
680) cell_string = trim(region%name)
681)
682) x_string = BestFloat(region%coordinates(ONE_INTEGER)%x,1.d4,1.d-2)
683) y_string = BestFloat(region%coordinates(ONE_INTEGER)%y,1.d4,1.d-2)
684) z_string = BestFloat(region%coordinates(ONE_INTEGER)%z,1.d4,1.d-2)
685) cell_string = trim(cell_string) // ' (' // trim(adjustl(x_string)) // ' ' // &
686) trim(adjustl(y_string)) // ' ' // &
687) trim(adjustl(z_string)) // ')'
688)
689) call WriteObservationHeaderSec(fid,realization_base,cell_string, &
690) print_secondary_data,icolumn)
691)
692) end subroutine WriteObservationHeaderForCoordSec
693)
694) ! ************************************************************************** !
695)
696) subroutine WriteObservationHeaderSec(fid,realization_base,cell_string, &
697) print_secondary_data,icolumn)
698) !
699) ! Print a header for secondary continuum data
700) !
701) ! Author: Satish Karra, LANL
702) ! Date: 10/27/13
703) !
704)
705) use Realization_Base_class, only : realization_base_type
706) use Option_module
707) use Reaction_Aux_module
708)
709) implicit none
710)
711) PetscInt :: fid
712) class(realization_base_type) :: realization_base
713) type(reaction_type), pointer :: reaction
714) PetscBool :: print_secondary_data(5)
715) character(len=MAXSTRINGLENGTH) :: cell_string
716) PetscInt :: icolumn
717)
718) PetscInt :: i,j
719) character(len=MAXSTRINGLENGTH) :: string
720) type(option_type), pointer :: option
721) type(output_option_type), pointer :: output_option
722)
723) option => realization_base%option
724) output_option => realization_base%output_option
725)
726) ! add secondary temperature to header
727) if (print_secondary_data(1)) then
728) select case (option%iflowmode)
729) case (TH_MODE, MPH_MODE)
730) do i = 1, option%nsec_cells
731) write(string,'(i2)') i
732) string = 'T(' // trim(adjustl(string)) // ')'
733) call OutputWriteToHeader(fid,string,'C',cell_string,icolumn)
734) enddo
735) case default
736) end select
737) endif
738)
739) ! add secondary concentrations to header
740) if (option%ntrandof > 0) then
741) reaction => realization_base%reaction
742) if (print_secondary_data(2)) then
743) do j = 1, reaction%naqcomp
744) do i = 1, option%nsec_cells
745) write(string,'(i2)') i
746) string = 'C(' // trim(adjustl(string)) // ') ' &
747) // trim(reaction%primary_species_names(j))
748) call OutputWriteToHeader(fid,string,'molal',cell_string, &
749) icolumn)
750) enddo
751) enddo
752) endif
753)
754) ! add secondary mineral volume fractions to header
755) if (print_secondary_data(3)) then
756) do j = 1, reaction%mineral%nkinmnrl
757) do i = 1, option%nsec_cells
758) write(string,'(i2)') i
759) string = 'VF(' // trim(adjustl(string)) // ') ' &
760) // trim(reaction%mineral%mineral_names(j))
761) call OutputWriteToHeader(fid,string,'',cell_string,icolumn)
762) enddo
763) enddo
764) endif
765)
766) ! add secondary mineral rates to header
767) if (print_secondary_data(4)) then
768) do j = 1, reaction%mineral%nkinmnrl
769) do i = 1, option%nsec_cells
770) write(string,'(i2)') i
771) string = 'Rate(' // trim(adjustl(string)) // ') ' &
772) // trim(reaction%mineral%mineral_names(j))
773) call OutputWriteToHeader(fid,string,'',cell_string,icolumn)
774) enddo
775) enddo
776) endif
777)
778) ! add secondary mineral volume fractions to header
779) if (print_secondary_data(5)) then
780) do j = 1, reaction%mineral%nkinmnrl
781) do i = 1, option%nsec_cells
782) write(string,'(i2)') i
783) string = 'SI(' // trim(adjustl(string)) // ') ' &
784) // trim(reaction%mineral%mineral_names(j))
785) call OutputWriteToHeader(fid,string,'',cell_string,icolumn)
786) enddo
787) enddo
788) endif
789)
790) endif
791)
792) end subroutine WriteObservationHeaderSec
793)
794) ! ************************************************************************** !
795)
796) subroutine WriteObservationHeaderForBC(fid,realization_base,coupler_name)
797) !
798) ! Print a header for data over a region
799) !
800) ! Author: Glenn Hammond
801) ! Date: 12/18/08
802) !
803)
804) use Realization_Base_class, only : realization_base_type
805) use Option_module
806) use Reaction_Aux_module
807)
808) implicit none
809)
810) PetscInt :: fid
811) class(realization_base_type) :: realization_base
812) character(len=MAXWORDLENGTH) :: coupler_name
813)
814) PetscInt :: i
815) character(len=MAXSTRINGLENGTH) :: string
816) type(option_type), pointer :: option
817) type(reaction_type), pointer :: reaction
818)
819) option => realization_base%option
820) reaction => realization_base%reaction
821)
822) select case(option%iflowmode)
823) case(FLASH2_MODE)
824) case(MPH_MODE)
825) case(IMS_MODE)
826) case(TH_MODE)
827) case(MIS_MODE)
828) case(RICHARDS_MODE)
829) string = ',"Darcy flux ' // trim(coupler_name) // &
830) ' [m^3/' // trim(realization_base%output_option%tunit) // ']"'
831) case default
832) end select
833) write(fid,'(a)',advance="no") trim(string)
834)
835) if (associated(reaction)) then
836) do i=1, reaction%naqcomp
837) ! may need to modify for molality vs molarity, but I believe molarity is correct
838) write(fid,'(a)',advance="no") ',"' // &
839) trim(reaction%primary_species_names(i)) // ' ' // &
840) trim(coupler_name) // &
841) ' [mol/' // trim(realization_base%output_option%tunit) // ']"'
842) enddo
843) endif
844)
845) end subroutine WriteObservationHeaderForBC
846)
847) ! ************************************************************************** !
848)
849) subroutine WriteObservationDataForCell(fid,realization_base,local_id)
850) !
851) ! Print data for data at a cell
852) !
853) ! Author: Glenn Hammond
854) ! Date: 02/11/08
855) !
856)
857) use Realization_Base_class, only : realization_base_type, &
858) RealizGetVariableValueAtCell
859) use Option_module
860) use Grid_module
861) use Field_module
862) use Patch_module
863) use Reaction_Aux_module
864) use Variables_module
865)
866) implicit none
867)
868) PetscInt :: fid, i
869) class(realization_base_type) :: realization_base
870) PetscInt :: local_id
871) PetscInt :: ghosted_id
872) PetscReal :: temp_real
873) type(option_type), pointer :: option
874) type(grid_type), pointer :: grid
875) type(field_type), pointer :: field
876) type(patch_type), pointer :: patch
877) type(reaction_type), pointer :: reaction
878) type(output_option_type), pointer :: output_option
879) type(output_variable_type), pointer :: cur_variable
880)
881) option => realization_base%option
882) patch => realization_base%patch
883) grid => patch%grid
884) field => realization_base%field
885) output_option => realization_base%output_option
886)
887) 100 format(es14.6)
888) 101 format(i2)
889) 110 format(es14.6)
890) 111 format(i2)
891)
892) ghosted_id = grid%nL2G(local_id)
893)
894) ! loop over observation variables and write to file
895) cur_variable => output_option%output_obs_variable_list%first
896) do
897) if (.not.associated(cur_variable)) exit
898) if (cur_variable%plot_only) then
899) cur_variable => cur_variable%next
900) cycle
901) endif
902) temp_real = RealizGetVariableValueAtCell(realization_base, &
903) cur_variable%ivar, &
904) cur_variable%isubvar,ghosted_id)
905) if (cur_variable%iformat == 0) then ! real
906) write(fid,110,advance="no") temp_real
907) else ! integer
908) write(fid,111,advance="no") int(temp_real + 1.d-5)
909) endif
910) cur_variable => cur_variable%next
911) enddo
912)
913) end subroutine WriteObservationDataForCell
914)
915) ! ************************************************************************** !
916)
917) subroutine WriteObservationDataForCoord(fid,realization_base,region)
918) !
919) ! Print data for data at a coordinate
920) !
921) ! Author: Glenn Hammond
922) ! Date: 04/11/08
923) !
924)
925) use Realization_Base_class, only : realization_base_type
926) use Option_module
927) use Region_module
928) use Grid_module
929) use Field_module
930) use Patch_module
931) use Reaction_Aux_module
932) use Variables_module
933)
934) use Grid_Structured_module
935)
936) implicit none
937)
938) PetscInt :: fid
939) class(realization_base_type) :: realization_base
940) type(region_type) :: region
941)
942) PetscInt :: local_id
943) PetscInt :: ghosted_id
944) PetscReal :: temp_real
945) type(option_type), pointer :: option
946) type(grid_type), pointer :: grid
947) type(field_type), pointer :: field
948) type(patch_type), pointer :: patch
949) type(reaction_type), pointer :: reaction
950) type(output_option_type), pointer :: output_option
951) type(output_variable_type), pointer :: cur_variable
952)
953) PetscInt :: ghosted_ids(8)
954) PetscInt :: count
955) PetscInt :: i, j, k
956) PetscInt :: istart, iend, jstart, jend, kstart, kend
957)
958) option => realization_base%option
959) patch => realization_base%patch
960) grid => patch%grid
961) field => realization_base%field
962) output_option => realization_base%output_option
963)
964) 100 format(es14.6)
965) 101 format(i2)
966) 110 format(es14.6)
967) 111 format(i2)
968)
969) count = 0
970) local_id = region%cell_ids(1)
971) ghosted_id = grid%nL2G(local_id)
972) call StructGridGetIJKFromGhostedID(grid%structured_grid,ghosted_id,i,j,k)
973) istart = i
974) iend = i
975) jstart = j
976) jend = j
977) kstart = k
978) kend = k
979) ! find the neighboring cells, between which to interpolate
980) if (grid%x(ghosted_id) > region%coordinates(ONE_INTEGER)%x) then
981) if (i > 1) then
982) istart = i-1
983) endif
984) else
985) if (i < grid%structured_grid%ngx) then
986) iend = i+1
987) endif
988) endif
989) if (grid%y(ghosted_id) > region%coordinates(ONE_INTEGER)%y) then
990) if (j > 1) then
991) jstart = j-1
992) endif
993) else
994) if (j < grid%structured_grid%ngy) then
995) jend = j+1
996) endif
997) endif
998) if (grid%z(ghosted_id) > region%coordinates(ONE_INTEGER)%z) then
999) if (k > 1) then
1000) kstart = k-1
1001) endif
1002) else
1003) if (k < grid%structured_grid%ngz) then
1004) kend = k+1
1005) endif
1006) endif
1007) count = 0
1008) do k=kstart,kend
1009) do j=jstart,jend
1010) do i=istart,iend
1011) count = count + 1
1012) ghosted_ids(count) = i + (j-1)*grid%structured_grid%ngx + &
1013) (k-1)*grid%structured_grid%ngxy
1014) enddo
1015) enddo
1016) enddo
1017)
1018) ! loop over observation variables and write to file
1019) cur_variable => output_option%output_obs_variable_list%first
1020) do
1021) if (.not.associated(cur_variable)) exit
1022) if (cur_variable%plot_only) then
1023) cur_variable => cur_variable%next
1024) cycle
1025) endif
1026) temp_real = OutputGetVarFromArrayAtCoord(realization_base, &
1027) cur_variable%ivar, &
1028) cur_variable%isubvar, &
1029) region%coordinates(ONE_INTEGER)%x, &
1030) region%coordinates(ONE_INTEGER)%y, &
1031) region%coordinates(ONE_INTEGER)%z, &
1032) count,ghosted_ids)
1033) if (cur_variable%iformat == 0) then ! real
1034) write(fid,110,advance="no") temp_real
1035) else ! integer
1036) write(fid,111,advance="no") int(temp_real + 1.d-5)
1037) endif
1038) cur_variable => cur_variable%next
1039) enddo
1040)
1041) end subroutine WriteObservationDataForCoord
1042)
1043) ! ************************************************************************** !
1044)
1045) subroutine WriteObservationDataForBC(fid,realization_base,patch,connection_set)
1046) !
1047) ! Print flux data for a boundary condition
1048) !
1049) ! Author: Glenn Hammond
1050) ! Date: 12/18/08
1051) !
1052)
1053) use Realization_Base_class, only : realization_base_type
1054) use Option_module
1055) use Connection_module
1056) use Patch_module
1057) use Reaction_Aux_module
1058)
1059) implicit none
1060)
1061) PetscInt :: fid
1062) class(realization_base_type) :: realization_base
1063) type(patch_type), pointer :: patch
1064) type(connection_set_type), pointer :: connection_set
1065)
1066) PetscInt :: i
1067) PetscInt :: iconn
1068) PetscInt :: offset
1069) PetscInt :: iphase
1070) PetscMPIInt :: int_mpi
1071) PetscReal :: sum_volumetric_flux(realization_base%option%nphase)
1072) PetscReal :: sum_volumetric_flux_global(realization_base%option%nphase)
1073) PetscReal :: sum_solute_flux(realization_base%option%ntrandof)
1074) PetscReal :: sum_solute_flux_global(realization_base%option%ntrandof)
1075) type(option_type), pointer :: option
1076) type(reaction_type), pointer :: reaction
1077) PetscErrorCode :: ierr
1078)
1079) option => realization_base%option
1080) reaction => realization_base%reaction
1081)
1082) 100 format(es14.6)
1083) !100 format(es16.9)
1084) 101 format(i2)
1085) 110 format(es14.6)
1086) !110 format(',',es16.9)
1087) 111 format(i2)
1088)
1089) iphase = 1
1090)
1091) ! sum up fluxes across region
1092) if (associated(connection_set)) then
1093) offset = connection_set%offset
1094) select case(option%iflowmode)
1095) case(MPH_MODE,TH_MODE,IMS_MODE,FLASH2_MODE,G_MODE)
1096) case(MIS_MODE)
1097) case(RICHARDS_MODE)
1098) sum_volumetric_flux = 0.d0
1099) if (associated(connection_set)) then
1100) do iconn = 1, connection_set%num_connections
1101) sum_volumetric_flux(:) = sum_volumetric_flux(:) + &
1102) patch%boundary_velocities(iphase,offset+iconn)* &
1103) connection_set%area(iconn)
1104) enddo
1105) endif
1106) int_mpi = option%nphase
1107) call MPI_Reduce(sum_volumetric_flux,sum_volumetric_flux_global, &
1108) int_mpi,MPI_DOUBLE_PRECISION,MPI_SUM, &
1109) option%io_rank,option%mycomm,ierr)
1110) if (option%myrank == option%io_rank) then
1111) do i = 1, option%nphase
1112) write(fid,110,advance="no") sum_volumetric_flux_global(i)
1113) enddo
1114) endif
1115) end select
1116)
1117) if (associated(reaction)) then
1118) sum_solute_flux = 0.d0
1119) if (associated(connection_set)) then
1120) do iconn = 1, connection_set%num_connections
1121) sum_solute_flux(:) = sum_solute_flux(:) + &
1122) patch%boundary_tran_fluxes(:,offset+iconn)* &
1123) connection_set%area(iconn)
1124) enddo
1125) endif
1126) int_mpi = option%ntrandof
1127) call MPI_Reduce(sum_solute_flux,sum_solute_flux_global, &
1128) int_mpi,MPI_DOUBLE_PRECISION,MPI_SUM, &
1129) option%io_rank,option%mycomm,ierr)
1130) if (option%myrank == option%io_rank) then
1131) !we currently only print the aqueous components
1132) do i = 1, reaction%naqcomp
1133) write(fid,110,advance="no") sum_solute_flux_global(i)
1134) enddo
1135) endif
1136) endif
1137)
1138) endif
1139)
1140) end subroutine WriteObservationDataForBC
1141)
1142) ! ************************************************************************** !
1143)
1144) subroutine WriteVelocityAtCell(fid,realization_base,local_id)
1145) !
1146) ! Computes velocities at a grid cell
1147) ! note: limited to structured grids
1148) !
1149) ! Author: Glenn Hammond
1150) ! Date: 03/20/08
1151) !
1152)
1153) use Realization_Base_class, only : realization_base_type
1154) use Option_module
1155)
1156) implicit none
1157)
1158) PetscInt :: fid
1159) class(realization_base_type) :: realization_base
1160) type(option_type), pointer :: option
1161) PetscInt :: local_id
1162) PetscInt :: iphase
1163)
1164) PetscReal :: velocity(1:3)
1165) option => realization_base%option
1166)
1167) 200 format(3(es14.6))
1168)
1169) iphase = 1
1170) velocity = GetVelocityAtCell(fid,realization_base,local_id,iphase)
1171) write(fid,200,advance="no") velocity(1:3)* &
1172) realization_base%output_option%tconv
1173)
1174) if (option%nphase > 1) then
1175) iphase = 2
1176) velocity = GetVelocityAtCell(fid,realization_base,local_id,iphase)
1177) write(fid,200,advance="no") velocity(1:3)* &
1178) realization_base%output_option%tconv
1179) endif
1180)
1181) end subroutine WriteVelocityAtCell
1182)
1183) ! ************************************************************************** !
1184)
1185) subroutine WriteVelocityAtCell2(fid,realization_base,local_id,velocities)
1186) !
1187) ! Writes the velocity previoiusly calculated and stored in vecs at cell
1188) !
1189) ! Author: Glenn Hammond
1190) ! Date: 07/08/16
1191) !
1192) use Realization_Base_class, only : realization_base_type
1193) use Option_module
1194)
1195) implicit none
1196)
1197) PetscInt :: fid
1198) class(realization_base_type) :: realization_base
1199) type(option_type), pointer :: option
1200) PetscInt :: local_id
1201) PetscReal :: velocities(:,:,:)
1202)
1203) PetscReal :: velocity(3)
1204) option => realization_base%option
1205)
1206) 200 format(3(es14.6))
1207)
1208) velocity = velocities(:,local_id,1)
1209) write(fid,200,advance="no") velocity(1:3)* &
1210) realization_base%output_option%tconv
1211)
1212) if (option%nphase > 1) then
1213) velocity = velocities(:,local_id,2)
1214) write(fid,200,advance="no") velocity(1:3)* &
1215) realization_base%output_option%tconv
1216) endif
1217)
1218) end subroutine WriteVelocityAtCell2
1219)
1220) ! ************************************************************************** !
1221)
1222) function GetVelocityAtCell(fid,realization_base,local_id,iphase)
1223) !
1224) ! Computes velocities at a grid cell
1225) ! note: limited to structured grids
1226) !
1227) ! Author: Glenn Hammond
1228) ! Date: 03/20/08
1229) !
1230)
1231) use Realization_Base_class, only : realization_base_type
1232) use Option_module
1233) use Grid_module
1234) use Field_module
1235) use Patch_module
1236) use Connection_module
1237) use Coupler_module
1238)
1239) implicit none
1240)
1241) PetscReal :: GetVelocityAtCell(3)
1242) PetscInt :: fid
1243) class(realization_base_type) :: realization_base
1244) PetscInt :: local_id
1245)
1246) PetscInt :: ghosted_id
1247) type(option_type), pointer :: option
1248) type(grid_type), pointer :: grid
1249) type(field_type), pointer :: field
1250) type(patch_type), pointer :: patch
1251) type(connection_set_list_type), pointer :: connection_set_list
1252) type(coupler_type), pointer :: boundary_condition
1253) type(connection_set_type), pointer :: cur_connection_set
1254) PetscInt :: iconn, sum_connection
1255) PetscInt :: local_id_up, local_id_dn
1256) PetscInt :: direction, iphase
1257) PetscReal :: area
1258) PetscReal :: sum_velocity(1:3), sum_area(1:3), velocity(1:3)
1259)
1260) option => realization_base%option
1261) patch => realization_base%patch
1262) grid => patch%grid
1263) field => realization_base%field
1264)
1265) sum_velocity = 0.d0
1266) sum_area = 0.d0
1267) ! iphase = 1
1268)
1269) ! interior velocities
1270) connection_set_list => grid%internal_connection_set_list
1271) cur_connection_set => connection_set_list%first
1272) sum_connection = 0
1273) do
1274) if (.not.associated(cur_connection_set)) exit
1275) do iconn = 1, cur_connection_set%num_connections
1276) sum_connection = sum_connection + 1
1277) local_id_up = grid%nG2L(cur_connection_set%id_up(iconn)) ! = zero for ghost nodes
1278) local_id_dn = grid%nG2L(cur_connection_set%id_dn(iconn)) ! = zero for ghost nodes
1279) if (local_id_up == local_id .or. local_id_dn == local_id) then
1280) do direction=1,3
1281) area = cur_connection_set%area(iconn)* &
1282) !geh: no dabs() here
1283) cur_connection_set%dist(direction,iconn)
1284) sum_velocity(direction) = sum_velocity(direction) + &
1285) patch%internal_velocities(iphase,sum_connection)* &
1286) area
1287) sum_area(direction) = sum_area(direction) + dabs(area)
1288) enddo
1289) endif
1290) enddo
1291) cur_connection_set => cur_connection_set%next
1292) enddo
1293)
1294) ! boundary velocities
1295) boundary_condition => patch%boundary_condition_list%first
1296) sum_connection = 0
1297) do
1298) if (.not.associated(boundary_condition)) exit
1299) cur_connection_set => boundary_condition%connection_set
1300) do iconn = 1, cur_connection_set%num_connections
1301) sum_connection = sum_connection + 1
1302) if (cur_connection_set%id_dn(iconn) == local_id) then
1303) do direction=1,3
1304) area = cur_connection_set%area(iconn)* &
1305) !geh: no dabs() here
1306) cur_connection_set%dist(direction,iconn)
1307) sum_velocity(direction) = sum_velocity(direction) + &
1308) patch%boundary_velocities(iphase,sum_connection)* &
1309) area
1310) sum_area(direction) = sum_area(direction) + dabs(area)
1311) enddo
1312) endif
1313) enddo
1314) boundary_condition => boundary_condition%next
1315) enddo
1316)
1317) velocity = 0.d0
1318) do direction = 1,3
1319) if (abs(sum_area(direction)) > 1.d-40) &
1320) velocity(direction) = sum_velocity(direction)/sum_area(direction)
1321) enddo
1322)
1323) GetVelocityAtCell = velocity
1324)
1325) end function GetVelocityAtCell
1326)
1327) ! ************************************************************************** !
1328)
1329) subroutine WriteVelocityAtCoord(fid,realization_base,region)
1330) !
1331) ! Computes velocities at a coordinate
1332) ! note: limited to structured grids
1333) !
1334) ! Author: Glenn Hammond
1335) ! Date: 03/20/08
1336) !
1337)
1338) use Realization_Base_class, only : realization_base_type
1339) use Region_module
1340) use Option_module
1341)
1342) implicit none
1343)
1344) PetscInt :: fid
1345) class(realization_base_type) :: realization_base
1346) type(region_type) :: region
1347) type(option_type), pointer :: option
1348) PetscInt :: local_id
1349) PetscInt :: iphase
1350) PetscReal :: coordinate(3)
1351)
1352) PetscReal :: velocity(1:3)
1353)
1354) option => realization_base%option
1355)
1356) 200 format(3(es14.6))
1357)
1358) iphase = 1
1359) velocity = GetVelocityAtCoord(fid,realization_base,region%cell_ids(1), &
1360) region%coordinates(ONE_INTEGER)%x, &
1361) region%coordinates(ONE_INTEGER)%y, &
1362) region%coordinates(ONE_INTEGER)%z,iphase)
1363) write(fid,200,advance="no") velocity(1:3)*realization_base%output_option%tconv
1364)
1365) if (option%nphase > 1) then
1366) iphase = 2
1367) velocity = GetVelocityAtCoord(fid,realization_base,region%cell_ids(1), &
1368) region%coordinates(ONE_INTEGER)%x, &
1369) region%coordinates(ONE_INTEGER)%y, &
1370) region%coordinates(ONE_INTEGER)%z,iphase)
1371) write(fid,200,advance="no") velocity(1:3)*realization_base%output_option%tconv
1372) endif
1373)
1374) end subroutine WriteVelocityAtCoord
1375)
1376) ! ************************************************************************** !
1377)
1378) function GetVelocityAtCoord(fid,realization_base,local_id,x,y,z,iphase)
1379) !
1380) ! Computes velocities at a coordinate
1381) ! note: limited to structured grids
1382) !
1383) ! Author: Glenn Hammond
1384) ! Date: 03/20/08
1385) !
1386) use Realization_Base_class, only : realization_base_type
1387) use Option_module
1388) use Grid_module
1389) use Field_module
1390) use Patch_module
1391) use Connection_module
1392) use Coupler_module
1393)
1394) implicit none
1395)
1396) PetscReal :: GetVelocityAtCoord(3)
1397) PetscInt :: fid
1398) class(realization_base_type) :: realization_base
1399) PetscInt :: local_id
1400) PetscReal :: x, y, z
1401)
1402) PetscInt :: ghosted_id
1403) type(option_type), pointer :: option
1404) type(grid_type), pointer :: grid
1405) type(field_type), pointer :: field
1406) type(patch_type), pointer :: patch
1407) type(connection_set_list_type), pointer :: connection_set_list
1408) type(coupler_type), pointer :: boundary_condition
1409) type(connection_set_type), pointer :: cur_connection_set
1410) PetscInt :: iconn, sum_connection
1411) PetscInt :: local_id_up, local_id_dn
1412) PetscReal :: cell_coord(3), face_coord
1413) PetscReal :: coordinate(3)
1414) PetscInt :: direction, iphase
1415) PetscReal :: area, weight, distance
1416) PetscReal :: sum_velocity(1:3), velocity(1:3)
1417) PetscReal :: sum_weight(1:3)
1418)
1419) option => realization_base%option
1420) patch => realization_base%patch
1421) grid => patch%grid
1422) field => realization_base%field
1423)
1424) sum_velocity = 0.d0
1425) sum_weight = 0.d0
1426) ! iphase = 1
1427)
1428) ghosted_id = grid%nL2G(local_id)
1429)
1430) coordinate(X_DIRECTION) = x
1431) coordinate(Y_DIRECTION) = y
1432) coordinate(Z_DIRECTION) = z
1433)
1434) cell_coord(X_DIRECTION) = grid%x(ghosted_id)
1435) cell_coord(Y_DIRECTION) = grid%y(ghosted_id)
1436) cell_coord(Z_DIRECTION) = grid%z(ghosted_id)
1437)
1438) ! interior velocities
1439) connection_set_list => grid%internal_connection_set_list
1440) cur_connection_set => connection_set_list%first
1441) sum_connection = 0
1442) do
1443) if (.not.associated(cur_connection_set)) exit
1444) do iconn = 1, cur_connection_set%num_connections
1445) sum_connection = sum_connection + 1
1446) local_id_up = grid%nG2L(cur_connection_set%id_up(iconn)) ! = zero for ghost nodes
1447) local_id_dn = grid%nG2L(cur_connection_set%id_dn(iconn)) ! = zero for ghost nodes
1448) if (local_id_up == local_id .or. local_id_dn == local_id) then
1449) do direction=1,3
1450) if (local_id_up == local_id) then
1451) face_coord = cell_coord(direction) + &
1452) cur_connection_set%dist(-1,iconn)* &
1453) cur_connection_set%dist(0,iconn)* &
1454) cur_connection_set%dist(direction,iconn)
1455) else
1456) face_coord = cell_coord(direction) - &
1457) (1.d0-cur_connection_set%dist(-1,iconn))* &
1458) cur_connection_set%dist(0,iconn)* &
1459) cur_connection_set%dist(direction,iconn)
1460) endif
1461) distance = dabs(face_coord-coordinate(direction))
1462) if (distance < 1.d-40) distance = 1.d-40
1463) weight = cur_connection_set%area(iconn)* &
1464) dabs(cur_connection_set%dist(direction,iconn))/ &
1465) distance
1466)
1467) sum_velocity(direction) = sum_velocity(direction) + &
1468) cur_connection_set%dist(direction,iconn)* &
1469) patch%internal_velocities(iphase,sum_connection)* &
1470) weight
1471) sum_weight(direction) = sum_weight(direction) + weight
1472) enddo
1473) endif
1474) enddo
1475) cur_connection_set => cur_connection_set%next
1476) enddo
1477)
1478) ! boundary velocities
1479) boundary_condition => patch%boundary_condition_list%first
1480) sum_connection = 0
1481) do
1482) if (.not.associated(boundary_condition)) exit
1483) cur_connection_set => boundary_condition%connection_set
1484) do iconn = 1, cur_connection_set%num_connections
1485) sum_connection = sum_connection + 1
1486) if (cur_connection_set%id_dn(iconn) == local_id) then
1487) do direction=1,3
1488) face_coord = cell_coord(direction) - &
1489) ! (1.d0-cur_connection_set%dist(-1,iconn))* & ! fraction upwind is always 0.d0
1490) cur_connection_set%dist(0,iconn)* &
1491) cur_connection_set%dist(direction,iconn)
1492) distance = dabs(face_coord-coordinate(direction))
1493) if (distance < 1.d-40) distance = 1.d-40
1494) weight = cur_connection_set%area(iconn)* &
1495) dabs(cur_connection_set%dist(direction,iconn))/ &
1496) distance
1497) sum_velocity(direction) = sum_velocity(direction) + &
1498) cur_connection_set%dist(direction,iconn)* &
1499) patch%boundary_velocities(iphase,sum_connection)* &
1500) weight
1501) sum_weight(direction) = sum_weight(direction) + weight
1502) enddo
1503) endif
1504) enddo
1505) boundary_condition => boundary_condition%next
1506) enddo
1507)
1508) velocity = 0.d0
1509) do direction = 1,3
1510) if (abs(sum_weight(direction)) > 1.d-40) &
1511) velocity(direction) = sum_velocity(direction)/sum_weight(direction)
1512) enddo
1513)
1514) GetVelocityAtCoord = velocity
1515)
1516) end function GetVelocityAtCoord
1517)
1518) ! ************************************************************************** !
1519)
1520) subroutine WriteObservationSecondaryDataAtCell(fid,realization_base,local_id,ivar)
1521) !
1522) ! Print data for data at a cell
1523) !
1524) ! Author: Satish Karra, LANL
1525) ! Date: 10/4/12
1526) !
1527)
1528) use Realization_Base_class, only : realization_base_type, &
1529) RealizGetVariableValueAtCell
1530) use Option_module
1531) use Grid_module
1532) use Field_module
1533) use Patch_module
1534) use Variables_module
1535) use Reaction_Aux_module
1536)
1537) implicit none
1538)
1539) PetscInt :: fid,i,naqcomp,nkinmnrl
1540) class(realization_base_type) :: realization_base
1541) PetscInt :: local_id
1542) PetscInt :: ghosted_id
1543) type(option_type), pointer :: option
1544) type(grid_type), pointer :: grid
1545) type(field_type), pointer :: field
1546) type(patch_type), pointer :: patch
1547) type(output_option_type), pointer :: output_option
1548) type(reaction_type), pointer :: reaction
1549) PetscInt :: ivar
1550)
1551) option => realization_base%option
1552) patch => realization_base%patch
1553) grid => patch%grid
1554) field => realization_base%field
1555) output_option => realization_base%output_option
1556)
1557) 100 format(es14.6)
1558) 101 format(i2)
1559) 110 format(es14.6)
1560) 111 format(i2)
1561)
1562) ghosted_id = grid%nL2G(local_id)
1563)
1564) if (option%nsec_cells > 0) then
1565) if (ivar == PRINT_SEC_TEMP) then
1566) select case(option%iflowmode)
1567) case(MPH_MODE,TH_MODE)
1568) do i = 1, option%nsec_cells
1569) write(fid,110,advance="no") &
1570) RealizGetVariableValueAtCell(realization_base,SECONDARY_TEMPERATURE,i, &
1571) ghosted_id)
1572) enddo
1573) end select
1574) endif
1575) if (option%ntrandof > 0) then
1576) reaction => realization_base%reaction
1577) if (ivar == PRINT_SEC_CONC) then
1578) do naqcomp = 1, reaction%naqcomp
1579) do i = 1, option%nsec_cells
1580) write(fid,110,advance="no") &
1581) RealizGetVariableValueAtCell(realization_base,SECONDARY_CONCENTRATION,i, &
1582) ghosted_id,naqcomp)
1583) enddo
1584) enddo
1585) endif
1586) if (ivar == PRINT_SEC_MIN_VOLFRAC) then
1587) do nkinmnrl = 1, reaction%mineral%nkinmnrl
1588) do i = 1, option%nsec_cells
1589) write(fid,110,advance="no") &
1590) RealizGetVariableValueAtCell(realization_base,SEC_MIN_VOLFRAC,i, &
1591) ghosted_id,nkinmnrl)
1592) enddo
1593) enddo
1594) endif
1595) if (ivar == PRINT_SEC_MIN_RATE) then
1596) do nkinmnrl = 1, reaction%mineral%nkinmnrl
1597) do i = 1, option%nsec_cells
1598) write(fid,110,advance="no") &
1599) RealizGetVariableValueAtCell(realization_base,SEC_MIN_RATE,i, &
1600) ghosted_id,nkinmnrl)
1601) enddo
1602) enddo
1603) endif
1604) if (ivar == PRINT_SEC_MIN_SI) then
1605) do nkinmnrl = 1, reaction%mineral%nkinmnrl
1606) do i = 1, option%nsec_cells
1607) write(fid,110,advance="no") &
1608) RealizGetVariableValueAtCell(realization_base,SEC_MIN_SI,i, &
1609) ghosted_id,nkinmnrl)
1610) enddo
1611) enddo
1612) endif
1613) endif
1614) endif
1615)
1616) end subroutine WriteObservationSecondaryDataAtCell
1617)
1618) ! ************************************************************************** !
1619)
1620) subroutine OutputIntegralFlux(realization_base)
1621) !
1622) ! Print integral fluxes to Tecplot POINT format
1623) !
1624) ! Author: Glenn Hammond
1625) ! Date: 10/21/14
1626) !
1627)
1628) use Realization_Subsurface_class, only : realization_subsurface_type
1629) use Realization_Base_class, only : realization_base_type
1630) use Option_module
1631) use Grid_module
1632) use Patch_module
1633) use Output_Aux_module
1634) use Reaction_Aux_module
1635) use Integral_Flux_module
1636) use Utility_module
1637)
1638) implicit none
1639)
1640) class(realization_base_type), target :: realization_base
1641)
1642) type(option_type), pointer :: option
1643) type(patch_type), pointer :: patch
1644) type(grid_type), pointer :: grid
1645) type(output_option_type), pointer :: output_option
1646) type(reaction_type), pointer :: reaction
1647)
1648) character(len=MAXSTRINGLENGTH) :: filename
1649) character(len=MAXWORDLENGTH) :: word, units
1650) character(len=MAXSTRINGLENGTH) :: string
1651) type(integral_flux_type), pointer :: integral_flux
1652) PetscReal :: flow_dof_scale(10)
1653) PetscReal, allocatable :: array(:,:)
1654) PetscReal, allocatable :: array_global(:,:)
1655) PetscReal, allocatable :: instantaneous_array(:)
1656) PetscInt, parameter :: fid = 86
1657) PetscInt :: i, j
1658) PetscInt :: istart, iend
1659) PetscInt :: icol
1660) PetscMPIInt :: int_mpi
1661) PetscReal :: tempreal
1662) PetscErrorCode :: ierr
1663)
1664) patch => realization_base%patch
1665) grid => patch%grid
1666) option => realization_base%option
1667) output_option => realization_base%output_option
1668) reaction => realization_base%reaction
1669)
1670) if (.not.associated(patch%integral_flux_list%first)) return
1671)
1672) flow_dof_scale = 1.d0
1673) select case(option%iflowmode)
1674) case(RICHARDS_MODE)
1675) flow_dof_scale(1) = FMWH2O
1676) case(TH_MODE)
1677) flow_dof_scale(1) = FMWH2O
1678) case(MIS_MODE)
1679) flow_dof_scale(1) = FMWH2O
1680) flow_dof_scale(2) = FMWGLYC
1681) case(G_MODE)
1682) flow_dof_scale(1) = FMWH2O
1683) flow_dof_scale(2) = FMWAIR
1684) case(MPH_MODE,FLASH2_MODE,IMS_MODE)
1685) flow_dof_scale(1) = FMWH2O
1686) flow_dof_scale(2) = FMWCO2
1687) end select
1688)
1689) if (len_trim(output_option%plot_name) > 2) then
1690) filename = trim(output_option%plot_name) // '-int.dat'
1691) else
1692) filename = trim(option%global_prefix) // trim(option%group_prefix) // &
1693) '-int.dat'
1694) endif
1695)
1696) ! open file
1697) if (option%myrank == option%io_rank) then
1698)
1699) !geh option%io_buffer = '--> write tecplot mass balance file: ' // trim(filename)
1700) !geh call printMsg(option)
1701)
1702) if (output_option%print_column_ids) then
1703) icol = 1
1704) else
1705) icol = -1
1706) endif
1707)
1708) if (integral_flux_first .or. .not.FileExists(filename)) then
1709) open(unit=fid,file=filename,action="write",status="replace")
1710)
1711) ! write header
1712) write(fid,'(a)',advance="no") ' "Time [' // trim(output_option%tunit) // &
1713) ']"'
1714)
1715) if (option%iflowmode > 0) then
1716) call OutputWriteToHeader(fid,'dt_flow',output_option%tunit,'',icol)
1717) endif
1718)
1719) if (option%ntrandof > 0) then
1720) call OutputWriteToHeader(fid,'dt_tran',output_option%tunit,'',icol)
1721) endif
1722)
1723) integral_flux => patch%integral_flux_list%first
1724) do
1725) if (.not.associated(integral_flux)) exit
1726) select case(option%iflowmode)
1727) case(RICHARDS_MODE,TH_MODE,MIS_MODE,G_MODE,MPH_MODE,FLASH2_MODE, &
1728) IMS_MODE)
1729) string = trim(integral_flux%name) // ' Water'
1730) call OutputWriteToHeader(fid,string,'kg','',icol)
1731) units = 'kg/' // trim(output_option%tunit) // ''
1732) string = trim(integral_flux%name) // ' Water'
1733) call OutputWriteToHeader(fid,string,units,'',icol)
1734) end select
1735) select case(option%iflowmode)
1736) case(MIS_MODE)
1737) string = trim(integral_flux%name) // ' Glycol'
1738) call OutputWriteToHeader(fid,string,'kg','',icol)
1739) units = 'kg/' // trim(output_option%tunit) // ''
1740) string = trim(integral_flux%name) // ' Glycol'
1741) call OutputWriteToHeader(fid,string,units,'',icol)
1742) case(G_MODE)
1743) string = trim(integral_flux%name) // ' Air'
1744) call OutputWriteToHeader(fid,string,'kg','',icol)
1745) units = 'kg/' // trim(output_option%tunit) // ''
1746) string = trim(integral_flux%name) // ' Air'
1747) call OutputWriteToHeader(fid,string,units,'',icol)
1748) case(MPH_MODE,FLASH2_MODE,IMS_MODE)
1749) string = trim(integral_flux%name) // ' CO2'
1750) call OutputWriteToHeader(fid,string,'kg','',icol)
1751) units = 'kg/' // trim(output_option%tunit) // ''
1752) string = trim(integral_flux%name) // ' CO2'
1753) call OutputWriteToHeader(fid,string,units,'',icol)
1754) end select
1755) select case(option%iflowmode)
1756) case(TH_MODE,MIS_MODE,G_MODE,MPH_MODE,FLASH2_MODE,IMS_MODE)
1757) string = trim(integral_flux%name) // ' Energy'
1758) call OutputWriteToHeader(fid,string,'MJ','',icol)
1759) units = 'MJ/' // trim(output_option%tunit) // ''
1760) string = trim(integral_flux%name) // ' Energy'
1761) call OutputWriteToHeader(fid,string,units,'',icol)
1762) end select
1763)
1764) if (option%ntrandof > 0) then
1765) units = 'mol/' // trim(output_option%tunit) // ''
1766) do i=1,reaction%naqcomp
1767) if (reaction%primary_species_print(i)) then
1768) string = trim(integral_flux%name) // ' ' // &
1769) trim(reaction%primary_species_names(i))
1770) call OutputWriteToHeader(fid,string,'mol','',icol)
1771) string = trim(integral_flux%name) // ' ' // &
1772) trim(reaction%primary_species_names(i))
1773) call OutputWriteToHeader(fid,string,units,'',icol)
1774) endif
1775) enddo
1776) endif
1777) integral_flux => integral_flux%next
1778) enddo
1779) write(fid,'(a)') ''
1780) else
1781) open(unit=fid,file=filename,action="write",status="old",position="append")
1782) endif
1783) endif
1784)
1785) 100 format(100es17.8)
1786) 110 format(100es17.8)
1787) 120 format(100es17.8e3)
1788)
1789) ! write time
1790) if (option%myrank == option%io_rank) then
1791) write(fid,100,advance="no") option%time/output_option%tconv
1792) endif
1793)
1794) if (option%nflowdof > 0) then
1795) if (option%myrank == option%io_rank) &
1796) write(fid,100,advance="no") option%flow_dt/output_option%tconv
1797) endif
1798) if (option%ntrandof > 0) then
1799) if (option%myrank == option%io_rank) &
1800) write(fid,100,advance="no") option%tran_dt/output_option%tconv
1801) endif
1802)
1803) allocate(array(option%nflowdof + option%ntrandof,2))
1804) allocate(array_global(option%nflowdof + option%ntrandof,2))
1805) allocate(instantaneous_array(max(option%nflowdof,option%ntrandof)))
1806) integral_flux => patch%integral_flux_list%first
1807) do
1808) if (.not.associated(integral_flux)) exit
1809) array = 0.d0
1810) array_global = 0.d0
1811) if (option%nflowdof > 0) then
1812) istart = 1
1813) iend = option%nflowdof
1814) instantaneous_array = 0.d0
1815) call IntegralFluxGetInstantaneous(integral_flux, &
1816) patch%internal_flow_fluxes, &
1817) patch%boundary_flow_fluxes, &
1818) option%nflowdof, &
1819) instantaneous_array,option)
1820) array(istart:iend,1) = &
1821) integral_flux%integral_value(istart:iend)
1822) array(istart:iend,2) = &
1823) instantaneous_array(1:option%nflowdof)
1824) endif
1825) if (option%ntrandof > 0) then
1826) istart = option%nflowdof+1
1827) iend = option%nflowdof+option%ntrandof
1828) instantaneous_array = 0.d0
1829) call IntegralFluxGetInstantaneous(integral_flux, &
1830) patch%internal_tran_fluxes, &
1831) patch%boundary_tran_fluxes, &
1832) option%ntrandof, &
1833) instantaneous_array,option)
1834) array(istart:iend,1) = &
1835) integral_flux%integral_value(istart:iend)
1836) array(istart:iend,2) = &
1837) instantaneous_array(1:option%ntrandof)
1838) endif
1839) int_mpi = size(array)
1840) call MPI_Reduce(array,array_global, &
1841) int_mpi,MPI_DOUBLE_PRECISION,MPI_SUM, &
1842) option%io_rank,option%mycomm,ierr)
1843) ! time units conversion
1844) array_global(:,2) = array_global(:,2) * output_option%tconv
1845) if (option%myrank == option%io_rank) then
1846) if (option%nflowdof > 0) then
1847) do i = 1, option%nflowdof
1848) do j = 1, 2 ! 1 = integral, 2 = instantaneous
1849) tempreal = array_global(i,j)*flow_dof_scale(i)
1850) if (dabs(tempreal) > 0.d0 .and. dabs(tempreal) < 1.d-99) then
1851) write(fid,120,advance="no") tempreal
1852) else
1853) write(fid,110,advance="no") tempreal
1854) endif
1855) enddo
1856) enddo
1857) endif
1858) if (option%ntrandof > 0) then
1859) istart = option%nflowdof
1860) do i=1,reaction%naqcomp
1861) do j = 1, 2 ! 1 = integral, 2 = instantaneous
1862) if (reaction%primary_species_print(i)) then
1863) tempreal = array_global(istart+i,j)
1864) if (dabs(tempreal) > 0.d0 .and. dabs(tempreal) < 1.d-99) then
1865) write(fid,120,advance="no") tempreal
1866) else
1867) write(fid,110,advance="no") tempreal
1868) endif
1869) endif
1870) enddo
1871) enddo
1872) endif
1873) endif
1874) integral_flux => integral_flux%next
1875) enddo
1876) deallocate(array)
1877) deallocate(array_global)
1878) deallocate(instantaneous_array)
1879)
1880) if (option%myrank == option%io_rank) then
1881) write(fid,'(a)') ''
1882) close(fid)
1883) endif
1884)
1885) integral_flux_first = PETSC_FALSE
1886)
1887) end subroutine OutputIntegralFlux
1888)
1889) ! ************************************************************************** !
1890)
1891) subroutine OutputMassBalance(realization_base)
1892) !
1893) ! Print to Tecplot POINT format
1894) !
1895) ! Author: Glenn Hammond
1896) ! Date: 06/18/08
1897) !
1898)
1899) use Realization_Subsurface_class, only : realization_subsurface_type
1900) use Realization_Base_class, only : realization_base_type
1901) use Patch_module
1902) use Grid_module
1903) use Option_module
1904) use Coupler_module
1905) use Utility_module
1906) use Output_Aux_module
1907)
1908) use Richards_module, only : RichardsComputeMassBalance
1909) use Mphase_module, only : MphaseComputeMassBalance
1910) use Flash2_module, only : Flash2ComputeMassBalance
1911) use Immis_module, only : ImmisComputeMassBalance
1912) use Miscible_module, only : MiscibleComputeMassBalance
1913) use TH_module, only : THComputeMassBalance
1914) use Reactive_Transport_module, only : RTComputeMassBalance
1915) use General_module, only : GeneralComputeMassBalance
1916) use TOilIms_module, only : TOilImsComputeMassBalance
1917) use TOilIms_Aux_module ! for formula weights
1918)
1919) use Global_Aux_module
1920) use Reactive_Transport_Aux_module
1921) use Reaction_Aux_module
1922) use Material_Aux_class
1923)
1924) implicit none
1925)
1926) class(realization_base_type), target :: realization_base
1927)
1928) type(option_type), pointer :: option
1929) type(patch_type), pointer :: patch
1930) type(grid_type), pointer :: grid
1931) type(output_option_type), pointer :: output_option
1932) type(coupler_type), pointer :: coupler
1933) type(mass_balance_region_type), pointer :: cur_mbr
1934) type(global_auxvar_type), pointer :: global_auxvars_bc_or_ss(:)
1935) type(reactive_transport_auxvar_type), pointer :: rt_auxvars_bc_or_ss(:)
1936) type(reactive_transport_auxvar_type), pointer :: rt_auxvars(:)
1937)
1938) class(material_auxvar_type), pointer :: material_auxvars(:)
1939) type(reaction_type), pointer :: reaction
1940)
1941) character(len=MAXSTRINGLENGTH) :: filename
1942) character(len=MAXWORDLENGTH) :: word, units
1943) character(len=MAXSTRINGLENGTH) :: string
1944) PetscInt :: fid = 86
1945) PetscInt :: ios
1946) PetscInt :: i,icol
1947) PetscInt :: k, j
1948) PetscInt :: local_id
1949) PetscInt :: ghosted_id
1950) PetscInt :: iconn
1951) PetscInt :: offset
1952) PetscInt :: iphase, ispec
1953) PetscInt :: icomp, imnrl
1954) PetscReal :: sum_area(4)
1955) PetscReal :: sum_area_global(4)
1956) PetscReal :: sum_kg(realization_base%option%nflowspec, &
1957) realization_base%option%nphase)
1958) PetscReal :: sum_kg_global(realization_base%option%nflowspec, &
1959) realization_base%option%nphase)
1960) PetscReal :: sum_mol(realization_base%option%ntrandof, &
1961) realization_base%option%nphase)
1962) PetscReal :: sum_mol_global(realization_base%option%ntrandof, &
1963) realization_base%option%nphase)
1964)
1965) PetscReal, allocatable :: sum_mol_mnrl(:)
1966) PetscReal, allocatable :: sum_mol_mnrl_global(:)
1967)
1968) PetscReal :: global_total_mass
1969)
1970) PetscReal :: sum_trapped(realization_base%option%nphase)
1971) PetscReal :: sum_trapped_global(realization_base%option%nphase)
1972) PetscReal :: sum_mol_ye(3), sum_mol_global_ye(3)
1973)
1974) PetscMPIInt :: int_mpi
1975) PetscBool :: bcs_done
1976) PetscErrorCode :: ierr
1977)
1978) patch => realization_base%patch
1979) grid => patch%grid
1980) option => realization_base%option
1981) reaction => realization_base%reaction
1982) output_option => realization_base%output_option
1983)
1984) if (option%ntrandof > 0) then
1985) rt_auxvars => patch%aux%RT%auxvars
1986) material_auxvars => patch%aux%Material%auxvars
1987) endif
1988)
1989) if (len_trim(output_option%plot_name) > 2) then
1990) filename = trim(output_option%plot_name) // '-mas.dat'
1991) else
1992) filename = trim(option%global_prefix) // trim(option%group_prefix) // &
1993) '-mas.dat'
1994) endif
1995)
1996) ! open file
1997) if (option%myrank == option%io_rank) then
1998)
1999) !geh option%io_buffer = '--> write tecplot mass balance file: ' // trim(filename)
2000) !geh call printMsg(option)
2001)
2002) if (output_option%print_column_ids) then
2003) icol = 1
2004) else
2005) icol = -1
2006) endif
2007)
2008) if (mass_balance_first .or. .not.FileExists(filename)) then
2009) open(unit=fid,file=filename,action="write",status="replace")
2010)
2011) ! write header
2012) write(fid,'(a)',advance="no") ' "Time [' // trim(output_option%tunit) // &
2013) ']"'
2014)
2015) if (option%iflowmode > 0) then
2016) call OutputWriteToHeader(fid,'dt_flow',output_option%tunit,'',icol)
2017) endif
2018)
2019) if (option%ntrandof > 0) then
2020) call OutputWriteToHeader(fid,'dt_tran',output_option%tunit,'',icol)
2021) endif
2022)
2023) select case(option%iflowmode)
2024) case(RICHARDS_MODE)
2025) call OutputWriteToHeader(fid,'Global Water Mass','kg','',icol)
2026)
2027) case(TH_MODE)
2028) call OutputWriteToHeader(fid,'Global Water Mass in Liquid Phase', &
2029) 'kg','',icol)
2030) case(G_MODE)
2031) call OutputWriteToHeader(fid,'Global Water Mass in Liquid Phase', &
2032) 'kg','',icol)
2033) call OutputWriteToHeader(fid,'Global Air Mass in Liquid Phase', &
2034) 'kg','',icol)
2035) call OutputWriteToHeader(fid,'Global Water Mass in Gas Phase', &
2036) 'kg','',icol)
2037) call OutputWriteToHeader(fid,'Global Air Mass in Gas Phase', &
2038) 'kg','',icol)
2039) case(TOIL_IMS_MODE)
2040) call OutputWriteToHeader(fid,'Global Water Mass', &
2041) 'kg','',icol)
2042) call OutputWriteToHeader(fid,'Global Oil Mass', &
2043) 'kg','',icol)
2044) case(MPH_MODE,FLASH2_MODE)
2045) call OutputWriteToHeader(fid,'Global Water Mass in Water Phase', &
2046) 'kmol','',icol)
2047) call OutputWriteToHeader(fid,'Global CO2 Mass in Water Phase', &
2048) 'kmol','',icol)
2049) call OutputWriteToHeader(fid,'Trapped CO2 Mass in Water Phase', &
2050) 'kmol','',icol)
2051) call OutputWriteToHeader(fid,'Global Water Mass in Gas Phase', &
2052) 'kmol','',icol)
2053) call OutputWriteToHeader(fid,'Global CO2 Mass in Gas Phase', &
2054) 'kmol','',icol)
2055) call OutputWriteToHeader(fid,'Trapped CO2 Mass in Gas Phase', &
2056) 'kmol','',icol)
2057) case(IMS_MODE)
2058) call OutputWriteToHeader(fid,'Global Water Mass in Water Phase', &
2059) 'kmol','',icol)
2060) call OutputWriteToHeader(fid,'Global CO2 Mass in Gas Phase', &
2061) 'kmol','',icol)
2062) case(MIS_MODE)
2063) call OutputWriteToHeader(fid,'Global Water Mass in Liquid Phase', &
2064) 'kg','',icol)
2065) call OutputWriteToHeader(fid,'Global Glycol Mass in Liquid Phase', &
2066) 'kg','',icol)
2067) end select
2068)
2069) if (option%ntrandof > 0) then
2070) do i=1,reaction%naqcomp
2071) if (reaction%primary_species_print(i)) then
2072) string = 'Global ' // trim(reaction%primary_species_names(i))
2073) call OutputWriteToHeader(fid,string,'mol','',icol)
2074) endif
2075) enddo
2076)
2077) do i=1,reaction%nimcomp
2078) if (reaction%immobile%print_me(i)) then
2079) string = 'Global ' // trim(reaction%immobile%names(i))
2080) call OutputWriteToHeader(fid,string,'mol','',icol)
2081) endif
2082) enddo
2083)
2084) if (option%mass_bal_detailed) then
2085) do i=1,reaction%mineral%nkinmnrl
2086) if (reaction%mineral%kinmnrl_print(i)) then
2087) string = 'Global ' // trim(reaction%mineral%kinmnrl_names(i))
2088) call OutputWriteToHeader(fid,string,'mol','',icol)
2089) endif
2090) enddo
2091) endif
2092) endif
2093)
2094) coupler => patch%boundary_condition_list%first
2095) bcs_done = PETSC_FALSE
2096) do
2097) if (.not.associated(coupler)) then
2098) if (bcs_done) then
2099) exit
2100) else
2101) bcs_done = PETSC_TRUE
2102) if (associated(patch%source_sink_list)) then
2103) coupler => patch%source_sink_list%first
2104) if (.not.associated(coupler)) exit
2105) else
2106) exit
2107) endif
2108) endif
2109) endif
2110)
2111) select case(option%iflowmode)
2112) case(RICHARDS_MODE)
2113) string = trim(coupler%name) // ' Water Mass'
2114) call OutputWriteToHeader(fid,string,'kg','',icol)
2115)
2116) units = 'kg/' // trim(output_option%tunit) // ''
2117) string = trim(coupler%name) // ' Water Mass'
2118) call OutputWriteToHeader(fid,string,units,'',icol)
2119) case(TH_MODE)
2120) string = trim(coupler%name) // ' Water Mass'
2121) call OutputWriteToHeader(fid,string,'kg','',icol)
2122)
2123) units = 'kg/' // trim(output_option%tunit) // ''
2124) string = trim(coupler%name) // ' Water Mass'
2125) call OutputWriteToHeader(fid,string,units,'',icol)
2126) case(MIS_MODE)
2127) string = trim(coupler%name) // ' Water Mass'
2128) call OutputWriteToHeader(fid,string,'kg','',icol)
2129) string = trim(coupler%name) // ' Glycol Mass'
2130) call OutputWriteToHeader(fid,string,'kg','',icol)
2131)
2132) units = 'kg/' // trim(output_option%tunit) // ''
2133) string = trim(coupler%name) // ' Water Mass'
2134) call OutputWriteToHeader(fid,string,units,'',icol)
2135) string = trim(coupler%name) // ' Glycol Mass'
2136) call OutputWriteToHeader(fid,string,units,'',icol)
2137) case(G_MODE)
2138) string = trim(coupler%name) // ' Water Mass'
2139) call OutputWriteToHeader(fid,string,'kg','',icol)
2140) string = trim(coupler%name) // ' Air Mass'
2141) call OutputWriteToHeader(fid,string,'kg','',icol)
2142)
2143) units = 'kg/' // trim(output_option%tunit) // ''
2144) string = trim(coupler%name) // ' Water Mass'
2145) call OutputWriteToHeader(fid,string,units,'',icol)
2146) string = trim(coupler%name) // ' Air Mass'
2147) call OutputWriteToHeader(fid,string,units,'',icol)
2148) case(TOIL_IMS_MODE)
2149) string = trim(coupler%name) // ' Water Mass'
2150) call OutputWriteToHeader(fid,string,'kg','',icol)
2151) string = trim(coupler%name) // ' Oil Mass'
2152) call OutputWriteToHeader(fid,string,'kg','',icol)
2153)
2154) units = 'kg/' // trim(output_option%tunit) // ''
2155) string = trim(coupler%name) // ' Water Mass'
2156) call OutputWriteToHeader(fid,string,units,'',icol)
2157) string = trim(coupler%name) // ' Oil Mass'
2158) call OutputWriteToHeader(fid,string,units,'',icol)
2159) case(MPH_MODE,FLASH2_MODE,IMS_MODE)
2160) string = trim(coupler%name) // ' Water Mass'
2161) call OutputWriteToHeader(fid,string,'kmol','',icol)
2162) string = trim(coupler%name) // ' CO2 Mass'
2163) call OutputWriteToHeader(fid,string,'kmol','',icol)
2164)
2165) units = 'kmol/' // trim(output_option%tunit) // ''
2166) string = trim(coupler%name) // ' Water Mass'
2167) call OutputWriteToHeader(fid,string,units,'',icol)
2168) string = trim(coupler%name) // ' CO2 Mass'
2169) call OutputWriteToHeader(fid,string,units,'',icol)
2170) end select
2171)
2172) if (option%ntrandof > 0) then
2173) do i=1,reaction%naqcomp
2174) if (reaction%primary_species_print(i)) then
2175) ! option%io_buffer = 'Check OutputObservation to ensure that ' // &
2176) ! 'reactive transport species units are really kmol.'
2177) ! call printErrMsg(option)
2178) string = trim(coupler%name) // ' ' // &
2179) trim(reaction%primary_species_names(i))
2180) call OutputWriteToHeader(fid,string,'mol','',icol)
2181) endif
2182) enddo
2183)
2184) units = 'mol/' // trim(output_option%tunit) // ''
2185) do i=1,reaction%naqcomp
2186) if (reaction%primary_species_print(i)) then
2187) string = trim(coupler%name) // ' ' // &
2188) trim(reaction%primary_species_names(i))
2189) call OutputWriteToHeader(fid,string,units,'',icol)
2190) endif
2191) enddo
2192) endif
2193) coupler => coupler%next
2194)
2195) enddo
2196)
2197) ! Print the mass [mol] in the specified regions (header)
2198) if (associated(output_option%mass_balance_region_list)) then
2199) cur_mbr => output_option%mass_balance_region_list
2200) do
2201) if (.not.associated(cur_mbr)) exit
2202) string = 'Region ' // trim(cur_mbr%region_name) // ' Total Mass'
2203) call OutputWriteToHeader(fid,string,'mol','',icol)
2204) cur_mbr => cur_mbr%next
2205) enddo
2206) endif
2207)
2208) #ifdef YE_FLUX
2209) !geh do offset = 1, 4
2210) !geh write(word,'(i6)') offset*100
2211) select case(option%iflowmode)
2212) case(FLASH2_MODE,MPH_MODE)
2213) write(fid,'(a)',advance="no") ',"' // &
2214) 'Plane Water Flux [mol/s]","Plane CO2 Flux [mol/s]",' // &
2215) '"Plane Energy Flux [MJ/s]"'
2216) case(RICHARDS_MODE)
2217) write(fid,'(a)',advance="no") ',"' // &
2218) 'Plane Water Flux [mol/s]"'
2219) case(TH_MODE)
2220) write(fid,'(a)',advance="no") ',"' // &
2221) trim(adjustl(word)) // 'm Water Mass [kg]"'
2222) end select
2223)
2224) if (option%ntrandof > 0) then
2225) do i=1,reaction%naqcomp
2226) if (reaction%primary_species_print(i)) then
2227) write(fid,'(a)',advance="no") ',"' // &
2228) trim(adjustl(word)) // 'm ' // &
2229) trim(reaction%primary_species_names(i)) // ' [mol]"'
2230) endif
2231) enddo
2232) endif
2233) !geh enddo
2234) #endif
2235) write(fid,'(a)') ''
2236) else
2237) open(unit=fid,file=filename,action="write",status="old",position="append")
2238) endif
2239)
2240) endif
2241)
2242) 100 format(100es16.8)
2243) 110 format(100es16.8)
2244)
2245) ! write time
2246) if (option%myrank == option%io_rank) then
2247) write(fid,100,advance="no") option%time/output_option%tconv
2248) endif
2249)
2250) if (option%nflowdof > 0) then
2251) if (option%myrank == option%io_rank) &
2252) write(fid,100,advance="no") option%flow_dt/output_option%tconv
2253) endif
2254) if (option%ntrandof > 0) then
2255) if (option%myrank == option%io_rank) &
2256) write(fid,100,advance="no") option%tran_dt/output_option%tconv
2257) endif
2258)
2259) ! print out global mass balance
2260)
2261) if (option%nflowdof > 0) then
2262) sum_kg = 0.d0
2263) sum_trapped = 0.d0
2264) select type(realization_base)
2265) class is(realization_subsurface_type)
2266) select case(option%iflowmode)
2267) case(RICHARDS_MODE)
2268) call RichardsComputeMassBalance(realization_base,sum_kg(1,:))
2269) case(TH_MODE)
2270) call THComputeMassBalance(realization_base,sum_kg(1,:))
2271) case(MIS_MODE)
2272) call MiscibleComputeMassBalance(realization_base,sum_kg(:,1))
2273) case(MPH_MODE)
2274) call MphaseComputeMassBalance(realization_base,sum_kg(:,:),sum_trapped(:))
2275) case(FLASH2_MODE)
2276) call Flash2ComputeMassBalance(realization_base,sum_kg(:,:),sum_trapped(:))
2277) case(IMS_MODE)
2278) call ImmisComputeMassBalance(realization_base,sum_kg(:,1))
2279) case(G_MODE)
2280) call GeneralComputeMassBalance(realization_base,sum_kg(:,:))
2281) case(TOIL_IMS_MODE)
2282) call TOilImsComputeMassBalance(realization_base,sum_kg(:,:))
2283) end select
2284) class default
2285) option%io_buffer = 'Unrecognized realization class in MassBalance().'
2286) call printErrMsg(option)
2287) end select
2288)
2289) int_mpi = option%nflowspec*option%nphase
2290) call MPI_Reduce(sum_kg,sum_kg_global, &
2291) int_mpi,MPI_DOUBLE_PRECISION,MPI_SUM, &
2292) option%io_rank,option%mycomm,ierr)
2293)
2294) if (option%iflowmode == MPH_MODE .or. option%iflowmode == FLASH2_MODE) then
2295) ! call MPI_Barrier(option%mycomm,ierr)
2296) int_mpi = option%nphase
2297) call MPI_Reduce(sum_trapped,sum_trapped_global, &
2298) int_mpi,MPI_DOUBLE_PRECISION,MPI_SUM, &
2299) option%io_rank,option%mycomm,ierr)
2300) endif
2301)
2302) if (option%myrank == option%io_rank) then
2303) select case(option%iflowmode)
2304) case(RICHARDS_MODE,IMS_MODE,MIS_MODE,G_MODE, &
2305) TH_MODE)
2306) do iphase = 1, option%nphase
2307) do ispec = 1, option%nflowspec
2308) write(fid,110,advance="no") sum_kg_global(ispec,iphase)
2309) enddo
2310) enddo
2311) case(TOIL_IMS_MODE)
2312) do iphase = 1, option%nphase
2313) write(fid,110,advance="no") sum_kg_global(iphase,1)
2314) enddo
2315) case(MPH_MODE,FLASH2_MODE)
2316) do iphase = 1, option%nphase
2317) do ispec = 1, option%nflowspec
2318) write(fid,110,advance="no") sum_kg_global(ispec,iphase)
2319) enddo
2320) write(fid,110,advance="no") sum_trapped_global(iphase)
2321) enddo
2322) end select
2323) endif
2324) endif
2325)
2326) if (option%ntrandof > 0) then
2327) sum_mol = 0.d0
2328) select type(realization_base)
2329) class is(realization_subsurface_type)
2330) call RTComputeMassBalance(realization_base,sum_mol)
2331) class default
2332) option%io_buffer = 'Unrecognized realization class in MassBalance().'
2333) call printErrMsg(option)
2334) end select
2335) int_mpi = option%nphase*option%ntrandof
2336) call MPI_Reduce(sum_mol,sum_mol_global,int_mpi, &
2337) MPI_DOUBLE_PRECISION,MPI_SUM, &
2338) option%io_rank,option%mycomm,ierr)
2339)
2340) if (option%myrank == option%io_rank) then
2341) ! do iphase = 1, option%nphase
2342) iphase = 1
2343) do icomp = 1, reaction%naqcomp
2344) if (reaction%primary_species_print(icomp)) then
2345) write(fid,110,advance="no") sum_mol_global(icomp,iphase)
2346) endif
2347) enddo
2348) ! enddo
2349) ! immobile species
2350) do icomp = 1, reaction%nimcomp
2351) if (reaction%immobile%print_me(icomp)) then
2352) write(fid,110,advance="no") &
2353) sum_mol_global(reaction%offset_immobile+icomp,iphase)
2354) endif
2355) enddo
2356) endif
2357)
2358) ! print out mineral contribution to mass balance
2359) if (option%mass_bal_detailed) then
2360) allocate(sum_mol_mnrl(realization_base%reaction%mineral%nkinmnrl))
2361) allocate(sum_mol_mnrl_global(realization_base%reaction%mineral%nkinmnrl))
2362)
2363) ! store integral over mineral volume fractions
2364) do imnrl = 1, reaction%mineral%nkinmnrl
2365) sum_mol_mnrl(imnrl) = 0.d0
2366) do local_id = 1, grid%nlmax
2367) ghosted_id = grid%nL2G(local_id)
2368) if (patch%imat(ghosted_id) <= 0) cycle
2369) sum_mol_mnrl(imnrl) = sum_mol_mnrl(imnrl) &
2370) + rt_auxvars(ghosted_id)%mnrl_volfrac(imnrl) &
2371) * material_auxvars(ghosted_id)%volume &
2372) / reaction%mineral%kinmnrl_molar_vol(imnrl)
2373) enddo
2374) enddo
2375)
2376) int_mpi = reaction%mineral%nkinmnrl
2377) call MPI_Reduce(sum_mol_mnrl,sum_mol_mnrl_global,int_mpi, &
2378) MPI_DOUBLE_PRECISION,MPI_SUM, &
2379) option%io_rank,option%mycomm,ierr)
2380) if (option%myrank == option%io_rank) then
2381) do imnrl = 1, reaction%mineral%nkinmnrl
2382) if (reaction%mineral%kinmnrl_print(imnrl)) then
2383) write(fid,110,advance="no") sum_mol_mnrl_global(imnrl)
2384) endif
2385) enddo
2386) endif
2387) deallocate(sum_mol_mnrl)
2388) deallocate(sum_mol_mnrl_global)
2389) endif
2390) endif
2391)
2392) coupler => patch%boundary_condition_list%first
2393) global_auxvars_bc_or_ss => patch%aux%Global%auxvars_bc
2394) if (option%ntrandof > 0) then
2395) rt_auxvars_bc_or_ss => patch%aux%RT%auxvars_bc
2396) endif
2397) bcs_done = PETSC_FALSE
2398) do
2399) if (.not.associated(coupler)) then
2400) if (bcs_done) then
2401) exit
2402) else
2403) bcs_done = PETSC_TRUE
2404) if (associated(patch%source_sink_list)) then
2405) coupler => patch%source_sink_list%first
2406) if (.not.associated(coupler)) exit
2407) global_auxvars_bc_or_ss => patch%aux%Global%auxvars_ss
2408) if (option%ntrandof > 0) then
2409) rt_auxvars_bc_or_ss => patch%aux%RT%auxvars_ss
2410) endif
2411) else
2412) exit
2413) endif
2414) endif
2415) endif
2416)
2417) offset = coupler%connection_set%offset
2418)
2419) if (option%nflowdof > 0) then
2420)
2421) #if 0
2422) ! compute the total area of the boundary condition
2423) if (.not.bcs_done) then
2424) sum_area = 0.d0
2425) do iconn = 1, coupler%connection_set%num_connections
2426) sum_area(1) = sum_area(1) + &
2427) coupler%connection_set%area(iconn)
2428) if (global_auxvars_bc_or_ss(offset+iconn)%sat(1) >= 0.5d0) then
2429) sum_area(2) = sum_area(2) + &
2430) coupler%connection_set%area(iconn)
2431) endif
2432) if (global_auxvars_bc_or_ss(offset+iconn)%sat(1) > 0.99d0) then
2433) sum_area(3) = sum_area(3) + &
2434) coupler%connection_set%area(iconn)
2435) endif
2436) sum_area(4) = sum_area(4) + &
2437) coupler%connection_set%area(iconn)* &
2438) global_auxvars_bc_or_ss(offset+iconn)%sat(1)
2439) enddo
2440)
2441) call MPI_Reduce(sum_area,sum_area_global, &
2442) FOUR_INTEGER_MPI,MPI_DOUBLE_PRECISION,MPI_SUM, &
2443) option%io_rank,option%mycomm,ierr)
2444)
2445) if (option%myrank == option%io_rank) then
2446) print *
2447) write(word,'(es16.6)') sum_area_global(1)
2448) print *, 'Total area in ' // trim(coupler%name) // &
2449) ' boundary condition: ' // trim(adjustl(word)) // ' m^2'
2450) write(word,'(es16.6)') sum_area_global(2)
2451) print *, 'Total half-saturated area in '// &
2452) trim(coupler%name) // &
2453) ' boundary condition: ' // trim(adjustl(word)) // ' m^2'
2454) write(word,'(es16.6)') sum_area_global(3)
2455) print *, 'Total saturated area in '// trim(coupler%name) // &
2456) ' boundary condition: ' // trim(adjustl(word)) // ' m^2'
2457) write(word,'(es16.6)') sum_area_global(4)
2458) print *, 'Total saturation-weighted area [=sum(saturation*area)] in '//&
2459) trim(coupler%name) // &
2460) ' boundary condition: ' // trim(adjustl(word)) // ' m^2'
2461) print *
2462) endif
2463) endif
2464) #endif
2465)
2466) select case(option%iflowmode)
2467) case(RICHARDS_MODE)
2468) ! print out cumulative H2O flux
2469) sum_kg = 0.d0
2470) do iconn = 1, coupler%connection_set%num_connections
2471) sum_kg = sum_kg + global_auxvars_bc_or_ss(offset+iconn)%mass_balance
2472) enddo
2473)
2474) int_mpi = option%nphase
2475) call MPI_Reduce(sum_kg,sum_kg_global, &
2476) int_mpi,MPI_DOUBLE_PRECISION,MPI_SUM, &
2477) option%io_rank,option%mycomm,ierr)
2478)
2479) if (option%myrank == option%io_rank) then
2480) ! change sign for positive in / negative out
2481) write(fid,110,advance="no") -sum_kg_global
2482) endif
2483)
2484) ! print out H2O flux
2485) sum_kg = 0.d0
2486) do iconn = 1, coupler%connection_set%num_connections
2487) sum_kg = sum_kg + global_auxvars_bc_or_ss(offset+iconn)%mass_balance_delta
2488) enddo
2489)
2490) ! mass_balance_delta units = delta kmol h2o; must convert to delta kg h2o
2491) sum_kg = sum_kg*FMWH2O
2492)
2493) int_mpi = option%nphase
2494) call MPI_Reduce(sum_kg,sum_kg_global, &
2495) int_mpi,MPI_DOUBLE_PRECISION,MPI_SUM, &
2496) option%io_rank,option%mycomm,ierr)
2497)
2498) if (option%myrank == option%io_rank) then
2499) ! change sign for positive in / negative out
2500) write(fid,110,advance="no") -sum_kg_global*output_option%tconv
2501) endif
2502)
2503) case(TH_MODE)
2504) ! print out cumulative H2O flux
2505) sum_kg = 0.d0
2506) do iconn = 1, coupler%connection_set%num_connections
2507) sum_kg = sum_kg + global_auxvars_bc_or_ss(offset+iconn)%mass_balance
2508) enddo
2509)
2510) int_mpi = option%nphase
2511) call MPI_Reduce(sum_kg,sum_kg_global, &
2512) int_mpi,MPI_DOUBLE_PRECISION,MPI_SUM, &
2513) option%io_rank,option%mycomm,ierr)
2514)
2515) if (option%myrank == option%io_rank) then
2516) ! change sign for positive in / negative out
2517) write(fid,110,advance="no") -sum_kg_global
2518) endif
2519)
2520) ! print out H2O flux
2521) sum_kg = 0.d0
2522) do iconn = 1, coupler%connection_set%num_connections
2523) sum_kg = sum_kg + global_auxvars_bc_or_ss(offset+iconn)%mass_balance_delta
2524) enddo
2525) ! mass_balance_delta units = delta kmol h2o; must convert to delta kg h2o
2526) sum_kg = sum_kg*FMWH2O
2527)
2528) int_mpi = option%nphase
2529) call MPI_Reduce(sum_kg,sum_kg_global, &
2530) int_mpi,MPI_DOUBLE_PRECISION,MPI_SUM, &
2531) option%io_rank,option%mycomm,ierr)
2532)
2533) if (option%myrank == option%io_rank) then
2534) ! change sign for positive in / negative out
2535) write(fid,110,advance="no") -sum_kg_global*output_option%tconv
2536) endif
2537)
2538) case(MIS_MODE)
2539) ! print out cumulative mixture flux
2540) sum_kg = 0.d0
2541) do icomp = 1, option%nflowspec
2542) do iconn = 1, coupler%connection_set%num_connections
2543) sum_kg(icomp,1) = sum_kg(icomp,1) + &
2544) global_auxvars_bc_or_ss(offset+iconn)%mass_balance(icomp,1)
2545) enddo
2546)
2547) if (icomp == 1) then
2548) sum_kg(icomp,1) = sum_kg(icomp,1)*FMWH2O
2549) else
2550) sum_kg(icomp,1) = sum_kg(icomp,1)*FMWGLYC
2551) endif
2552)
2553) int_mpi = option%nphase
2554) call MPI_Reduce(sum_kg(icomp,1),sum_kg_global(icomp,1), &
2555) int_mpi,MPI_DOUBLE_PRECISION,MPI_SUM, &
2556) option%io_rank,option%mycomm,ierr)
2557)
2558) if (option%myrank == option%io_rank) then
2559) ! change sign for positive in / negative out
2560) write(fid,110,advance="no") -sum_kg_global(icomp,1)
2561) endif
2562) enddo
2563)
2564) ! print out mixture flux
2565) sum_kg = 0.d0
2566) do icomp = 1, option%nflowspec
2567) do iconn = 1, coupler%connection_set%num_connections
2568) sum_kg(icomp,1) = sum_kg(icomp,1) + &
2569) global_auxvars_bc_or_ss(offset+iconn)%mass_balance_delta(icomp,1)
2570) enddo
2571)
2572) ! mass_balance_delta units = delta kmol h2o; must convert to delta kg h2o/glycol
2573) if (icomp == 1) then
2574) sum_kg(icomp,1) = sum_kg(icomp,1)*FMWH2O
2575) else
2576) sum_kg(icomp,1) = sum_kg(icomp,1)*FMWGLYC
2577) endif
2578)
2579) int_mpi = option%nphase
2580) call MPI_Reduce(sum_kg(icomp,1),sum_kg_global(icomp,1), &
2581) int_mpi,MPI_DOUBLE_PRECISION,MPI_SUM, &
2582) option%io_rank,option%mycomm,ierr)
2583)
2584) if (option%myrank == option%io_rank) then
2585) ! change sign for positive in / negative out
2586) write(fid,110,advance="no") -sum_kg_global(icomp,1)*output_option%tconv
2587) endif
2588) enddo
2589)
2590) case(MPH_MODE,FLASH2_MODE)
2591) ! print out cumulative H2O & CO2 fluxes in kmol and kmol/time
2592) sum_kg = 0.d0
2593) do icomp = 1, option%nflowspec
2594) do iconn = 1, coupler%connection_set%num_connections
2595) sum_kg(icomp,1) = sum_kg(icomp,1) + &
2596) global_auxvars_bc_or_ss(offset+iconn)%mass_balance(icomp,1)
2597) enddo
2598) int_mpi = option%nphase
2599) call MPI_Reduce(sum_kg(icomp,1),sum_kg_global(icomp,1), &
2600) int_mpi,MPI_DOUBLE_PRECISION,MPI_SUM, &
2601) option%io_rank,option%mycomm,ierr)
2602)
2603) if (option%myrank == option%io_rank) then
2604) ! change sign for positive in / negative out
2605) write(fid,110,advance="no") -sum_kg_global(icomp,1)
2606) endif
2607) enddo
2608)
2609) ! print out H2O & CO2 fluxes in kmol and kmol/time
2610) sum_kg = 0.d0
2611) do icomp = 1, option%nflowspec
2612) do iconn = 1, coupler%connection_set%num_connections
2613) sum_kg(icomp,1) = sum_kg(icomp,1) + &
2614) global_auxvars_bc_or_ss(offset+iconn)%mass_balance_delta(icomp,1)
2615) enddo
2616)
2617) ! mass_balance_delta units = delta kmol h2o; must convert to delta kg h2o
2618) ! sum_kg(icomp,1) = sum_kg(icomp,1)*FMWH2O ! <<---fix for multiphase!
2619)
2620) int_mpi = option%nphase
2621) call MPI_Reduce(sum_kg(icomp,1),sum_kg_global(icomp,1), &
2622) int_mpi,MPI_DOUBLE_PRECISION,MPI_SUM, &
2623) option%io_rank,option%mycomm,ierr)
2624)
2625) if (option%myrank == option%io_rank) then
2626) ! change sign for positive in / negative out
2627) write(fid,110,advance="no") -sum_kg_global(icomp,1)*output_option%tconv
2628) endif
2629) enddo
2630)
2631) case(IMS_MODE)
2632) ! print out cumulative H2O & CO2 fluxes
2633) sum_kg = 0.d0
2634) do icomp = 1, option%nflowspec
2635) do iconn = 1, coupler%connection_set%num_connections
2636) sum_kg(icomp,1) = sum_kg(icomp,1) + &
2637) global_auxvars_bc_or_ss(offset+iconn)%mass_balance(icomp,1)
2638) enddo
2639) int_mpi = option%nphase
2640) call MPI_Reduce(sum_kg(icomp,1),sum_kg_global(icomp,1), &
2641) int_mpi,MPI_DOUBLE_PRECISION,MPI_SUM, &
2642) option%io_rank,option%mycomm,ierr)
2643)
2644) if (option%myrank == option%io_rank) then
2645) ! change sign for positive in / negative out
2646) write(fid,110,advance="no") -sum_kg_global(icomp,1)
2647) endif
2648) enddo
2649)
2650) ! print out H2O & CO2 fluxes
2651) sum_kg = 0.d0
2652) do icomp = 1, option%nflowspec
2653) do iconn = 1, coupler%connection_set%num_connections
2654) sum_kg(icomp,1) = sum_kg(icomp,1) + &
2655) global_auxvars_bc_or_ss(offset+iconn)%mass_balance_delta(icomp,1)
2656) enddo
2657)
2658) ! mass_balance_delta units = delta kmol h2o; must convert to delta kg h2o
2659) ! sum_kg(icomp,1) = sum_kg(icomp,1)*FMWH2O ! <<---fix for multiphase!
2660)
2661) int_mpi = option%nphase
2662) call MPI_Reduce(sum_kg(icomp,1),sum_kg_global(icomp,1), &
2663) int_mpi,MPI_DOUBLE_PRECISION,MPI_SUM, &
2664) option%io_rank,option%mycomm,ierr)
2665)
2666) if (option%myrank == option%io_rank) then
2667) ! change sign for positive in / negative out
2668) write(fid,110,advance="no") -sum_kg_global(icomp,1)*output_option%tconv
2669) endif
2670) enddo
2671) case(G_MODE)
2672) ! print out cumulative H2O flux
2673) sum_kg = 0.d0
2674) do iconn = 1, coupler%connection_set%num_connections
2675) sum_kg = sum_kg + global_auxvars_bc_or_ss(offset+iconn)%mass_balance
2676) enddo
2677)
2678) int_mpi = option%nphase
2679) call MPI_Reduce(sum_kg,sum_kg_global, &
2680) int_mpi,MPI_DOUBLE_PRECISION,MPI_SUM, &
2681) option%io_rank,option%mycomm,ierr)
2682)
2683) if (option%myrank == option%io_rank) then
2684) ! change sign for positive in / negative out
2685) write(fid,110,advance="no") -sum_kg_global(:,1)
2686) endif
2687)
2688) ! print out H2O flux
2689) sum_kg = 0.d0
2690) do iconn = 1, coupler%connection_set%num_connections
2691) sum_kg(:,1) = sum_kg(:,1) + &
2692) global_auxvars_bc_or_ss(offset+iconn)%mass_balance_delta(:,1)
2693) enddo
2694) sum_kg(1,1) = sum_kg(1,1)*FMWH2O
2695) sum_kg(2,1) = sum_kg(2,1)*FMWAIR
2696)
2697) int_mpi = option%nphase
2698) call MPI_Reduce(sum_kg,sum_kg_global, &
2699) int_mpi,MPI_DOUBLE_PRECISION,MPI_SUM, &
2700) option%io_rank,option%mycomm,ierr)
2701)
2702) if (option%myrank == option%io_rank) then
2703) ! change sign for positive in / negative out
2704) write(fid,110,advance="no") -sum_kg_global(:,1)*output_option%tconv
2705) endif
2706) case(TOIL_IMS_MODE)
2707) ! print out cumulative H2O and Oil fluxes
2708) sum_kg = 0.d0
2709) do iconn = 1, coupler%connection_set%num_connections
2710) sum_kg = sum_kg + global_auxvars_bc_or_ss(offset+iconn)%mass_balance
2711) enddo
2712)
2713) int_mpi = option%nphase
2714) call MPI_Reduce(sum_kg,sum_kg_global, &
2715) int_mpi,MPI_DOUBLE_PRECISION,MPI_SUM, &
2716) option%io_rank,option%mycomm,ierr)
2717)
2718) if (option%myrank == option%io_rank) then
2719) ! change sign for positive in / negative out
2720) write(fid,110,advance="no") -sum_kg_global(:,1)
2721) endif
2722)
2723) ! print out H2O and oil fluxes
2724) sum_kg = 0.d0
2725) do iconn = 1, coupler%connection_set%num_connections
2726) sum_kg(:,1) = sum_kg(:,1) + &
2727) global_auxvars_bc_or_ss(offset+iconn)%mass_balance_delta(:,1)
2728) enddo
2729) sum_kg(1,1) = sum_kg(1,1)*toil_ims_fmw_comp(1)
2730) sum_kg(2,1) = sum_kg(2,1)*toil_ims_fmw_comp(2)
2731)
2732) int_mpi = option%nphase
2733) call MPI_Reduce(sum_kg,sum_kg_global, &
2734) int_mpi,MPI_DOUBLE_PRECISION,MPI_SUM, &
2735) option%io_rank,option%mycomm,ierr)
2736)
2737) if (option%myrank == option%io_rank) then
2738) ! change sign for positive in / negative out
2739) write(fid,110,advance="no") -sum_kg_global(:,1)*output_option%tconv
2740) endif
2741) end select
2742) endif
2743)
2744) if (option%ntrandof > 0) then
2745)
2746) ! print out cumulative boundary flux
2747) sum_mol = 0.d0
2748) do iconn = 1, coupler%connection_set%num_connections
2749) sum_mol = sum_mol + rt_auxvars_bc_or_ss(offset+iconn)%mass_balance
2750) enddo
2751)
2752) int_mpi = option%nphase*option%ntrandof
2753) call MPI_Reduce(sum_mol,sum_mol_global,int_mpi, &
2754) MPI_DOUBLE_PRECISION,MPI_SUM, &
2755) option%io_rank,option%mycomm,ierr)
2756)
2757) if (option%myrank == option%io_rank) then
2758) ! change sign for positive in / negative out
2759) do icomp = 1, reaction%naqcomp
2760) if (reaction%primary_species_print(icomp)) then
2761) write(fid,110,advance="no") -sum_mol_global(icomp,1)
2762) endif
2763) enddo
2764) endif
2765)
2766) ! print out boundary flux
2767) sum_mol = 0.d0
2768) do iconn = 1, coupler%connection_set%num_connections
2769) sum_mol = sum_mol + rt_auxvars_bc_or_ss(offset+iconn)%mass_balance_delta
2770) enddo
2771)
2772) int_mpi = option%nphase*option%ntrandof
2773) call MPI_Reduce(sum_mol,sum_mol_global,int_mpi, &
2774) MPI_DOUBLE_PRECISION,MPI_SUM, &
2775) option%io_rank,option%mycomm,ierr)
2776)
2777) if (option%myrank == option%io_rank) then
2778) ! change sign for positive in / negative out
2779) do icomp = 1, reaction%naqcomp
2780) if (reaction%primary_species_print(icomp)) then
2781) write(fid,110,advance="no") -sum_mol_global(icomp,1)* &
2782) output_option%tconv
2783) endif
2784) enddo
2785) endif
2786) endif
2787)
2788) coupler => coupler%next
2789) enddo
2790)
2791) ! Print the total mass in the specified regions (data)
2792) if (associated(output_option%mass_balance_region_list)) then
2793) cur_mbr => output_option%mass_balance_region_list
2794) do
2795) if (.not.associated(cur_mbr)) exit
2796) call PatchGetCompMassInRegion(cur_mbr%region_cell_ids, &
2797) cur_mbr%num_cells,patch,option,global_total_mass)
2798) write(fid,110,advance="no") global_total_mass
2799) cur_mbr => cur_mbr%next
2800) enddo
2801) endif
2802)
2803) #ifdef YE_FLUX
2804)
2805) !geh do offset = 1, 4
2806) !geh iconn = offset*20-1
2807)
2808) !TODO(ye): The flux will be calculated at the plane intersecting the top
2809) ! of the kth cell in the z-direction. You need to update this.
2810) k = 30
2811)
2812) if (option%nflowdof > 0) then
2813) ! really summation of moles, but we are hijacking the variable
2814) sum_mol_ye = 0.d0
2815) if (k-1 >= grid%structured_grid%lzs .and. &
2816) k-1 < grid%structured_grid%lze) then
2817) offset = (grid%structured_grid%ngx-1)*grid%structured_grid%nlyz + &
2818) (grid%structured_grid%ngy-1)*grid%structured_grid%nlxz
2819) do j = grid%structured_grid%lys, grid%structured_grid%lye-1
2820) do i = grid%structured_grid%lxs, grid%structured_grid%lxe-1
2821) iconn = offset + (i-grid%structured_grid%lxs)* &
2822) (grid%structured_grid%ngz-1) + &
2823) (j-grid%structured_grid%lys)* &
2824) grid%structured_grid%nlx* &
2825) (grid%structured_grid%ngz-1) + &
2826) k-grid%structured_grid%lzs+1
2827) !gehprint *, option%myrank, grid%nG2A(grid%internal_connection_set_list%first%id_up(iconn)), &
2828) !gehpatch%internal_fluxes(1:option%nflowdof,1,iconn), 'sum_mol_by_conn'
2829) sum_mol_ye(1:option%nflowdof) = sum_mol_ye(1:option%nflowdof) + &
2830) patch%internal_flow_fluxes(1:option%nflowdof,iconn)
2831) enddo
2832) enddo
2833) endif
2834) !geh int_mpi = option%nphase
2835) int_mpi = option%nflowdof
2836) !gehprint *, option%myrank, sum_mol_ye(1,1), 'sum_mol_ye'
2837) call MPI_Reduce(sum_mol_ye,sum_mol_global_ye, &
2838) int_mpi,MPI_DOUBLE_PRECISION,MPI_SUM, &
2839) option%io_rank,option%mycomm,ierr)
2840)
2841) if (option%myrank == option%io_rank) then
2842) ! change sign for positive in / negative out
2843) ! for mphase use:
2844) write(fid,110,advance="no") -sum_mol_global_ye(1:option%nflowdof)/option%flow_dt
2845)
2846) ! for Richards eqn. use:
2847) ! write(fid,110,advance="no") -sum_mol_global_ye(1:option%nflowdof)
2848) endif
2849) endif
2850)
2851) if (option%ntrandof > 0) then
2852)
2853) sum_mol = 0.d0
2854) sum_mol = sum_mol + patch%aux%RT%auxvars(iconn)%mass_balance
2855)
2856) int_mpi = option%nphase*option%ntrandof
2857) call MPI_Reduce(sum_mol,sum_mol_global,int_mpi, &
2858) MPI_DOUBLE_PRECISION,MPI_SUM, &
2859) option%io_rank,option%mycomm,ierr)
2860)
2861) if (option%myrank == option%io_rank) then
2862) do iphase = 1, option%nphase
2863) do icomp = 1, reaction%naqcomp
2864) if (reaction%primary_species_print(icomp)) then
2865) ! change sign for positive in / negative out
2866) write(fid,110,advance="no") -sum_mol_global(icomp,iphase)
2867) endif
2868) enddo
2869) enddo
2870) endif
2871) endif
2872) !geh enddo
2873) #endif
2874)
2875) if (option%myrank == option%io_rank) then
2876) write(fid,'(a)') ''
2877) close(fid)
2878) endif
2879)
2880) mass_balance_first = PETSC_FALSE
2881)
2882) end subroutine OutputMassBalance
2883)
2884) end module Output_Observation_module