factory_subsurface.F90 coverage: 100.00 %func 73.64 %block
1) module Factory_Subsurface_module
2)
3) use Simulation_Subsurface_class
4)
5) use PFLOTRAN_Constants_module
6)
7) implicit none
8)
9) private
10)
11) #include "petsc/finclude/petscsys.h"
12)
13) public :: SubsurfaceInitialize, &
14) SubsurfaceInitializePostPETSc, &
15) SubsurfaceJumpStart, &
16) ! move to init_subsurface
17) SubsurfaceReadFlowPM, &
18) SubsurfaceReadRTPM, &
19) SubsurfaceReadWasteFormPM, &
20) SubsurfaceReadUFDDecayPM
21)
22) contains
23)
24) ! ************************************************************************** !
25)
26) subroutine SubsurfaceInitialize(simulation)
27) !
28) ! Sets up PFLOTRAN subsurface simulation
29) !
30) ! Author: Glenn Hammond
31) ! Date: 06/10/13
32) !
33)
34) use WIPP_module
35) use Klinkenberg_module
36)
37) implicit none
38)
39) class(simulation_subsurface_type) :: simulation
40)
41) ! Modules that must be initialized
42) call WIPPInit()
43) call KlinkenbergInit()
44)
45) ! NOTE: PETSc must already have been initialized here!
46) call SubsurfaceInitializePostPetsc(simulation)
47)
48) end subroutine SubsurfaceInitialize
49)
50) ! ************************************************************************** !
51)
52) subroutine SubsurfaceInitializePostPetsc(simulation)
53) !
54) ! Sets up PFLOTRAN subsurface simulation
55) ! framework after to PETSc initialization
56) !
57) ! Author: Glenn Hammond
58) ! Date: 06/07/13
59) !
60)
61) use Option_module
62) use PM_Subsurface_Flow_class
63) use PM_Base_class
64) use PM_RT_class
65) use PM_Waste_Form_class
66) use PM_UFD_Decay_class
67) use PM_Auxiliary_class
68) use PMC_Subsurface_class
69) use PMC_Auxiliary_class
70) use PMC_Third_Party_class
71) use Timestepper_BE_class
72) use Realization_Subsurface_class
73) use Logging_module
74) use Simulation_Subsurface_class
75) use Solver_module
76) use Waypoint_module
77) use Init_Common_module
78) use Init_Subsurface_module
79) use Input_Aux_module
80) use String_module
81) use Checkpoint_module
82)
83) implicit none
84)
85) class(simulation_subsurface_type) :: simulation
86)
87) type(option_type), pointer :: option
88) class(pmc_subsurface_type), pointer :: pmc_subsurface
89) class(pmc_auxiliary_type), pointer :: auxiliary_process_model_coupler
90) class(pmc_third_party_type), pointer :: pmc_third_party
91) class(pm_subsurface_flow_type), pointer :: pm_flow
92) class(pm_rt_type), pointer :: pm_rt
93) class(pm_waste_form_type), pointer :: pm_waste_form
94) class(pm_ufd_decay_type), pointer :: pm_ufd_decay
95) class(pm_auxiliary_type), pointer :: pm_auxiliary
96) class(pm_base_type), pointer :: cur_pm, prev_pm
97) class(realization_subsurface_type), pointer :: realization
98) class(timestepper_BE_type), pointer :: timestepper
99) type(waypoint_list_type), pointer :: sync_waypoint_list
100) character(len=MAXSTRINGLENGTH) :: string
101)
102) option => simulation%option
103) ! process command line arguments specific to subsurface
104) call SubsurfInitCommandLineSettings(option)
105) nullify(pm_flow)
106) nullify(pm_rt)
107) nullify(pm_waste_form)
108) nullify(pm_ufd_decay)
109) nullify(pm_auxiliary)
110) cur_pm => simulation%process_model_list
111) do
112) if (.not.associated(cur_pm)) exit
113) select type(cur_pm)
114) class is(pm_subsurface_flow_type)
115) pm_flow => cur_pm
116) class is(pm_rt_type)
117) pm_rt => cur_pm
118) class is(pm_waste_form_type)
119) pm_waste_form => cur_pm
120) class is(pm_ufd_decay_type)
121) pm_ufd_decay => cur_pm
122) class is(pm_auxiliary_type)
123) pm_auxiliary => cur_pm
124) class default
125) option%io_buffer = &
126) 'PM Class unrecognized in SubsurfaceInitializePostPetsc.'
127) call printErrMsg(option)
128) end select
129) prev_pm => cur_pm
130) cur_pm => cur_pm%next
131) ! we must destroy the linkage between pms so that they are in independent
132) ! lists among pmcs
133) nullify(prev_pm%next)
134) enddo
135) call SubsurfaceSetFlowMode(pm_flow,option)
136) realization => RealizationCreate(option)
137) simulation%realization => realization
138) realization%output_option => simulation%output_option
139) simulation%waypoint_list_subsurface => WaypointListCreate()
140) if (associated(pm_flow)) then
141) pmc_subsurface => PMCSubsurfaceCreate()
142) pmc_subsurface%option => option
143) pmc_subsurface%checkpoint_option => simulation%checkpoint_option
144) pmc_subsurface%waypoint_list => simulation%waypoint_list_subsurface
145) pmc_subsurface%pm_list => pm_flow
146) pmc_subsurface%pm_ptr%pm => pm_flow
147) pmc_subsurface%realization => realization
148) ! set up logging stage
149) string = trim(pm_flow%name) // 'Flow'
150) call LoggingCreateStage(string,pmc_subsurface%stage)
151) ! timestepper => TimestepperBECreate()
152) ! timestepper%solver => SolverCreate()
153) ! simulation%flow_process_model_coupler%timestepper => timestepper
154) simulation%flow_process_model_coupler => pmc_subsurface
155) simulation%process_model_coupler_list => simulation%flow_process_model_coupler
156) nullify(pmc_subsurface)
157) endif
158) if (associated(pm_rt)) then
159) pmc_subsurface => PMCSubsurfaceCreate()
160) pmc_subsurface%name = 'PMCSubsurfaceTransport'
161) pmc_subsurface%option => option
162) pmc_subsurface%checkpoint_option => simulation%checkpoint_option
163) pmc_subsurface%waypoint_list => simulation%waypoint_list_subsurface
164) pmc_subsurface%pm_list => pm_rt
165) pmc_subsurface%pm_ptr%pm => pm_rt
166) pmc_subsurface%realization => realization
167) ! set up logging stage
168) string = 'Reactive Transport'
169) call LoggingCreateStage(string,pmc_subsurface%stage)
170) ! timestepper => TimestepperBECreate()
171) ! timestepper%solver => SolverCreate()
172) ! simulation%rt_process_model_coupler%timestepper => timestepper
173) simulation%rt_process_model_coupler => pmc_subsurface
174) if (.not.associated(simulation%process_model_coupler_list)) then
175) simulation%process_model_coupler_list => pmc_subsurface
176) else
177) simulation%flow_process_model_coupler%child => pmc_subsurface
178) endif
179) nullify(pmc_subsurface)
180) endif
181)
182) realization%input => InputCreate(IN_UNIT,option%input_filename,option)
183) call SubsurfaceReadRequiredCards(simulation)
184) call SubsurfaceReadInput(simulation)
185) if (associated(pm_waste_form)) then
186) string = 'WASTE_FORM_GENERAL'
187) call InputFindStringInFile(realization%input,option,string)
188) call InputFindStringErrorMsg(realization%input,option,string)
189) call pm_waste_form%Read(realization%input)
190) endif
191) if (associated(pm_ufd_decay)) then
192) string = 'UFD_DECAY'
193) call InputFindStringInFile(realization%input,option,string)
194) call InputFindStringErrorMsg(realization%input,option,string)
195) call pm_ufd_decay%Read(realization%input)
196) endif
197) call InputDestroy(realization%input)
198)
199) if (associated(pm_waste_form)) then
200) if (.not.associated(simulation%rt_process_model_coupler)) then
201) option%io_buffer = 'The Waste Form process model requires ' // &
202) 'reactive transport.'
203) call printErrMsg(option)
204) endif
205) pmc_third_party => PMCThirdPartyCreate()
206) pmc_third_party%option => option
207) pmc_third_party%checkpoint_option => simulation%checkpoint_option
208) pmc_third_party%pm_list => pm_waste_form
209) pmc_third_party%pm_ptr%pm => pm_waste_form
210) pmc_third_party%realization => realization
211) ! set up logging stage
212) string = 'WASTE_FORM_GENERAL'
213) call LoggingCreateStage(string,pmc_third_party%stage)
214) if (associated(simulation%rt_process_model_coupler%child)) then
215) simulation%rt_process_model_coupler%child%peer => pmc_third_party
216) else
217) simulation%rt_process_model_coupler%child => pmc_third_party
218) endif
219) nullify(pmc_third_party)
220) endif
221)
222) if (associated(pm_ufd_decay)) then
223) if (.not.associated(simulation%rt_process_model_coupler)) then
224) option%io_buffer = 'The UFD Decay process model requires reactive ' // &
225) 'transport.'
226) call printErrMsg(option)
227) endif
228) pmc_third_party => PMCThirdPartyCreate()
229) pmc_third_party%option => option
230) pmc_third_party%checkpoint_option => simulation%checkpoint_option
231) pmc_third_party%waypoint_list => simulation%waypoint_list_subsurface
232) pmc_third_party%pm_list => pm_ufd_decay
233) pmc_third_party%pm_ptr%pm => pm_ufd_decay
234) pmc_third_party%realization => realization
235) ! set up logging stage
236) string = 'UFD_DECAY'
237) call LoggingCreateStage(string,pmc_third_party%stage)
238) if (associated(simulation%rt_process_model_coupler%child)) then
239) simulation%rt_process_model_coupler%child%peer => pmc_third_party
240) else
241) simulation%rt_process_model_coupler%child => pmc_third_party
242) endif
243) nullify(pmc_third_party)
244) endif
245)
246) if (associated(pm_auxiliary)) then
247) string = 'salinity'
248) if (StringCompareIgnoreCase(pm_auxiliary%ctype,string)) then
249) if (associated(simulation%rt_process_model_coupler)) then
250) auxiliary_process_model_coupler => PMCAuxiliaryCreate()
251) simulation%rt_process_model_coupler%peer => auxiliary_process_model_coupler
252) pm_auxiliary%realization => realization
253) auxiliary_process_model_coupler%pm_list => pm_auxiliary
254) auxiliary_process_model_coupler%pm_aux => pm_auxiliary
255) auxiliary_process_model_coupler%option => option
256) else
257) option%io_buffer = 'Reactive transport must be included in the &
258) &SIMULATION block in order to use the SALINITY process model.'
259) call printErrMsg(option)
260) endif
261) endif
262) endif
263)
264) ! SubsurfaceInitSimulation() must be called after pmc linkages are set above.
265) call SubsurfaceInitSimulation(simulation)
266)
267) ! create sync waypoint list to be used a few lines below
268) sync_waypoint_list => &
269) WaypointCreateSyncWaypointList(simulation%waypoint_list_subsurface)
270) ! merge in outer waypoints (e.g. checkpoint times)
271) call WaypointListCopyAndMerge(simulation%waypoint_list_subsurface, &
272) simulation%waypoint_list_outer,option)
273) ! add sync waypoints into outer list
274) call WaypointListMerge(simulation%waypoint_list_outer,sync_waypoint_list, &
275) option)
276) ! add in periodic time waypoints for checkpointing. these will not appear
277) ! in the outer list
278) call CheckpointPeriodicTimeWaypoints(simulation%checkpoint_option, &
279) simulation%waypoint_list_subsurface)
280)
281) ! clean up waypoints
282) if (.not.option%steady_state) then
283) ! fill in holes in waypoint data
284) call WaypointListFillIn(simulation%waypoint_list_subsurface,option)
285) call WaypointListRemoveExtraWaypnts(simulation%waypoint_list_subsurface, &
286) option)
287) endif
288)
289) ! debugging output
290) if (realization%debug%print_couplers) then
291) call InitCommonVerifyAllCouplers(realization)
292) endif
293) if (realization%debug%print_waypoints) then
294) call WaypointListPrint(simulation%waypoint_list_subsurface,option, &
295) realization%output_option)
296) endif
297)
298) call SubsurfaceJumpStart(simulation)
299) ! set first process model coupler as the master
300) simulation%process_model_coupler_list%is_master = PETSC_TRUE
301)
302) end subroutine SubsurfaceInitializePostPetsc
303)
304) ! ************************************************************************** !
305)
306) subroutine SubsurfInitCommandLineSettings(option)
307) !
308) ! Initializes PFLTORAN subsurface output
309) ! filenames, etc.
310) !
311) ! Author: Glenn Hammond
312) ! Date: 06/06/13
313) !
314)
315) use Option_module
316) use Input_Aux_module
317)
318) implicit none
319)
320) type(option_type) :: option
321)
322) character(len=MAXSTRINGLENGTH) :: string
323) PetscBool :: option_found
324) PetscBool :: bool_flag
325)
326) string = '-multisimulation'
327) call InputGetCommandLineTruth(string,bool_flag,option_found,option)
328) if (option_found) then
329) option%subsurface_simulation_type = MULTISIMULATION_SIM_TYPE
330) endif
331)
332) string = '-stochastic'
333) call InputGetCommandLineTruth(string,bool_flag,option_found,option)
334) if (option_found) then
335) option%subsurface_simulation_type = STOCHASTIC_SIM_TYPE
336) endif
337)
338) end subroutine SubsurfInitCommandLineSettings
339)
340) ! ************************************************************************** !
341)
342) subroutine SubsurfaceSetFlowMode(pm_flow,option)
343) !
344) ! Sets the flow mode (richards, vadose, mph, etc.)
345) !
346) ! Author: Glenn Hammond
347) ! Date: 10/26/07
348) !
349)
350) use Option_module
351) use PM_Subsurface_Flow_class
352) use PM_Base_class
353) use PM_Flash2_class
354) use PM_General_class
355) use PM_Immis_class
356) use PM_Miscible_class
357) use PM_Mphase_class
358) use PM_Richards_class
359) use PM_TH_class
360) use PM_TOilIms_class
361)
362) implicit none
363)
364) type(option_type) :: option
365) class(pm_subsurface_flow_type), pointer :: pm_flow
366)
367) if (.not.associated(pm_flow)) then
368) option%nphase = 1
369) option%liquid_phase = 1
370) ! assume default isothermal when only transport
371) option%use_isothermal = PETSC_TRUE
372) return
373) endif
374)
375) select type(pm_flow)
376) class is (pm_flash2_type)
377) option%iflowmode = FLASH2_MODE
378) option%nphase = 2
379) option%liquid_phase = 1
380) option%gas_phase = 2
381) option%nflowdof = 3
382) option%nflowspec = 2
383) option%itable = 2
384) option%use_isothermal = PETSC_FALSE
385) option%water_id = 1
386) option%air_id = 2
387) class is (pm_general_type)
388) option%iflowmode = G_MODE
389) option%nphase = 2
390) option%liquid_phase = 1 ! liquid_pressure
391) option%gas_phase = 2 ! gas_pressure
392)
393) option%air_pressure_id = 3
394) option%capillary_pressure_id = 4
395) option%vapor_pressure_id = 5
396) option%saturation_pressure_id = 6
397)
398) option%water_id = 1
399) option%air_id = 2
400) option%energy_id = 3
401)
402) option%nflowdof = 3
403) option%nflowspec = 2
404) option%use_isothermal = PETSC_FALSE
405) class is (pm_toil_ims_type)
406) option%iflowmode = TOIL_IMS_MODE
407) option%nphase = 2
408) option%liquid_phase = 1 ! liquid_pressure
409) option%oil_phase = 2 ! oil_pressure
410)
411) option%capillary_pressure_id = 3 ! capillary pressure
412)
413) option%nflowdof = 3
414) !two species (H2O,OIL): each present only in its own rich phase
415) option%nflowspec = 2
416)
417) option%water_id = 1
418) option%oil_id = 2
419) option%energy_id = 3
420)
421) option%use_isothermal = PETSC_FALSE
422) class is (pm_immis_type)
423) option%iflowmode = IMS_MODE
424) option%nphase = 2
425) option%liquid_phase = 1
426) option%gas_phase = 2
427) option%nflowdof = 3
428) option%nflowspec = 2
429) option%itable = 2
430) option%io_buffer = 'Material AuxVars must be refactored for IMMIS.'
431) call printErrMsg(option)
432) class is (pm_miscible_type)
433) option%iflowmode = MIS_MODE
434) option%nphase = 1
435) option%liquid_phase = 1
436) option%gas_phase = 2
437) option%nflowdof = 2
438) option%nflowspec = 2
439) option%io_buffer = 'Material AuxVars must be refactored for MISCIBLE.'
440) call printErrMsg(option)
441) class is (pm_mphase_type)
442) option%iflowmode = MPH_MODE
443) option%nphase = 2
444) option%liquid_phase = 1
445) option%gas_phase = 2
446) option%nflowdof = 3
447) option%nflowspec = 2
448) option%itable = 2 ! read CO2DATA0.dat
449) ! option%itable = 1 ! create CO2 database: co2data.dat
450) option%use_isothermal = PETSC_FALSE
451) option%water_id = 1
452) option%air_id = 2
453) class is (pm_richards_type)
454) option%iflowmode = RICHARDS_MODE
455) option%nphase = 1
456) option%liquid_phase = 1
457) option%nflowdof = 1
458) option%nflowspec = 1
459) option%use_isothermal = PETSC_TRUE
460) class is (pm_th_type)
461) option%iflowmode = TH_MODE
462) option%nphase = 1
463) option%liquid_phase = 1
464) option%gas_phase = 2
465) option%nflowdof = 2
466) option%nflowspec = 1
467) option%use_isothermal = PETSC_FALSE
468) option%flow%store_fluxes = PETSC_TRUE
469) class default
470) end select
471)
472) end subroutine SubsurfaceSetFlowMode
473)
474) ! ************************************************************************** !
475)
476) subroutine SubsurfaceReadFlowPM(input, option, pm)
477) !
478) ! Author: Glenn Hammond
479) ! Date: 06/11/13
480) !
481) use Input_Aux_module
482) use Option_module
483) use String_module
484)
485) use PMC_Base_class
486)
487) use PM_Base_class
488) use PM_Flash2_class
489) use PM_General_class
490) use PM_Immis_class
491) use PM_Miscible_class
492) use PM_Mphase_class
493) use PM_Richards_class
494) use PM_TH_class
495) use PM_TOilIms_class
496)
497) use Init_Common_module
498)
499) use General_module
500)
501) implicit none
502)
503) type(input_type), pointer :: input
504) type(option_type), pointer :: option
505) class(pm_base_type), pointer :: pm
506)
507) character(len=MAXWORDLENGTH) :: word
508) character(len=MAXSTRINGLENGTH) :: error_string
509)
510) error_string = 'SIMULATION,PROCESS_MODELS,SUBSURFACE_FLOW'
511)
512) nullify(pm)
513) word = ''
514) do
515) call InputReadPflotranString(input,option)
516) if (InputCheckExit(input,option)) exit
517) call InputReadWord(input,option,word,PETSC_FALSE)
518) call StringToUpper(word)
519) select case(word)
520) case('MODE')
521) call InputReadWord(input,option,word,PETSC_FALSE)
522) call InputErrorMsg(input,option,'mode',error_string)
523) call StringToUpper(word)
524) select case(word)
525) case('GENERAL')
526) pm => PMGeneralCreate()
527) case('MPHASE')
528) pm => PMMphaseCreate()
529) case('FLASH2')
530) pm => PMFlash2Create()
531) case('IMS','IMMIS','THS')
532) pm => PMImmisCreate()
533) case('MIS','MISCIBLE')
534) pm => PMMiscibleCreate()
535) case('RICHARDS')
536) pm => PMRichardsCreate()
537) case('TH')
538) pm => PMTHCreate()
539) case('TOIL_IMS')
540) pm => PMToilImsCreate()
541) case default
542) error_string = trim(error_string) // ',MODE'
543) call InputKeywordUnrecognized(word,error_string,option)
544) end select
545) pm%option => option
546) case('OPTIONS')
547) if (.not.associated(pm)) then
548) option%io_buffer = 'MODE keyword must be read first under ' // &
549) trim(error_string)
550) call printErrMsg(option)
551) endif
552) call pm%Read(input)
553) case default
554) error_string = trim(error_string) // ',SUBSURFACE_FLOW'
555) call InputKeywordUnrecognized(word,error_string,option)
556) end select
557) enddo
558)
559) if (.not.associated(pm)) then
560) option%io_buffer = 'A flow MODE (card) must be included in the ' // &
561) 'SUBSURFACE_FLOW block in ' // trim(error_string) // '.'
562) call printErrMsg(option)
563) endif
564)
565) end subroutine SubsurfaceReadFlowPM
566)
567) ! ************************************************************************** !
568)
569) subroutine SubsurfaceReadRTPM(input, option, pm)
570) !
571) ! Author: Glenn Hammond
572) ! Date: 06/11/13
573) !
574) use Input_Aux_module
575) use Option_module
576) use String_module
577)
578) use PMC_Base_class
579) use PM_Base_class
580) use PM_RT_class
581)
582) use Init_Common_module
583)
584) implicit none
585)
586) type(input_type), pointer :: input
587) type(option_type), pointer :: option
588) class(pm_base_type), pointer :: pm
589)
590) character(len=MAXWORDLENGTH) :: word
591) character(len=MAXSTRINGLENGTH) :: error_string
592)
593) error_string = 'SIMULATION,PROCESS_MODELS,SUBSURFACE_TRANSPORT'
594)
595) pm => PMRTCreate()
596) pm%option => option
597)
598) call pm%Read(input)
599)
600) end subroutine SubsurfaceReadRTPM
601)
602) ! ************************************************************************** !
603)
604) subroutine SubsurfaceReadWasteFormPM(input, option, pm)
605) !
606) ! Author: Glenn Hammond
607) ! Date: 06/11/13
608) !
609) use Input_Aux_module
610) use Option_module
611) use String_module
612)
613) use PM_Base_class
614) use PM_Waste_Form_class
615)
616) implicit none
617)
618) type(input_type), pointer :: input
619) type(option_type), pointer :: option
620) class(pm_base_type), pointer :: pm
621)
622) character(len=MAXWORDLENGTH) :: word
623) character(len=MAXSTRINGLENGTH) :: error_string
624)
625) error_string = 'SIMULATION,PROCESS_MODELS,WASTE_FORM'
626)
627) word = ''
628) do
629) call InputReadPflotranString(input,option)
630) if (InputCheckExit(input,option)) exit
631) call InputReadWord(input,option,word,PETSC_FALSE)
632) call StringToUpper(word)
633) select case(word)
634) case('TYPE')
635) call InputReadWord(input,option,word,PETSC_FALSE)
636) call InputErrorMsg(input,option,'mode',error_string)
637) call StringToUpper(word)
638) select case(word)
639) case('GENERAL')
640) pm => PMWFCreate()
641) case default
642) option%io_buffer = 'WASTE FORM type ' // trim(word) // &
643) ' not recognized. Only TYPE GENERAL currently supported. &
644) & TYPE GLASS or TYPE FMDM no longer supported.'
645) call printErrMsg(option)
646) end select
647) case default
648) option%io_buffer = 'Keyword ' // trim(word) // &
649) ' not recognized for the ' // trim(error_string) // ' block.'
650) call printErrMsg(option)
651) end select
652) enddo
653)
654) if (.not.associated(pm)) then
655) option%io_buffer = 'TYPE card missing in ' // trim(error_string)
656) call printErrMsg(option)
657) endif
658)
659) pm%option => option
660)
661) end subroutine SubsurfaceReadWasteFormPM
662)
663) ! ************************************************************************** !
664)
665) subroutine SubsurfaceReadUFDDecayPM(input, option, pm)
666) !
667) ! Author: Glenn Hammond
668) ! Date: 06/11/13
669) !
670) use Input_Aux_module
671) use Option_module
672) use String_module
673)
674) use PM_Base_class
675) use PM_UFD_Decay_class
676)
677) implicit none
678)
679) type(input_type), pointer :: input
680) type(option_type), pointer :: option
681) class(pm_base_type), pointer :: pm
682)
683) character(len=MAXWORDLENGTH) :: word
684) character(len=MAXSTRINGLENGTH) :: error_string
685)
686) error_string = 'SIMULATION,PROCESS_MODELS,UFD_DECAY'
687)
688) pm => PMUFDDecayCreate()
689) pm%option => option
690)
691) word = ''
692) do
693) call InputReadPflotranString(input,option)
694) if (InputCheckExit(input,option)) exit
695) call InputReadWord(input,option,word,PETSC_FALSE)
696) call StringToUpper(word)
697) select case(word)
698) case default
699) option%io_buffer = 'Keyword ' // trim(word) // &
700) ' not recognized for the ' // trim(error_string) // ' block.'
701) call printErrMsg(option)
702) end select
703) enddo
704)
705) end subroutine SubsurfaceReadUFDDecayPM
706)
707) ! ************************************************************************** !
708)
709) subroutine SubsurfaceInitSimulation(simulation)
710) !
711) ! Author: Glenn Hammond
712) ! Date: 06/11/13
713) !
714)
715) use Realization_Subsurface_class
716) use Realization_Base_class
717) use Discretization_module
718) use Option_module
719) use Output_module, only : Output
720) use Output_Aux_module
721)
722)
723) use Global_module
724) use Init_Subsurface_module
725) use Init_Subsurface_Flow_module
726) use Init_Subsurface_Tran_module
727) use Init_Common_module
728) use Waypoint_module
729) use Strata_module
730) use Regression_module
731)
732) use PMC_Subsurface_class
733) use PMC_Auxiliary_class
734) use PMC_Base_class
735) use PM_Base_class
736) use PM_Base_Pointer_module
737) use PM_Subsurface_Flow_class
738) use PM_Auxiliary_class
739)
740) use PM_General_class
741) use PM_Richards_class
742) use PM_TH_class
743) use PM_RT_class
744) use PM_Waste_Form_class
745) use PM_UFD_Decay_class
746) use PM_TOilIms_class
747)
748) use Timestepper_BE_class
749)
750) implicit none
751)
752) #include "petsc/finclude/petscsnes.h"
753)
754) class(simulation_subsurface_type) :: simulation
755)
756) class(pmc_subsurface_type), pointer :: flow_process_model_coupler
757) class(pmc_subsurface_type), pointer :: tran_process_model_coupler
758) class(pmc_auxiliary_type), pointer :: auxiliary_process_model_coupler
759) class(pmc_base_type), pointer :: cur_process_model_coupler
760) class(pmc_base_type), pointer :: cur_process_model_coupler_top
761) class(pm_base_type), pointer :: cur_process_model
762) class(pm_auxiliary_type), pointer :: pm_aux
763)
764) class(realization_subsurface_type), pointer :: realization
765) type(option_type), pointer :: option
766) character(len=MAXSTRINGLENGTH) :: string
767) SNESLineSearch :: linesearch
768) PetscInt :: ndof
769) PetscBool, allocatable :: dof_is_active(:)
770) PetscErrorCode :: ierr
771)
772) realization => simulation%realization
773) option => realization%option
774)
775) ! begin from old Init()
776) call SubsurfaceSetupRealization(simulation)
777) call InitCommonAddOutputWaypoints(option,simulation%output_option, &
778) simulation%waypoint_list_subsurface)
779)
780) !TODO(geh): refactor
781) if (associated(simulation%flow_process_model_coupler)) then
782) if (associated(simulation%flow_process_model_coupler%timestepper)) then
783) simulation%flow_process_model_coupler%timestepper%cur_waypoint => &
784) simulation%waypoint_list_subsurface%first
785) endif
786) endif
787) if (associated(simulation%rt_process_model_coupler)) then
788) if (associated(simulation%rt_process_model_coupler%timestepper)) then
789) simulation%rt_process_model_coupler%timestepper%cur_waypoint => &
790) simulation%waypoint_list_subsurface%first
791) endif
792) endif
793)
794) !TODO(geh): refactor
795) ! initialize global auxiliary variable object
796) call GlobalSetup(realization)
797)
798) ! always call the flow side since a velocity field still has to be
799) ! set if no flow exists
800) call InitSubsurfFlowSetupRealization(realization)
801) if (option%ntrandof > 0) then
802) call InitSubsurfTranSetupRealization(realization)
803) endif
804) ! InitSubsurfaceSetupZeroArray must come after InitSubsurfaceXXXRealization
805) call InitSubsurfaceSetupZeroArrays(realization)
806) call OutputVariableAppendDefaults(realization%output_option% &
807) output_snap_variable_list,option)
808)
809) if (option%nflowdof > 0) then
810) select type(ts => simulation%flow_process_model_coupler%timestepper)
811) class is (timestepper_BE_type)
812) call InitSubsurfFlowSetupSolvers(realization,ts%convergence_context, &
813) ts%solver)
814) end select
815) endif
816) if (option%ntrandof > 0) then
817) select type(ts => simulation%rt_process_model_coupler%timestepper)
818) class is (timestepper_BE_type)
819) call InitSubsurfTranSetupSolvers(realization,ts%convergence_context, &
820) ts%solver)
821) end select
822) endif
823) call RegressionCreateMapping(simulation%regression,realization)
824) ! end from old Init()
825)
826) call DiscretizationPrintInfo(realization%discretization, &
827) realization%patch%grid,option)
828)
829) !----------------------------------------------------------------------------!
830) ! This section for setting up new process model approach
831) !----------------------------------------------------------------------------!
832)
833) if (StrataEvolves(realization%patch%strata_list)) then
834) auxiliary_process_model_coupler => PMCAuxiliaryCreate()
835) allocate(pm_aux)
836) call PMAuxiliaryInit(pm_aux)
837) string = 'EVOLVING_STRATA'
838) call PMAuxiliarySetFunctionPointer(pm_aux,string)
839) pm_aux%realization => realization
840) pm_aux%option => option
841) auxiliary_process_model_coupler%pm_list => pm_aux
842) auxiliary_process_model_coupler%pm_aux => pm_aux
843) auxiliary_process_model_coupler%option => option
844) ! place the material process model as %peer for the top pmc
845) simulation%process_model_coupler_list%peer => &
846) auxiliary_process_model_coupler
847) endif
848)
849) ! For each ProcessModel, set:
850) ! - realization (subsurface or surface),
851) ! - stepper (flow/trans/surf_flow),
852) ! - SNES functions (Residual/Jacobain), or TS function (RHSFunction)
853)
854) cur_process_model_coupler_top => simulation%process_model_coupler_list
855) ! the following recursive subroutine will also call each pmc child
856) ! and each pms's peers
857) if (associated(cur_process_model_coupler_top)) then
858) call SetUpPMApproach(cur_process_model_coupler_top,simulation)
859) endif
860)
861) ! point the top process model coupler to Output
862) simulation%process_model_coupler_list%Output => Output
863)
864) end subroutine SubsurfaceInitSimulation
865)
866) ! ************************************************************************** !
867)
868) recursive subroutine SetUpPMApproach(pmc,simulation)
869) !
870) ! Loops through all of the PMC's recursively and sets their realization,
871) ! timestepper, and solver.
872) !
873) ! Author: Jenn Frederick, SNL
874) ! Date: 04/04/2016
875) !
876) use PMC_Base_class
877) use PMC_Subsurface_class
878) use PM_Base_Pointer_module
879) use PM_Base_class
880) use PM_Subsurface_Flow_class
881) use PM_General_class
882) use PM_Richards_class
883) use PM_TH_class
884) use PM_RT_class
885) use PM_Waste_Form_class
886) use PM_UFD_Decay_class
887) use PM_TOilIms_class
888) use Option_module
889) use Simulation_Subsurface_class
890) use Realization_Subsurface_class
891) use Timestepper_BE_class
892)
893) implicit none
894)
895) #include "petsc/finclude/petscsnes.h"
896)
897) class(pmc_base_type), pointer :: pmc
898) class(simulation_subsurface_type) :: simulation
899)
900) class(realization_subsurface_type), pointer :: realization
901) class(pm_base_type), pointer :: cur_pm
902) type(option_type), pointer :: option
903) SNESLineSearch :: linesearch
904) PetscErrorCode :: ierr
905)
906) realization => simulation%realization
907) option => realization%option
908)
909) if (.not.associated(pmc)) return
910)
911) pmc%waypoint_list => simulation%waypoint_list_subsurface
912)
913) ! loop through this pmc's process models:
914) cur_pm => pmc%pm_list
915) do
916) if (.not.associated(cur_pm)) exit
917) ! set realization
918) select type(cur_pm)
919) !-----------------------------------
920) class is(pm_rt_type)
921) if (.not.associated(realization%reaction)) then
922) option%io_buffer = 'SUBSURFACE_TRANSPORT specified as a ' // &
923) 'process model without a corresponding CHEMISTRY block.'
924) call printErrMsg(option)
925) endif
926) call cur_pm%PMRTSetRealization(realization)
927) !-----------------------------------
928) class is(pm_subsurface_flow_type)
929) call cur_pm%PMSubsurfaceFlowSetRealization(realization)
930) !-----------------------------------
931) class is(pm_waste_form_type)
932) call cur_pm%PMWFSetRealization(realization)
933) !-----------------------------------
934) class is(pm_ufd_decay_type)
935) call cur_pm%PMUFDDecaySetRealization(realization)
936) !-----------------------------------
937) end select
938) ! set time stepper
939) select type(cur_pm)
940) !-----------------------------------
941) class is(pm_subsurface_flow_type)
942) pmc%timestepper%dt = option%flow_dt
943) !-----------------------------------
944) class is(pm_rt_type)
945) pmc%timestepper%dt = option%tran_dt
946) !-----------------------------------
947) end select
948) cur_pm%output_option => simulation%output_option
949) call cur_pm%Setup()
950) if (associated(pmc%timestepper)) then
951) select type(ts => pmc%timestepper)
952) !----------------------------------
953) class is(timestepper_BE_type)
954) call SNESGetLineSearch(ts%solver%snes,linesearch, &
955) ierr);CHKERRQ(ierr)
956) ! Post
957) select type(cur_pm)
958) !-----------------------------------
959) class is(pm_subsurface_flow_type)
960) if (cur_pm%check_post_convergence) then
961) call SNESLineSearchSetPostCheck(linesearch,PMCheckUpdatePost, &
962) pmc%pm_ptr,ierr);CHKERRQ(ierr)
963) !geh: it is possible that the other side has not been set
964) cur_pm%check_post_convergence = PETSC_TRUE
965) endif
966) !------------------------------------
967) class is(pm_rt_type)
968) if (cur_pm%print_EKG .or. option%use_mc .or. &
969) cur_pm%check_post_convergence) then
970) call SNESLineSearchSetPostCheck(linesearch,PMCheckUpdatePost, &
971) pmc%pm_ptr,ierr);CHKERRQ(ierr)
972) if (cur_pm%print_EKG) then
973) cur_pm%check_post_convergence = PETSC_TRUE
974) endif
975) endif
976) !-------------------------------------
977) end select
978) ! Pre
979) select type(cur_pm)
980) !-------------------------------------
981) class is(pm_richards_type)
982) if (Initialized(cur_pm%pressure_dampening_factor) .or. &
983) Initialized(cur_pm%saturation_change_limit)) then
984) call SNESLineSearchSetPreCheck(linesearch,PMCheckUpdatePre, &
985) pmc%pm_ptr,ierr);CHKERRQ(ierr)
986) endif
987) !-------------------------------------
988) class is(pm_general_type)
989) call SNESLineSearchSetPreCheck(linesearch,PMCheckUpdatePre, &
990) pmc%pm_ptr,ierr);CHKERRQ(ierr)
991) !-------------------------------------
992) class is(pm_toil_ims_type)
993) call SNESLineSearchSetPreCheck(linesearch,PMCheckUpdatePre, &
994) pmc%pm_ptr,ierr);CHKERRQ(ierr)
995) !-------------------------------------
996) class is(pm_th_type)
997) if (Initialized(cur_pm%pressure_dampening_factor) .or. &
998) Initialized(cur_pm%pressure_change_limit) .or. &
999) Initialized(cur_pm%temperature_change_limit)) then
1000) call SNESLineSearchSetPreCheck(linesearch,PMCheckUpdatePre, &
1001) pmc%pm_ptr,ierr);CHKERRQ(ierr)
1002) endif
1003) !-------------------------------------
1004) class is(pm_rt_type)
1005) if (realization%reaction%check_update) then
1006) call SNESLineSearchSetPreCheck(linesearch,PMCheckUpdatePre, &
1007) pmc%pm_ptr,ierr);CHKERRQ(ierr)
1008) endif
1009) !-------------------------------------
1010) end select
1011) !----------------------------------
1012) end select
1013) endif ! associated(pmc%timestepper)
1014) cur_pm => cur_pm%next
1015) enddo
1016) call pmc%SetupSolvers()
1017)
1018) ! call this function for this pmc's child
1019) if (associated(pmc%child)) then
1020) call SetUpPMApproach(pmc%child,simulation)
1021) endif
1022)
1023) ! call this function for this pmc's peer
1024) if (associated(pmc%peer)) then
1025) call SetUpPMApproach(pmc%peer,simulation)
1026) endif
1027)
1028)
1029) end subroutine SetUpPMApproach
1030)
1031) ! ************************************************************************** !
1032)
1033) subroutine SubsurfaceSetupRealization(simulation)
1034) !
1035) ! Initializes material property data structres and assign them to the domain.
1036) !
1037) ! Author: Glenn Hammond
1038) ! Date: 12/04/14
1039) !
1040) use Init_Subsurface_module
1041) use Simulation_Subsurface_class
1042) use Realization_Subsurface_class
1043) use Option_module
1044) use Logging_module
1045) use Waypoint_module
1046) use Init_Common_module
1047) use Reaction_Aux_module, only : ACT_COEF_FREQUENCY_OFF
1048) use Reaction_Database_module
1049) use EOS_Water_module
1050) use Dataset_module
1051) use Patch_module
1052)
1053) implicit none
1054)
1055) class(simulation_subsurface_type) :: simulation
1056)
1057) class(realization_subsurface_type), pointer :: realization
1058) type(option_type), pointer :: option
1059) PetscReal :: dum1
1060) PetscErrorCode :: ierr
1061)
1062) realization => simulation%realization
1063) option => realization%option
1064)
1065) call PetscLogEventBegin(logging%event_setup,ierr);CHKERRQ(ierr)
1066)
1067) ! initialize reference density
1068) if (option%reference_water_density < 1.d-40) then
1069) call EOSWaterDensity(option%reference_temperature, &
1070) option%reference_pressure, &
1071) option%reference_water_density, &
1072) dum1,ierr)
1073) endif
1074)
1075) ! read reaction database
1076) if (associated(realization%reaction)) then
1077) if (realization%reaction%use_full_geochemistry) then
1078) call DatabaseRead(realization%reaction,option)
1079) call BasisInit(realization%reaction,option)
1080) else
1081) ! turn off activity coefficients since the database has not been read
1082) realization%reaction%act_coef_update_frequency = ACT_COEF_FREQUENCY_OFF
1083) allocate(realization%reaction%primary_species_print(option%ntrandof))
1084) realization%reaction%primary_species_print = PETSC_TRUE
1085) endif
1086) endif
1087)
1088) ! SK 09/30/13, Added to check if Mphase is called with OS
1089) if (option%transport%reactive_transport_coupling == OPERATOR_SPLIT .and. &
1090) option%iflowmode == MPH_MODE) then
1091) option%io_buffer = 'Operator split not implemented with MPHASE. ' // &
1092) 'Switching to Global Implicit.'
1093) !geh: We should force the user to switch without automatically switching
1094) ! call printWrnMsg(option)
1095) call printErrMsg(option)
1096) option%transport%reactive_transport_coupling = GLOBAL_IMPLICIT
1097) endif
1098)
1099) ! create grid and allocate vectors
1100) call RealizationCreateDiscretization(realization)
1101)
1102) ! read any regions provided in external files
1103) call InitCommonReadRegionFiles(realization)
1104) ! clip regions and set up boundary connectivity, distance
1105) call RealizationLocalizeRegions(realization)
1106) call RealizationPassPtrsToPatches(realization)
1107) call RealizationProcessDatasets(realization)
1108) if (realization%output_option%mass_balance_region_flag) then
1109) call PatchGetCompMassInRegionAssign(realization%patch%region_list, &
1110) realization%output_option%mass_balance_region_list,option)
1111) endif
1112) ! link conditions with regions through couplers and generate connectivity
1113) call RealProcessMatPropAndSatFunc(realization)
1114) ! must process conditions before couplers in order to determine dataset types
1115) call RealizationProcessConditions(realization)
1116) call RealizationProcessCouplers(realization)
1117) call SubsurfSandboxesSetup(realization)
1118) call RealProcessFluidProperties(realization)
1119) call SubsurfInitMaterialProperties(realization)
1120) ! assignVolumesToMaterialAuxVars() must be called after
1121) ! RealizInitMaterialProperties() where the Material object is created
1122) call SubsurfAssignVolsToMatAuxVars(realization)
1123) call RealizationInitAllCouplerAuxVars(realization)
1124) if (option%ntrandof > 0) then
1125) call printMsg(option," Setting up TRAN Realization ")
1126) call RealizationInitConstraints(realization)
1127) call printMsg(option," Finished setting up TRAN Realization ")
1128) endif
1129) call RealizationPrintCouplers(realization)
1130) if (.not.option%steady_state) then
1131) ! add waypoints associated with boundary conditions, source/sinks etc. to list
1132) call RealizationAddWaypointsToList(realization, &
1133) simulation%waypoint_list_subsurface)
1134) ! fill in holes in waypoint data
1135) endif
1136) call PetscLogEventEnd(logging%event_setup,ierr);CHKERRQ(ierr)
1137)
1138) #ifdef OS_STATISTICS
1139) call RealizationPrintGridStatistics(realization)
1140) #endif
1141)
1142) #if defined(PETSC_HAVE_HDF5)
1143) #if !defined(HDF5_BROADCAST)
1144) call printMsg(option,"Default HDF5 method is used in Initialization")
1145) #else
1146) call printMsg(option,"Glenn's HDF5 broadcast method is used in Initialization")
1147) #endif
1148) #endif
1149)
1150) end subroutine SubsurfaceSetupRealization
1151)
1152) ! ************************************************************************** !
1153)
1154) subroutine SubsurfaceJumpStart(simulation)
1155) !
1156) ! Author: Glenn Hammond
1157) ! Date: 06/11/13
1158) !
1159)
1160) use Realization_Subsurface_class
1161) use Option_module
1162) use Timestepper_Base_class
1163) use Timestepper_BE_class
1164) use Output_Aux_module
1165) use Output_module, only : Output, OutputInit, OutputPrintCouplers
1166) use Condition_Control_module
1167) use Reactive_Transport_module, only : RTJumpStartKineticSorption
1168)
1169) implicit none
1170)
1171) type(simulation_subsurface_type) :: simulation
1172)
1173) class(realization_subsurface_type), pointer :: realization
1174) class(timestepper_base_type), pointer :: master_timestepper
1175) class(timestepper_BE_type), pointer :: flow_timestepper
1176) class(timestepper_BE_type), pointer :: tran_timestepper
1177) type(option_type), pointer :: option
1178) type(output_option_type), pointer :: output_option
1179)
1180) character(len=MAXSTRINGLENGTH) :: string
1181) PetscBool :: snapshot_plot_flag, observation_plot_flag
1182) PetscBool :: massbal_plot_flag
1183) PetscBool :: activity_coefs_read
1184) PetscBool :: flow_read
1185) PetscBool :: transport_read
1186) PetscBool :: failure
1187) PetscErrorCode :: ierr
1188)
1189) realization => simulation%realization
1190)
1191) if (associated(simulation%flow_process_model_coupler)) then
1192) select type(ts => simulation%flow_process_model_coupler%timestepper)
1193) class is(timestepper_BE_type)
1194) flow_timestepper => ts
1195) end select
1196) else
1197) nullify(flow_timestepper)
1198) endif
1199) if (associated(simulation%rt_process_model_coupler)) then
1200) select type(ts => simulation%rt_process_model_coupler%timestepper)
1201) class is(timestepper_BE_type)
1202) tran_timestepper => ts
1203) end select
1204) else
1205) nullify(tran_timestepper)
1206) endif
1207) nullify(master_timestepper)
1208)
1209) option => realization%option
1210)
1211) call PetscOptionsHasName(PETSC_NULL_OBJECT, &
1212) PETSC_NULL_CHARACTER, "-vecload_block_size", &
1213) failure, ierr);CHKERRQ(ierr)
1214)
1215) #if 0
1216) if (option%steady_state) then
1217) option%io_buffer = 'Running in steady-state not yet supported in &
1218) &refactored code.'
1219) call printErrMsg(option)
1220) #if 0
1221) call StepperRunSteadyState(realization,flow_timestepper,tran_timestepper)
1222) #endif
1223) ! do not want to run through time stepper
1224) option%status = DONE
1225) return
1226) endif
1227) #endif
1228)
1229) if (associated(flow_timestepper)) then
1230) master_timestepper => flow_timestepper
1231) else
1232) master_timestepper => tran_timestepper
1233) endif
1234)
1235) snapshot_plot_flag = PETSC_FALSE
1236) observation_plot_flag = PETSC_FALSE
1237) massbal_plot_flag = PETSC_FALSE
1238) activity_coefs_read = PETSC_FALSE
1239) flow_read = PETSC_FALSE
1240) transport_read = PETSC_FALSE
1241) failure = PETSC_FALSE
1242)
1243) !geh: now performed in PMRTInitializeRun()
1244) ! if (associated(simulation%rt_process_model_coupler)) then
1245) ! call simulation%rt_process_model_coupler%UpdateSolution()
1246) ! endif
1247)
1248) if (option%transport%jumpstart_kinetic_sorption .and. &
1249) option%time < 1.d-40) then
1250) ! only user jumpstart for a restarted simulation
1251) if (.not. option%restart_flag) then
1252) option%io_buffer = 'Only use JUMPSTART_KINETIC_SORPTION on a ' // &
1253) 'restarted simulation. ReactionEquilibrateConstraint() will ' // &
1254) 'appropriately set sorbed initial concentrations for a normal ' // &
1255) '(non-restarted) simulation.'
1256) call printErrMsg(option)
1257) endif
1258) call RTJumpStartKineticSorption(realization)
1259) endif
1260)
1261) end subroutine SubsurfaceJumpStart
1262)
1263) ! ************************************************************************** !
1264)
1265) subroutine SubsurfaceReadRequiredCards(simulation)
1266) !
1267) ! Reads required cards from input file
1268) !
1269) ! Author: Glenn Hammond
1270) ! Date: 10/23/07, refactored 08/20/14, refactored 12/10/14
1271) !
1272)
1273) use Option_module
1274) use Discretization_module
1275) use Grid_module
1276) use Input_Aux_module
1277) use String_module
1278) use Patch_module
1279) use Realization_Subsurface_class
1280) use HDF5_Aux_module
1281)
1282) use Simulation_Subsurface_class
1283) use General_module
1284) use Reaction_module
1285) use Reaction_Aux_module
1286) use Init_Common_module
1287)
1288) implicit none
1289)
1290) class(simulation_subsurface_type) :: simulation
1291)
1292) character(len=MAXSTRINGLENGTH) :: string
1293) character(len=MAXWORDLENGTH) :: word
1294) character(len=MAXWORDLENGTH) :: card
1295) type(patch_type), pointer :: patch, patch2
1296) type(grid_type), pointer :: grid
1297) class(realization_subsurface_type), pointer :: realization
1298) type(discretization_type), pointer :: discretization
1299) type(option_type), pointer :: option
1300) type(input_type), pointer :: input
1301)
1302) realization => simulation%realization
1303) patch => realization%patch
1304) option => realization%option
1305) discretization => realization%discretization
1306)
1307) input => realization%input
1308)
1309) ! Read in select required cards
1310) !.........................................................................
1311)
1312) ! GRID information - GRID is a required card for every simulation
1313) string = "GRID"
1314) call InputFindStringInFile(input,option,string)
1315) call InputFindStringErrorMsg(input,option,string)
1316)
1317) call DiscretizationReadRequiredCards(discretization,input,option)
1318)
1319) select case(discretization%itype)
1320) case(STRUCTURED_GRID,UNSTRUCTURED_GRID)
1321) patch => PatchCreate()
1322) patch%grid => discretization%grid
1323) if (.not.associated(realization%patch_list)) then
1324) realization%patch_list => PatchCreateList()
1325) endif
1326) call PatchAddToList(patch,realization%patch_list)
1327) realization%patch => patch
1328) end select
1329)
1330) ! optional required cards - yes, an oxymoron, but we need to know if
1331) ! these exist before we can go any further.
1332) rewind(input%fid)
1333) do
1334) call InputReadPflotranString(input,option)
1335) if (InputError(input)) exit
1336)
1337) call InputReadWord(input,option,word,PETSC_FALSE)
1338) call StringToUpper(word)
1339) card = trim(word)
1340)
1341) select case(trim(card))
1342)
1343) !....................
1344) case('DBASE_FILENAME')
1345) call InputReadWord(input,option,word,PETSC_FALSE)
1346) call InputErrorMsg(input,option,'filename','DBASE_FILENAME')
1347) if (index(word,'.h5') > 0) then
1348) #if defined(PETSC_HAVE_HDF5)
1349) call HDF5ReadDbase(word,option)
1350) #endif
1351) else
1352) call InputReadASCIIDbase(word,option)
1353) endif
1354)
1355) !....................
1356) #if defined(SCORPIO)
1357) case('HDF5_WRITE_GROUP_SIZE')
1358) call InputReadInt(input,option,option%hdf5_write_group_size)
1359) call InputErrorMsg(input,option,'HDF5_WRITE_GROUP_SIZE','Group size')
1360) call InputSkipToEnd(input,option,'HDF5_WRITE_GROUP_SIZE')
1361)
1362) case('HDF5_READ_GROUP_SIZE')
1363) call InputReadInt(input,option,option%hdf5_read_group_size)
1364) call InputErrorMsg(input,option,'HDF5_READ_GROUP_SIZE','Group size')
1365) #endif
1366)
1367) !....................
1368) case('PROC')
1369) ! processor decomposition
1370) if (realization%discretization%itype == STRUCTURED_GRID) then
1371) grid => realization%patch%grid
1372) ! strip card from front of string
1373) call InputReadInt(input,option,grid%structured_grid%npx)
1374) call InputDefaultMsg(input,option,'npx')
1375) call InputReadInt(input,option,grid%structured_grid%npy)
1376) call InputDefaultMsg(input,option,'npy')
1377) call InputReadInt(input,option,grid%structured_grid%npz)
1378) call InputDefaultMsg(input,option,'npz')
1379)
1380) if (option%myrank == option%io_rank .and. &
1381) option%print_to_screen) then
1382) option%io_buffer = ' Processor Decomposition:'
1383) call printMsg(option)
1384) write(option%io_buffer,'(" npx = ",3x,i4)') &
1385) grid%structured_grid%npx
1386) call printMsg(option)
1387) write(option%io_buffer,'(" npy = ",3x,i4)') &
1388) grid%structured_grid%npy
1389) call printMsg(option)
1390) write(option%io_buffer,'(" npz = ",3x,i4)') &
1391) grid%structured_grid%npz
1392) call printMsg(option)
1393) endif
1394)
1395) if (option%mycommsize /= grid%structured_grid%npx * &
1396) grid%structured_grid%npy * &
1397) grid%structured_grid%npz) then
1398) write(option%io_buffer,*) 'Incorrect number of processors &
1399) &specified: ',grid%structured_grid%npx*grid%structured_grid%npy* &
1400) grid%structured_grid%npz,' commsize = ',option%mycommsize
1401) call printErrMsg(option)
1402) endif
1403) endif
1404)
1405) !....................
1406) case('CHEMISTRY')
1407) if (.not.associated(simulation%rt_process_model_coupler)) then
1408) option%io_buffer = 'CHEMISTRY card included when no ' // &
1409) 'SUBSURFACE_TRANSPORT process model included in SIMULATION block.'
1410) call printErrMsg(option)
1411) endif
1412) !geh: for some reason, we need this with CHEMISTRY read for
1413) ! multicontinuum
1414) ! option%use_mc = PETSC_TRUE
1415) call ReactionInit(realization%reaction,input,option)
1416) end select
1417) enddo
1418)
1419) #if defined(SCORPIO)
1420) call InitCommonCreateIOGroups(option)
1421) #endif
1422)
1423) end subroutine SubsurfaceReadRequiredCards
1424)
1425) ! ************************************************************************** !
1426)
1427) subroutine SubsurfaceReadInput(simulation)
1428) !
1429) ! Reads pflow input file
1430) !
1431) ! Author: Glenn Hammond
1432) ! Date: 10/23/07
1433) !
1434)
1435) use Option_module
1436) use Field_module
1437) use Grid_module
1438) use Grid_Unstructured_Aux_module
1439) use Grid_Structured_module
1440) use Solver_module
1441) use Material_module
1442) use Saturation_Function_module
1443) use Characteristic_Curves_module
1444) use Dataset_Base_class
1445) use Dataset_module
1446) use Dataset_Common_HDF5_class
1447) use Fluid_module
1448) use Realization_Subsurface_class
1449) use Realization_Base_class
1450) use Region_module
1451) use Condition_module
1452) use Transport_Constraint_module
1453) use Coupler_module
1454) use Strata_module
1455) use Observation_module
1456) use Integral_Flux_module
1457) use Waypoint_module
1458) use Debug_module
1459) use Patch_module
1460) use Reaction_module
1461) use Reaction_Aux_module
1462) use Discretization_module
1463) use Input_Aux_module
1464) use String_module
1465) use Units_module
1466) use Uniform_Velocity_module
1467) use Reaction_Mineral_module
1468) use Regression_module
1469) use Output_Aux_module
1470) use Output_module
1471) use Output_Tecplot_module
1472) use Data_Mediator_Dataset_class
1473) use EOS_module
1474) use EOS_Water_module
1475) use SrcSink_Sandbox_module
1476) use Klinkenberg_module
1477) use WIPP_module
1478) use Utility_module
1479) use Checkpoint_module
1480)
1481) use Simulation_Subsurface_class
1482) use PMC_Subsurface_class
1483) use Timestepper_BE_class
1484) use Timestepper_Steady_class
1485)
1486) #ifdef SOLID_SOLUTION
1487) use Reaction_Solid_Solution_module, only : SolidSolutionReadFromInputFile
1488) #endif
1489)
1490) implicit none
1491)
1492) class(simulation_subsurface_type) :: simulation
1493)
1494) PetscErrorCode :: ierr
1495) character(len=MAXWORDLENGTH) :: word
1496) character(len=MAXWORDLENGTH) :: card
1497) character(len=MAXSTRINGLENGTH) :: string, temp_string
1498) character(len=MAXWORDLENGTH) :: internal_units
1499)
1500) character(len=1) :: backslash
1501) PetscReal :: temp_real, temp_real2
1502) PetscReal, pointer :: temp_real_array(:)
1503) PetscInt :: temp_int
1504) PetscInt :: id
1505)
1506) PetscBool :: vel_cent
1507) PetscBool :: vel_face
1508) PetscBool :: fluxes
1509) PetscBool :: mass_flowrate
1510) PetscBool :: energy_flowrate
1511) PetscBool :: aveg_mass_flowrate
1512) PetscBool :: aveg_energy_flowrate
1513)
1514) type(region_type), pointer :: region
1515) type(flow_condition_type), pointer :: flow_condition
1516) type(tran_condition_type), pointer :: tran_condition
1517) type(tran_constraint_type), pointer :: tran_constraint
1518) type(tran_constraint_type), pointer :: sec_tran_constraint
1519) type(coupler_type), pointer :: coupler
1520) type(strata_type), pointer :: strata
1521) type(observation_type), pointer :: observation
1522) type(integral_flux_type), pointer :: integral_flux
1523)
1524) type(waypoint_type), pointer :: waypoint
1525)
1526) type(material_property_type), pointer :: material_property
1527) type(fluid_property_type), pointer :: fluid_property
1528) type(saturation_function_type), pointer :: saturation_function
1529) class(characteristic_curves_type), pointer :: characteristic_curves
1530)
1531) class(realization_subsurface_type), pointer :: realization
1532) type(grid_type), pointer :: grid
1533) type(option_type), pointer :: option
1534) type(field_type), pointer :: field
1535) type(patch_type), pointer :: patch
1536) type(reaction_type), pointer :: reaction
1537) type(output_option_type), pointer :: output_option
1538) type(uniform_velocity_dataset_type), pointer :: uniform_velocity_dataset
1539) class(dataset_base_type), pointer :: dataset
1540) class(data_mediator_dataset_type), pointer :: flow_data_mediator
1541) class(data_mediator_dataset_type), pointer :: rt_data_mediator
1542) type(waypoint_list_type), pointer :: waypoint_list
1543) type(input_type), pointer :: input, input_parent
1544)
1545) PetscReal :: dt_init
1546) PetscReal :: dt_min
1547) PetscReal :: units_conversion
1548)
1549) class(timestepper_BE_type), pointer :: flow_timestepper
1550) class(timestepper_BE_type), pointer :: tran_timestepper
1551)
1552) internal_units = 'not_assigned'
1553)
1554) realization => simulation%realization
1555) output_option => simulation%output_option
1556) waypoint_list => simulation%waypoint_list_subsurface
1557) patch => realization%patch
1558)
1559) if (associated(patch)) grid => patch%grid
1560)
1561) option => realization%option
1562) field => realization%field
1563) reaction => realization%reaction
1564) input => realization%input
1565)
1566) flow_timestepper => TimestepperBECreate()
1567) flow_timestepper%solver%itype = FLOW_CLASS
1568) tran_timestepper => TimestepperBECreate()
1569) tran_timestepper%solver%itype = TRANSPORT_CLASS
1570)
1571) backslash = achar(92) ! 92 = "\" Some compilers choke on \" thinking it
1572) ! is a double quote as in c/c++
1573)
1574) rewind(input%fid)
1575) string = 'SUBSURFACE'
1576) call InputFindStringInFile(input,option,string)
1577) call InputFindStringErrorMsg(input,option,string)
1578)
1579) do
1580) call InputReadPflotranString(input,option)
1581) if (InputError(input)) exit
1582)
1583) call InputReadWord(input,option,word,PETSC_FALSE)
1584) call StringToUpper(word)
1585) card = trim(word)
1586)
1587) option%io_buffer = 'pflotran card:: ' // trim(card)
1588) call printMsg(option)
1589)
1590) select case(trim(card))
1591)
1592) !....................
1593) case ('GRID')
1594) call DiscretizationRead(realization%discretization,input,option)
1595)
1596) !....................
1597) case ('CHEMISTRY')
1598) call ReactionReadPass2(reaction,input,option)
1599)
1600) !....................
1601) case ('NONUNIFORM_VELOCITY')
1602) call InputReadNChars(input,option, &
1603) realization%nonuniform_velocity_filename, &
1604) MAXSTRINGLENGTH,PETSC_TRUE)
1605) call InputErrorMsg(input,option,'filename','NONUNIFORM_VELOCITY')
1606)
1607) case ('UNIFORM_VELOCITY')
1608) uniform_velocity_dataset => UniformVelocityDatasetCreate()
1609) uniform_velocity_dataset%rank = 3
1610) uniform_velocity_dataset%interpolation_method = 1 ! 1 = STEP
1611) uniform_velocity_dataset%is_cyclic = PETSC_FALSE
1612) allocate(uniform_velocity_dataset%times(1))
1613) uniform_velocity_dataset%times = 0.d0
1614) allocate(uniform_velocity_dataset%values(3,1))
1615) uniform_velocity_dataset%values = 0.d0
1616) call InputReadDouble(input,option,uniform_velocity_dataset%values(1,1))
1617) call InputErrorMsg(input,option,'velx','UNIFORM_VELOCITY')
1618) call InputReadDouble(input,option,uniform_velocity_dataset%values(2,1))
1619) call InputErrorMsg(input,option,'vely','UNIFORM_VELOCITY')
1620) call InputReadDouble(input,option,uniform_velocity_dataset%values(3,1))
1621) call InputErrorMsg(input,option,'velz','UNIFORM_VELOCITY')
1622) ! read units, if present
1623) temp_real = 1.d0
1624) call InputReadAndConvertUnits(input,temp_real, &
1625) 'meter/sec','UNIFORM_VELOCITY,units',option)
1626) uniform_velocity_dataset%values(:,1) = &
1627) uniform_velocity_dataset%values(:,1) * temp_real
1628) call UniformVelocityDatasetVerify(option,uniform_velocity_dataset)
1629) realization%uniform_velocity_dataset => uniform_velocity_dataset
1630)
1631) case ('VELOCITY_DATASET')
1632) uniform_velocity_dataset => UniformVelocityDatasetCreate()
1633) call UniformVelocityDatasetRead(uniform_velocity_dataset,input,option)
1634) realization%uniform_velocity_dataset => uniform_velocity_dataset
1635)
1636) !....................
1637) case ('DEBUG')
1638) call DebugRead(realization%debug,input,option)
1639)
1640) !....................
1641) case ('PRINT_PRIMAL_GRID')
1642) option%print_explicit_primal_grid = PETSC_TRUE
1643)
1644) !....................
1645) case ('PRINT_DUAL_GRID')
1646) option%print_explicit_dual_grid = PETSC_TRUE
1647)
1648) !....................
1649) case ('PROC')
1650)
1651) !....................
1652) case ('REGION')
1653) region => RegionCreate()
1654) call InputReadWord(input,option,region%name,PETSC_TRUE)
1655) call InputErrorMsg(input,option,'name','REGION')
1656) call printMsg(option,region%name)
1657) call RegionRead(region,input,option)
1658) ! we don't copy regions down to patches quite yet, since we
1659) ! don't want to duplicate IO in reading the regions
1660) call RegionAddToList(region,realization%region_list)
1661) nullify(region)
1662)
1663) !....................
1664) case ('FLOW_CONDITION')
1665) flow_condition => FlowConditionCreate(option)
1666) call InputReadWord(input,option,flow_condition%name,PETSC_TRUE)
1667) call InputErrorMsg(input,option,'FLOW_CONDITION','name')
1668) call printMsg(option,flow_condition%name)
1669) if (option%iflowmode == G_MODE) then
1670) call FlowConditionGeneralRead(flow_condition,input,option)
1671) else if(option%iflowmode == TOIL_IMS_MODE) then
1672) call FlowConditionTOilImsRead(flow_condition,input,option)
1673) else
1674) call FlowConditionRead(flow_condition,input,option)
1675) endif
1676) call FlowConditionAddToList(flow_condition,realization%flow_conditions)
1677) nullify(flow_condition)
1678)
1679) !....................
1680) case ('TRANSPORT_CONDITION')
1681) if (.not.associated(reaction)) then
1682) option%io_buffer = 'TRANSPORT_CONDITIONs not supported without ' // &
1683) 'CHEMISTRY.'
1684) call printErrMsg(option)
1685) endif
1686) tran_condition => TranConditionCreate(option)
1687) call InputReadWord(input,option,tran_condition%name,PETSC_TRUE)
1688) call InputErrorMsg(input,option,'TRANSPORT_CONDITION','name')
1689) call printMsg(option,tran_condition%name)
1690) call TranConditionRead(tran_condition,realization%transport_constraints, &
1691) reaction,input,option)
1692) call TranConditionAddToList(tran_condition,realization%transport_conditions)
1693) nullify(tran_condition)
1694)
1695) !....................
1696) case('CONSTRAINT')
1697) if (.not.associated(reaction)) then
1698) option%io_buffer = 'CONSTRAINTs not supported without CHEMISTRY.'
1699) call printErrMsg(option)
1700) endif
1701) tran_constraint => TranConstraintCreate(option)
1702) call InputReadWord(input,option,tran_constraint%name,PETSC_TRUE)
1703) call InputErrorMsg(input,option,'constraint','name')
1704) call printMsg(option,tran_constraint%name)
1705) call TranConstraintRead(tran_constraint,reaction,input,option)
1706) call TranConstraintAddToList(tran_constraint,realization%transport_constraints)
1707) nullify(tran_constraint)
1708)
1709)
1710) !....................
1711) case ('BOUNDARY_CONDITION')
1712) coupler => CouplerCreate(BOUNDARY_COUPLER_TYPE)
1713) call InputReadWord(input,option,coupler%name,PETSC_TRUE)
1714) call InputDefaultMsg(input,option,'Boundary Condition name')
1715) call CouplerRead(coupler,input,option)
1716) call RealizationAddCoupler(realization,coupler)
1717) nullify(coupler)
1718)
1719) !....................
1720) case ('INITIAL_CONDITION')
1721) coupler => CouplerCreate(INITIAL_COUPLER_TYPE)
1722) call InputReadWord(input,option,coupler%name,PETSC_TRUE)
1723) call InputDefaultMsg(input,option,'Initial Condition name')
1724) call CouplerRead(coupler,input,option)
1725) call RealizationAddCoupler(realization,coupler)
1726) nullify(coupler)
1727)
1728) !....................
1729) case ('SOURCE_SINK')
1730) coupler => CouplerCreate(SRC_SINK_COUPLER_TYPE)
1731) call InputReadWord(input,option,coupler%name,PETSC_TRUE)
1732) call InputDefaultMsg(input,option,'Source Sink name')
1733) call CouplerRead(coupler,input,option)
1734) call RealizationAddCoupler(realization,coupler)
1735) nullify(coupler)
1736)
1737) !....................
1738) case ('SOURCE_SINK_SANDBOX')
1739) call SSSandboxInit(option)
1740) call SSSandboxRead(input,option)
1741)
1742) !....................
1743) case ('FLOW_MASS_TRANSFER')
1744) flow_data_mediator => DataMediatorDatasetCreate()
1745) call InputReadWord(input,option,flow_data_mediator%name,PETSC_TRUE)
1746) call InputDefaultMsg(input,option,'Flow Mass Transfer name')
1747) call DataMediatorDatasetRead(flow_data_mediator,input,option)
1748) call flow_data_mediator%AddToList(realization%flow_data_mediator_list)
1749) nullify(flow_data_mediator)
1750)
1751) !....................
1752) case ('RT_MASS_TRANSFER')
1753) rt_data_mediator => DataMediatorDatasetCreate()
1754) call InputReadWord(input,option,rt_data_mediator%name,PETSC_TRUE)
1755) call InputDefaultMsg(input,option,'RT Mass Transfer name')
1756) call DataMediatorDatasetRead(rt_data_mediator,input,option)
1757) call rt_data_mediator%AddToList(realization%tran_data_mediator_list)
1758) nullify(rt_data_mediator)
1759)
1760) !....................
1761) case ('STRATIGRAPHY','STRATA')
1762) strata => StrataCreate()
1763) call StrataRead(strata,input,option)
1764) call RealizationAddStrata(realization,strata)
1765) nullify(strata)
1766)
1767) !.....................
1768) case ('DATASET')
1769) nullify(dataset)
1770) call DatasetRead(input,dataset,option)
1771) call DatasetBaseAddToList(dataset,realization%datasets)
1772) nullify(dataset)
1773)
1774) !....................
1775)
1776) case('REFERENCE_PRESSURE')
1777) call InputReadStringErrorMsg(input,option,card)
1778) call InputReadDouble(input,option,option%reference_pressure)
1779) call InputDefaultMsg(input,option,'Reference Pressure')
1780)
1781) !....................
1782)
1783) case('REFERENCE_DENSITY')
1784) call InputReadStringErrorMsg(input,option,card)
1785) call InputReadDouble(input,option,option%reference_water_density)
1786) call InputDefaultMsg(input,option,'Reference Density')
1787)
1788) !....................
1789)
1790) case('MINIMUM_HYDROSTATIC_PRESSURE')
1791) call InputReadStringErrorMsg(input,option,card)
1792) call InputReadDouble(input,option,option%minimum_hydrostatic_pressure)
1793) call InputDefaultMsg(input,option,'Minimum Hydrostatic Pressure')
1794)
1795) !......................
1796)
1797) case('REFERENCE_TEMPERATURE')
1798) call InputReadStringErrorMsg(input,option,card)
1799) call InputReadDouble(input,option,option%reference_temperature)
1800) call InputDefaultMsg(input,option,'Reference Temperature')
1801)
1802) !......................
1803)
1804) case('REFERENCE_POROSITY')
1805) call InputReadStringErrorMsg(input,option,card)
1806) call InputReadDouble(input,option,option%reference_porosity)
1807) call InputDefaultMsg(input,option,'Reference Porosity')
1808)
1809) !......................
1810)
1811) case('REFERENCE_SATURATION')
1812) call InputReadStringErrorMsg(input,option,card)
1813) call InputReadDouble(input,option,option%reference_saturation)
1814) call InputDefaultMsg(input,option,'Reference Saturation')
1815)
1816) !......................
1817)
1818) case('NONISOTHERMAL')
1819) option%use_isothermal = PETSC_FALSE
1820)
1821) !......................
1822)
1823) case('ISOTHERMAL')
1824) option%use_isothermal = PETSC_TRUE
1825)
1826) !......................
1827)
1828) case('UPDATE_FLOW_PERMEABILITY')
1829) option%update_flow_perm = PETSC_TRUE
1830)
1831) !......................
1832)
1833) case('DFN')
1834) grid%unstructured_grid%grid_type = TWO_DIM_GRID
1835)
1836) !......................
1837)
1838) case("MULTIPLE_CONTINUUM")
1839) option%use_mc = PETSC_TRUE
1840)
1841) !......................
1842)
1843) case('SECONDARY_CONTINUUM_SOLVER')
1844) if (.not.option%use_mc) then
1845) option%io_buffer = 'SECONDARY_CONTINUUM_SOLVER can only be used ' // &
1846) 'with MULTIPLE_CONTINUUM keyword.'
1847) call printErrMsg(option)
1848) endif
1849) call InputReadWord(input,option,word,PETSC_FALSE)
1850) call StringToUpper(word)
1851) select case(word)
1852) case('KEARST')
1853) option%secondary_continuum_solver = 1
1854) case('HINDMARSH')
1855) option%secondary_continuum_solver = 2
1856) case('THOMAS')
1857) option%secondary_continuum_solver = 3
1858) case default
1859) option%io_buffer = 'SECONDARY_CONTINUUM_SOLVER can be only ' // &
1860) 'HINDMARSH or KEARST. For single component'// &
1861) 'chemistry THOMAS can be used.'
1862) call printErrMsg(option)
1863) end select
1864) !....................
1865)
1866) case('SECONDARY_CONSTRAINT')
1867) if (.not.option%use_mc) then
1868) option%io_buffer = 'SECONDARY_CONSTRAINT can only be used with ' // &
1869) 'MULTIPLE_CONTINUUM keyword.'
1870) call printErrMsg(option)
1871) endif
1872) if (.not.associated(reaction)) then
1873) option%io_buffer = 'SECONDARY_CONSTRAINT not supported without' // &
1874) 'CHEMISTRY.'
1875) call printErrMsg(option)
1876) endif
1877) sec_tran_constraint => TranConstraintCreate(option)
1878) call InputReadWord(input,option,sec_tran_constraint%name,PETSC_TRUE)
1879) call InputErrorMsg(input,option,'secondary constraint','name')
1880) call printMsg(option,sec_tran_constraint%name)
1881) call TranConstraintRead(sec_tran_constraint,reaction,input,option)
1882) realization%sec_transport_constraint => sec_tran_constraint
1883) nullify(sec_tran_constraint)
1884)
1885) !......................
1886)
1887) case('BRIN','BRINE')
1888) call InputReadStringErrorMsg(input,option,card)
1889) call InputReadDouble(input,option,option%m_nacl)
1890) call InputDefaultMsg(input,option,'NaCl Concentration')
1891)
1892) call InputReadWord(input,option,word,PETSC_FALSE)
1893) call StringToUpper(word)
1894) select case(word(1:len_trim(word)))
1895) case('MOLAL')
1896) case('MASS')
1897) option%m_nacl = option%m_nacl /FMWNACL/(1.D0-option%m_nacl)
1898) case('MOLE')
1899) option%m_nacl = option%m_nacl /FMWH2O/(1.D0-option%m_nacl)
1900) case default
1901) print *, 'Wrong unit: ', word(1:len_trim(word))
1902) stop
1903) end select
1904) if (OptionPrintToScreen(option)) print *, option%m_nacl
1905) !......................
1906)
1907) case ('RESTART')
1908) option%io_buffer = 'The RESTART card within SUBSURFACE block has &
1909) &been deprecated.'
1910) call printErrMsg(option)
1911) ! option%restart_flag = PETSC_TRUE
1912) !call InputReadNChars(input,option,option%restart_filename,MAXSTRINGLENGTH, &
1913) ! PETSC_TRUE)
1914) !call InputErrorMsg(input,option,'RESTART','Restart file name')
1915) !call InputReadDouble(input,option,option%restart_time)
1916) !if (input%ierr == 0) then
1917) ! call InputReadWord(input,option,word,PETSC_TRUE)
1918) ! if (input%ierr == 0) then
1919) ! internal_units = 'sec'
1920) ! option%restart_time = option%restart_time* &
1921) ! UnitsConvertToInternal(word,internal_units,option)
1922) ! else
1923) ! call InputDefaultMsg(input,option,'RESTART, time units')
1924) ! endif
1925) !endif
1926)
1927) !......................
1928)
1929) case ('CHECKPOINT')
1930) option%io_buffer = 'The CHECKPOINT card within SUBSURFACE block must &
1931) &be moved to the SIMULATION block.'
1932) call printErrMsg(option)
1933) ! call CheckpointRead(input,option,realization%checkpoint_option, &
1934) ! realization%waypoint_list)
1935)
1936) !......................
1937)
1938) case ('NUMERICAL_JACOBIAN_FLOW')
1939) option%io_buffer = 'The NUMERICAL_JACOBIAN_FLOW card within &
1940) &SUBSURFACE block must be listed under the SIMULATION/&
1941) &PROCESS_MODELS/SUBSURFACE_FLOW/OPTIONS block as NUMERICAL_JACOBIAN.'
1942) call printErrMsg(option)
1943)
1944) !......................
1945)
1946) case ('NUMERICAL_JACOBIAN_RXN')
1947) option%io_buffer = 'The NUMERICAL_JACOBIAN_FLOW card within &
1948) &SUBSURFACE block must be listed under the SIMULATION/&
1949) &PROCESS_MODELS/SUBSURFACE_TRANSPORT block as &
1950) &NUMERICAL_JACOBIAN.'
1951) call printErrMsg(option)
1952)
1953) !......................
1954)
1955) case ('NUMERICAL_JACOBIAN_MULTI_COUPLE')
1956) option%numerical_derivatives_multi_coupling = PETSC_TRUE
1957)
1958) !......................
1959)
1960) case ('COMPUTE_STATISTICS')
1961) option%compute_statistics = PETSC_TRUE
1962)
1963) !....................
1964)
1965) case ('CO2_DATABASE')
1966) call InputReadNChars(input,option,option%co2_database_filename, &
1967) MAXSTRINGLENGTH,PETSC_TRUE)
1968) call InputErrorMsg(input,option,'CO2_DATABASE','filename')
1969)
1970) !....................
1971)
1972) case ('TIMESTEPPER')
1973) call InputReadWord(input,option,word,PETSC_FALSE)
1974) call StringToUpper(word)
1975) select case(word)
1976) case('FLOW')
1977) call flow_timestepper%ReadInput(input,option)
1978) case('TRAN','TRANSPORT')
1979) call tran_timestepper%ReadInput(input,option)
1980) case default
1981) option%io_buffer = 'TIMESTEPPER must specify FLOW or TRANSPORT.'
1982) call printErrMsg(option)
1983) end select
1984)
1985) !....................
1986)
1987) case ('LINEAR_SOLVER')
1988) call InputReadWord(input,option,word,PETSC_FALSE)
1989) call StringToUpper(word)
1990) select case(word)
1991) case('FLOW')
1992) call SolverReadLinear(flow_timestepper%solver,input,option)
1993) case('TRAN','TRANSPORT')
1994) call SolverReadLinear(tran_timestepper%solver,input,option)
1995) case default
1996) option%io_buffer = 'LINEAR_SOLVER must specify FLOW or TRANSPORT.'
1997) call printErrMsg(option)
1998) end select
1999)
2000) !....................
2001)
2002) case ('NEWTON_SOLVER')
2003) call InputReadWord(input,option,word,PETSC_FALSE)
2004) call StringToUpper(word)
2005) select case(word)
2006) case('FLOW')
2007) call SolverReadNewton(flow_timestepper%solver,input,option)
2008) case('TRAN','TRANSPORT')
2009) call SolverReadNewton(tran_timestepper%solver,input,option)
2010) case default
2011) option%io_buffer = 'NEWTON_SOLVER must specify FLOW or TRANSPORT.'
2012) call printErrMsg(option)
2013) end select
2014) !....................
2015)
2016) case ('FLUID_PROPERTY')
2017)
2018) fluid_property => FluidPropertyCreate()
2019) call FluidPropertyRead(fluid_property,input,option)
2020) call FluidPropertyAddToList(fluid_property,realization%fluid_properties)
2021) nullify(fluid_property)
2022)
2023) !....................
2024)
2025) case ('EOS')
2026) call EOSRead(input,option)
2027)
2028) !....................
2029)
2030) case ('SATURATION_FUNCTION')
2031) if (option%iflowmode == RICHARDS_MODE .or. &
2032) option%iflowmode == TOIL_IMS_MODE .or. &
2033) option%iflowmode == G_MODE) then
2034) option%io_buffer = &
2035) 'Must compile with legacy_saturation_function=1 ' //&
2036) 'to use the SATURATION_FUNCTION keyword. Otherwise, use ' // &
2037) 'CHARACTERISTIC_CURVES.'
2038) call printErrMsg(option)
2039) endif
2040) saturation_function => SaturationFunctionCreate(option)
2041) call InputReadWord(input,option,saturation_function%name,PETSC_TRUE)
2042) call InputErrorMsg(input,option,'name','SATURATION_FUNCTION')
2043) call SaturationFunctionRead(saturation_function,input,option)
2044) call SatFunctionComputePolynomial(option,saturation_function)
2045) call PermFunctionComputePolynomial(option,saturation_function)
2046) call SaturationFunctionAddToList(saturation_function, &
2047) realization%saturation_functions)
2048) nullify(saturation_function)
2049)
2050) !....................
2051)
2052) case ('CHARACTERISTIC_CURVES')
2053)
2054) if (.not.(option%iflowmode == NULL_MODE .or. &
2055) option%iflowmode == RICHARDS_MODE .or. &
2056) option%iflowmode == TOIL_IMS_MODE .or. &
2057) option%iflowmode == G_MODE)) then
2058) option%io_buffer = 'CHARACTERISTIC_CURVES not supported in flow ' // &
2059) 'modes other than RICHARDS, TOIL_IMS, or GENERAL. Use ' // &
2060) 'SATURATION_FUNCTION.'
2061) call printErrMsg(option)
2062) endif
2063) characteristic_curves => CharacteristicCurvesCreate()
2064) call InputReadWord(input,option,characteristic_curves%name,PETSC_TRUE)
2065) call InputErrorMsg(input,option,'name','CHARACTERISTIC_CURVES')
2066) call CharacteristicCurvesRead(characteristic_curves,input,option)
2067) ! call SatFunctionComputePolynomial(option,saturation_function)
2068) ! call PermFunctionComputePolynomial(option,saturation_function)
2069) call CharacteristicCurvesAddToList(characteristic_curves, &
2070) realization%characteristic_curves)
2071) nullify(characteristic_curves)
2072)
2073) !....................
2074)
2075) case ('MATERIAL_PROPERTY')
2076)
2077) material_property => MaterialPropertyCreate()
2078) call InputReadWord(input,option,material_property%name,PETSC_TRUE)
2079) call InputErrorMsg(input,option,'name','MATERIAL_PROPERTY')
2080) call MaterialPropertyRead(material_property,input,option)
2081) call MaterialPropertyAddToList(material_property, &
2082) realization%material_properties)
2083) nullify(material_property)
2084)
2085) !....................
2086)
2087) case ('USE_TOUCH_OPTIONS')
2088) option%use_touch_options = PETSC_TRUE
2089)
2090) case ('MPI_IO')
2091) ! call PetscOptionsInsertString(PETSC_NULL_OBJECT, &
2092) ! '-viewer_binary_mpiio')
2093)
2094) case ('HANDSHAKE_IO')
2095) call InputReadInt(input,option,option%io_handshake_buffer_size)
2096) call InputErrorMsg(input,option,'io_handshake_buffer_size','HANDSHAKE_IO')
2097)
2098) case ('OVERWRITE_RESTART_TRANSPORT')
2099) option%overwrite_restart_transport = PETSC_TRUE
2100)
2101) case ('OVERWRITE_RESTART_FLOW_PARAMS')
2102) option%overwrite_restart_flow = PETSC_TRUE
2103)
2104) case ('INITIALIZE_FLOW_FROM_FILE')
2105) call InputReadNChars(input,option,option%initialize_flow_filename, &
2106) MAXSTRINGLENGTH,PETSC_TRUE)
2107) call InputErrorMsg(input,option,'filename','INITIALIZE_FLOW_FROM_FILE')
2108)
2109) case ('INITIALIZE_TRANSPORT_FROM_FILE')
2110) call InputReadNChars(input,option,option%initialize_transport_filename, &
2111) MAXSTRINGLENGTH,PETSC_TRUE)
2112) call InputErrorMsg(input,option,'filename','INITIALIZE_TRANSPORT_FROM_FILE')
2113)
2114) case ('CENTRAL_DIFFERENCE')
2115) option%use_upwinding = PETSC_FALSE
2116)
2117) !....................
2118) case ('OBSERVATION')
2119) observation => ObservationCreate()
2120) call ObservationRead(observation,input,option)
2121) call ObservationAddToList(observation, &
2122) realization%patch%observation_list)
2123) nullify(observation)
2124)
2125) !....................
2126) case ('INTEGRAL_FLUX')
2127) integral_flux => IntegralFluxCreate()
2128) call InputReadWord(input,option,integral_flux%name,PETSC_TRUE)
2129) call InputDefaultMsg(input,option,'Integral Flux name')
2130) call IntegralFluxRead(integral_flux,input,option)
2131) call IntegralFluxAddToList(integral_flux, &
2132) realization%patch%integral_flux_list)
2133) nullify(integral_flux)
2134)
2135) !.....................
2136) case ('WALLCLOCK_STOP')
2137) option%wallclock_stop_flag = PETSC_TRUE
2138) call InputReadDouble(input,option,option%wallclock_stop_time)
2139) call InputErrorMsg(input,option,'stop time','WALLCLOCK_STOP')
2140)
2141) call InputReadWord(input,option,word,PETSC_TRUE)
2142) if (input%ierr /= 0) word = 'h'
2143) call InputDefaultMsg(input,option,'WALLCLOCK_STOP time units')
2144) internal_units = 'sec'
2145) units_conversion = UnitsConvertToInternal(word,internal_units,option)
2146) ! convert from hrs to seconds and add to start_time
2147) option%wallclock_stop_time = option%start_time + &
2148) option%wallclock_stop_time* &
2149) units_conversion
2150)
2151) !....................
2152) case ('OUTPUT')
2153) vel_cent = PETSC_FALSE
2154) vel_face = PETSC_FALSE
2155) fluxes = PETSC_FALSE
2156) mass_flowrate = PETSC_FALSE
2157) energy_flowrate = PETSC_FALSE
2158) aveg_mass_flowrate = PETSC_FALSE
2159) aveg_energy_flowrate = PETSC_FALSE
2160) do
2161) call InputReadPflotranString(input,option)
2162) call InputReadStringErrorMsg(input,option,card)
2163) if (InputCheckExit(input,option)) exit
2164) call InputReadWord(input,option,word,PETSC_TRUE)
2165) call InputErrorMsg(input,option,'keyword','OUTPUT')
2166) call StringToUpper(word)
2167) !----------------------------------------------------------------------
2168) !----- NEW INPUT FORMAT: ----------------------------------------------
2169) !----------------------------------------------------------------------
2170) select case(trim(word))
2171) case('OBSERVATION_FILE')
2172) call OutputFileRead(realization,output_option, &
2173) waypoint_list,trim(word))
2174) case('SNAPSHOT_FILE')
2175) call OutputFileRead(realization,output_option, &
2176) waypoint_list,trim(word))
2177) case('MASS_BALANCE_FILE')
2178) call OutputFileRead(realization,output_option, &
2179) waypoint_list,trim(word))
2180) case('TIME_UNITS')
2181) call InputReadWord(input,option,word,PETSC_TRUE)
2182) call InputErrorMsg(input,option,'Output Time Units','OUTPUT')
2183) output_option%tunit = trim(word)
2184) internal_units = 'sec'
2185) output_option%tconv = &
2186) UnitsConvertToInternal(word,internal_units,option)
2187) case('VARIABLES')
2188) select case (option%iflowmode)
2189) case(FLASH2_MODE,MPH_MODE)
2190) option%io_buffer = 'A variable list cannot be specified for &
2191) &the CO2 flow modes. Variables are determined internally.'
2192) call printErrMsg(option)
2193) end select
2194) call OutputVariableRead(input,option, &
2195) output_option%output_variable_list)
2196) case('AVERAGE_VARIABLES')
2197) call OutputVariableRead(input,option, &
2198) output_option%aveg_output_variable_list)
2199) case('UNFILTER_NON_STATE_VARIABLES')
2200) output_option%filter_non_state_variables = PETSC_FALSE
2201)
2202)
2203) !----------------------------------------------------------------------
2204) !----- SUPPORT FOR OLD INPUT FORMAT: ----------------------------------
2205) !----------------------------------------------------------------------
2206) case('NO_FINAL','NO_PRINT_FINAL')
2207) output_option%print_final_obs = PETSC_FALSE
2208) output_option%print_final_snap = PETSC_FALSE
2209) output_option%print_final_massbal = PETSC_FALSE
2210) case('NO_INITIAL','NO_PRINT_INITIAL')
2211) output_option%print_initial_obs = PETSC_FALSE
2212) output_option%print_initial_snap = PETSC_FALSE
2213) output_option%print_initial_massbal = PETSC_FALSE
2214) case('PROCESSOR_ID')
2215) option%io_buffer = 'PROCESSOR_ID output must now be entered &
2216) &under OUTPUT/VARIABLES card as PROCESS_ID.'
2217) call printErrMsg(option)
2218) ! output_option%print_iproc = PETSC_TRUE
2219) case('PERMEABILITY')
2220) option%io_buffer = 'PERMEABILITY output must now be entered &
2221) &under OUTPUT/VARIABLES card.'
2222) call printErrMsg(option)
2223) ! output_option%print_permeability = PETSC_TRUE
2224) case('POROSITY')
2225) option%io_buffer = 'POROSITY output must now be entered under &
2226) &OUTPUT/VARIABLES card.'
2227) call printErrMsg(option)
2228) ! output_option%print_porosity = PETSC_TRUE
2229) case('TORTUOSITY')
2230) option%io_buffer = 'TORTUOSITY output must now be entered under &
2231) &OUTPUT/VARIABLES card.'
2232) call printErrMsg(option)
2233) ! output_option%print_tortuosity = PETSC_TRUE
2234) case('VOLUME')
2235) option%io_buffer = 'VOLUME output must now be entered under &
2236) &OUTPUT/VARIABLES card.'
2237) call printErrMsg(option)
2238) ! output_option%print_volume = PETSC_TRUE
2239) case('MASS_BALANCE')
2240) option%compute_mass_balance_new = PETSC_TRUE
2241) output_option%periodic_msbl_output_ts_imod = 1
2242) call InputReadWord(input,option,word,PETSC_TRUE)
2243) call InputDefaultMsg(input,option, &
2244) 'OUTPUT,MASS_BALANCE,DETAILED')
2245) if (len_trim(word) > 0) then
2246) call StringToUpper(word)
2247) select case(trim(word))
2248) case('DETAILED')
2249) option%mass_bal_detailed = PETSC_TRUE
2250) case default
2251) call InputKeywordUnrecognized(word, &
2252) 'OUTPUT,MASS_BALANCE',option)
2253) end select
2254) endif
2255) case('PRINT_COLUMN_IDS')
2256) output_option%print_column_ids = PETSC_TRUE
2257) case('TIMES')
2258) call InputReadWord(input,option,word,PETSC_TRUE)
2259) call InputErrorMsg(input,option,'units','OUTPUT,TIMES')
2260) internal_units = 'sec'
2261) units_conversion = &
2262) UnitsConvertToInternal(word,internal_units,option)
2263) string = 'OUTPUT,TIMES'
2264) nullify(temp_real_array)
2265) call UtilityReadArray(temp_real_array,NEG_ONE_INTEGER, &
2266) string,input,option)
2267) do temp_int = 1, size(temp_real_array)
2268) waypoint => WaypointCreate()
2269) waypoint%time = temp_real_array(temp_int)*units_conversion
2270) waypoint%print_snap_output = PETSC_TRUE
2271) call WaypointInsertInList(waypoint,waypoint_list)
2272) enddo
2273) call DeallocateArray(temp_real_array)
2274) case('OUTPUT_FILE')
2275) call InputReadWord(input,option,word,PETSC_TRUE)
2276) call InputErrorMsg(input,option,'time increment', &
2277) 'OUTPUT,OUTPUT_FILE')
2278) call StringToUpper(word)
2279) select case(trim(word))
2280) case('OFF')
2281) option%print_to_file = PETSC_FALSE
2282) case('PERIODIC')
2283) call InputReadInt(input,option,output_option%output_file_imod)
2284) call InputErrorMsg(input,option,'timestep increment', &
2285) 'OUTPUT,OUTPUT_FILE,PERIODIC')
2286) case default
2287) call InputKeywordUnrecognized(word, &
2288) 'OUTPUT,OUTPUT_FILE',option)
2289) end select
2290) case('SCREEN')
2291) call InputReadWord(input,option,word,PETSC_TRUE)
2292) call InputErrorMsg(input,option,'time increment','OUTPUT,SCREEN')
2293) call StringToUpper(word)
2294) select case(trim(word))
2295) case('OFF')
2296) option%print_to_screen = PETSC_FALSE
2297) case('PERIODIC')
2298) call InputReadInt(input,option,output_option%screen_imod)
2299) call InputErrorMsg(input,option,'timestep increment', &
2300) 'OUTPUT,PERIODIC,SCREEN')
2301) case default
2302) call InputKeywordUnrecognized(word, &
2303) 'OUTPUT,SCREEN',option)
2304) end select
2305) case('PERIODIC')
2306) call InputReadWord(input,option,word,PETSC_TRUE)
2307) call InputErrorMsg(input,option,'time increment', &
2308) 'OUTPUT,PERIODIC')
2309) call StringToUpper(word)
2310) select case(trim(word))
2311) case('TIME')
2312) call InputReadDouble(input,option,temp_real)
2313) call InputErrorMsg(input,option,'time increment', &
2314) 'OUTPUT,PERIODIC,TIME')
2315) call InputReadWord(input,option,word,PETSC_TRUE)
2316) call InputErrorMsg(input,option,'time increment units', &
2317) 'OUTPUT,PERIODIC,TIME')
2318) internal_units = 'sec'
2319) units_conversion = UnitsConvertToInternal(word, &
2320) internal_units,option)
2321) output_option%periodic_snap_output_time_incr = temp_real* &
2322) units_conversion
2323) call InputReadWord(input,option,word,PETSC_TRUE)
2324) if (input%ierr == 0) then
2325) if (StringCompareIgnoreCase(word,'between')) then
2326)
2327) call InputReadDouble(input,option,temp_real)
2328) call InputErrorMsg(input,option,'start time', &
2329) 'OUTPUT,PERIODIC,TIME')
2330) call InputReadWord(input,option,word,PETSC_TRUE)
2331) call InputErrorMsg(input,option,'start time units', &
2332) 'OUTPUT,PERIODIC,TIME')
2333) internal_units = 'sec'
2334) units_conversion = UnitsConvertToInternal(word, &
2335) internal_units,option)
2336) temp_real = temp_real * units_conversion
2337) call InputReadWord(input,option,word,PETSC_TRUE)
2338) if (.not.StringCompareIgnoreCase(word,'and')) then
2339) input%ierr = 1
2340) endif
2341) call InputErrorMsg(input,option,'and', &
2342) 'OUTPUT,PERIODIC,TIME"')
2343) call InputReadDouble(input,option,temp_real2)
2344) call InputErrorMsg(input,option,'end time', &
2345) 'OUTPUT,PERIODIC,TIME')
2346) call InputReadWord(input,option,word,PETSC_TRUE)
2347) call InputErrorMsg(input,option,'end time units', &
2348) 'OUTPUT,PERIODIC,TIME')
2349) temp_real2 = temp_real2 * units_conversion
2350) do
2351) waypoint => WaypointCreate()
2352) waypoint%time = temp_real
2353) waypoint%print_snap_output = PETSC_TRUE
2354) call WaypointInsertInList(waypoint,waypoint_list)
2355) temp_real = temp_real + &
2356) output_option%periodic_snap_output_time_incr
2357) if (temp_real > temp_real2) exit
2358) enddo
2359) output_option%periodic_snap_output_time_incr = 0.d0
2360) else
2361) input%ierr = 1
2362) call InputErrorMsg(input,option,'between', &
2363) 'OUTPUT,PERIODIC,TIME')
2364) endif
2365) endif
2366) case('TIMESTEP')
2367) call InputReadInt(input,option, &
2368) output_option%periodic_snap_output_ts_imod)
2369) call InputErrorMsg(input,option,'timestep increment', &
2370) 'OUTPUT,PERIODIC,TIMESTEP')
2371) case default
2372) call InputKeywordUnrecognized(word, &
2373) 'OUTPUT,PERIODIC',option)
2374) end select
2375) case('OBSERVATION_TIMES')
2376) output_option%print_observation = PETSC_TRUE
2377) call InputReadWord(input,option,word,PETSC_TRUE)
2378) call InputErrorMsg(input,option,'time units', &
2379) 'OUTPUT,OBSERVATION_TIMES')
2380) internal_units = 'sec'
2381) units_conversion = &
2382) UnitsConvertToInternal(word,internal_units,option)
2383) string = 'OBSERVATION_TIMES,TIMES'
2384) nullify(temp_real_array)
2385) call UtilityReadArray(temp_real_array,NEG_ONE_INTEGER, &
2386) string,input,option)
2387) do temp_int = 1, size(temp_real_array)
2388) waypoint => WaypointCreate()
2389) waypoint%time = temp_real_array(temp_int)*units_conversion
2390) waypoint%print_obs_output = PETSC_TRUE
2391) call WaypointInsertInList(waypoint,waypoint_list)
2392) waypoint => WaypointCreate()
2393) waypoint%time = temp_real_array(temp_int)*units_conversion
2394) waypoint%print_msbl_output = PETSC_TRUE
2395) call WaypointInsertInList(waypoint,waypoint_list)
2396) enddo
2397) call DeallocateArray(temp_real_array)
2398) case('PERIODIC_OBSERVATION')
2399) output_option%print_observation = PETSC_TRUE
2400) call InputReadWord(input,option,word,PETSC_TRUE)
2401) call InputErrorMsg(input,option,'time increment', &
2402) 'OUTPUT, PERIODIC_OBSERVATION')
2403) call StringToUpper(word)
2404) select case(trim(word))
2405) case('TIME')
2406) call InputReadDouble(input,option,temp_real)
2407) call InputErrorMsg(input,option,'time increment', &
2408) 'OUTPUT,PERIODIC_OBSERVATION,TIME')
2409) call InputReadWord(input,option,word,PETSC_TRUE)
2410) call InputErrorMsg(input,option,'time increment units', &
2411) 'OUTPUT,PERIODIC_OBSERVATION,TIME')
2412) internal_units = 'sec'
2413) units_conversion = UnitsConvertToInternal(word, &
2414) internal_units,option)
2415) output_option%periodic_obs_output_time_incr = temp_real* &
2416) units_conversion
2417) case('TIMESTEP')
2418) call InputReadInt(input,option, &
2419) output_option%periodic_obs_output_ts_imod)
2420) call InputErrorMsg(input,option,'timestep increment', &
2421) 'OUTPUT,PERIODIC_OBSERVATION,TIMESTEP')
2422) case default
2423) call InputKeywordUnrecognized(word, &
2424) 'OUTPUT,PERIODIC_OBSERVATION',option)
2425) end select
2426) case('FORMAT')
2427) call InputReadWord(input,option,word,PETSC_TRUE)
2428) call InputErrorMsg(input,option,'keyword','OUTPUT,FORMAT')
2429) call StringToUpper(word)
2430) select case(trim(word))
2431) case ('HDF5')
2432) output_option%print_hdf5 = PETSC_TRUE
2433) call InputReadWord(input,option,word,PETSC_TRUE)
2434) call InputDefaultMsg(input,option, &
2435) 'OUTPUT,FORMAT,HDF5,# FILES')
2436) if (len_trim(word) > 0) then
2437) call StringToUpper(word)
2438) select case(trim(word))
2439) case('SINGLE_FILE')
2440) output_option%print_single_h5_file = PETSC_TRUE
2441) case('MULTIPLE_FILES')
2442) output_option%print_single_h5_file = PETSC_FALSE
2443) output_option%times_per_h5_file = 1
2444) call InputReadWord(input,option,word,PETSC_TRUE)
2445) if (len_trim(word)>0) then
2446) select case(trim(word))
2447) case('TIMES_PER_FILE')
2448) call InputReadInt(input,option, &
2449) output_option%times_per_h5_file)
2450) call InputErrorMsg(input,option, &
2451) 'timestep increment', &
2452) 'OUTPUT,FORMAT,HDF5,MULTIPLE_FILES,TIMES_PER_FILE')
2453) case default
2454) call InputKeywordUnrecognized(word, &
2455) 'OUTPUT,FORMAT,HDF5,MULTIPLE_FILES',option)
2456) end select
2457) endif
2458) case default
2459) call InputKeywordUnrecognized(word, &
2460) 'OUTPUT,FORMAT,HDF5',option)
2461) end select
2462) endif
2463) case ('MAD')
2464) output_option%print_mad = PETSC_TRUE
2465) case ('TECPLOT')
2466) output_option%print_tecplot = PETSC_TRUE
2467) call InputReadWord(input,option,word,PETSC_TRUE)
2468) call InputErrorMsg(input,option,'TECPLOT','OUTPUT,FORMAT')
2469) call StringToUpper(word)
2470) select case(trim(word))
2471) case('POINT')
2472) output_option%tecplot_format = TECPLOT_POINT_FORMAT
2473) case('BLOCK')
2474) output_option%tecplot_format = TECPLOT_BLOCK_FORMAT
2475) case('FEBRICK')
2476) output_option%tecplot_format = TECPLOT_FEBRICK_FORMAT
2477) case default
2478) call InputKeywordUnrecognized(word, &
2479) 'OUTPUT,FORMAT,TECPLOT',option)
2480) end select
2481) if (output_option%tecplot_format == TECPLOT_POINT_FORMAT &
2482) .and. option%mycommsize > 1) then
2483) output_option%tecplot_format = TECPLOT_BLOCK_FORMAT
2484) endif
2485) if (grid%itype == IMPLICIT_UNSTRUCTURED_GRID) then
2486) output_option%tecplot_format = TECPLOT_FEBRICK_FORMAT
2487) endif
2488) case ('VTK')
2489) output_option%print_vtk = PETSC_TRUE
2490) case default
2491) call InputKeywordUnrecognized(word,'OUTPUT,FORMAT',option)
2492) end select
2493) case('VELOCITY_AT_CENTER')
2494) vel_cent = PETSC_TRUE
2495) case('VELOCITY_AT_FACE')
2496) vel_face = PETSC_TRUE
2497) case('FLUXES')
2498) fluxes = PETSC_TRUE
2499) case('FLOWRATES','FLOWRATE')
2500) mass_flowrate = PETSC_TRUE
2501) energy_flowrate = PETSC_TRUE
2502) case('MASS_FLOWRATE')
2503) mass_flowrate = PETSC_TRUE
2504) case('ENERGY_FLOWRATE')
2505) energy_flowrate = PETSC_TRUE
2506) case('AVERAGE_FLOWRATES','AVERAGE_FLOWRATE')
2507) aveg_mass_flowrate = PETSC_TRUE
2508) aveg_energy_flowrate = PETSC_TRUE
2509) case('AVERAGE_MASS_FLOWRATE')
2510) aveg_mass_flowrate = PETSC_TRUE
2511) case('AVERAGE_ENERGY_FLOWRATE')
2512) aveg_energy_flowrate = PETSC_TRUE
2513) case ('HDF5_WRITE_GROUP_SIZE')
2514) call InputReadInt(input,option,option%hdf5_write_group_size)
2515) call InputErrorMsg(input,option,'HDF5_WRITE_GROUP_SIZE','Group size')
2516) case default
2517) call InputKeywordUnrecognized(word,'OUTPUT',option)
2518) end select
2519)
2520) enddo
2521)
2522) ! If VARIABLES were not specified within the *_FILE blocks, point their
2523) ! variable lists to the master variable list, which can be specified within
2524) ! the OUTPUT block. If no VARIABLES are specified for the master list, the
2525) ! defaults will be populated.
2526) if (.not.associated(output_option%output_snap_variable_list%first)) &
2527) then
2528) call OutputVariableListDestroy( &
2529) output_option%output_snap_variable_list)
2530) output_option%output_snap_variable_list => &
2531) output_option%output_variable_list
2532) endif
2533) if (.not.associated(output_option%output_obs_variable_list%first)) &
2534) then
2535) call OutputVariableListDestroy( &
2536) output_option%output_obs_variable_list)
2537) output_option%output_obs_variable_list => &
2538) output_option%output_variable_list
2539) endif
2540)
2541) if (vel_cent) then
2542) if (output_option%print_tecplot) &
2543) output_option%print_tecplot_vel_cent = PETSC_TRUE
2544) if (output_option%print_hdf5) &
2545) output_option%print_hdf5_vel_cent = PETSC_TRUE
2546) if (output_option%print_vtk) &
2547) output_option%print_vtk_vel_cent = PETSC_TRUE
2548) endif
2549) if (vel_face) then
2550) if (output_option%print_tecplot) &
2551) output_option%print_tecplot_vel_face = PETSC_TRUE
2552) if (output_option%print_hdf5) &
2553) output_option%print_hdf5_vel_face = PETSC_TRUE
2554) endif
2555) if (fluxes) then
2556) output_option%print_fluxes = PETSC_TRUE
2557) endif
2558) if(output_option%aveg_output_variable_list%nvars>0) then
2559) if(output_option%periodic_snap_output_time_incr==0.d0) then
2560) option%io_buffer = 'Keyword: AVERAGE_VARIABLES defined without' // &
2561) ' PERIODIC TIME being set.'
2562) call printErrMsg(option)
2563) endif
2564) if(.not.output_option%print_hdf5) then
2565) option%io_buffer = 'Keyword: AVERAGE_VARIABLES only defined for FORMAT HDF5'
2566) call printErrMsg(option)
2567) endif
2568) endif
2569) if (mass_flowrate.or.energy_flowrate.or.aveg_mass_flowrate &
2570) .or.aveg_energy_flowrate) then
2571) if (output_option%print_hdf5) then
2572) output_option%print_hdf5_mass_flowrate = mass_flowrate
2573) output_option%print_hdf5_energy_flowrate = energy_flowrate
2574) output_option%print_hdf5_aveg_mass_flowrate = aveg_mass_flowrate
2575) output_option%print_hdf5_aveg_energy_flowrate = aveg_energy_flowrate
2576) if(aveg_mass_flowrate.or.aveg_energy_flowrate) then
2577) if(output_option%periodic_snap_output_time_incr==0.d0) then
2578) option%io_buffer = 'Keyword: AVEGRAGE_FLOWRATES/ ' // &
2579) 'AVEGRAGE_MASS_FLOWRATE/ENERGY_FLOWRATE defined without' // &
2580) ' PERIODIC TIME being set.'
2581) call printErrMsg(option)
2582) endif
2583) endif
2584) option%flow%store_fluxes = PETSC_TRUE
2585) endif
2586) if (associated(grid%unstructured_grid)) then
2587) if (associated(grid%unstructured_grid%explicit_grid)) then
2588) option%flow%store_fluxes = PETSC_TRUE
2589) output_option%print_explicit_flowrate = mass_flowrate
2590) endif
2591) endif
2592) endif
2593)
2594) !.....................
2595) case ('REGRESSION')
2596) call RegressionRead(simulation%regression,input,option)
2597)
2598) !.....................
2599) case ('TIME')
2600) ! dt_init = UNINITIALIZED_DOUBLE
2601) dt_init = 1.d0
2602) dt_min = UNINITIALIZED_DOUBLE
2603) do
2604) call InputReadPflotranString(input,option)
2605) call InputReadStringErrorMsg(input,option,card)
2606) if (InputCheckExit(input,option)) exit
2607) call InputReadWord(input,option,word,PETSC_TRUE)
2608) call InputErrorMsg(input,option,'word','TIME')
2609) select case(trim(word))
2610) case('STEADY_STATE')
2611) option%steady_state = PETSC_TRUE
2612) case('FINAL_TIME')
2613) call InputReadDouble(input,option,temp_real)
2614) call InputErrorMsg(input,option,'Final Time','TIME')
2615) call InputReadWord(input,option,word,PETSC_TRUE)
2616) call InputErrorMsg(input,option,'Final Time Units','TIME')
2617) internal_units = 'sec'
2618) temp_real2 = UnitsConvertToInternal(word,internal_units,option)
2619) if (len_trim(output_option%tunit) == 0) then
2620) output_option%tunit = trim(word)
2621) output_option%tconv = temp_real2
2622) endif
2623) waypoint => WaypointCreate()
2624) waypoint%final = PETSC_TRUE
2625) waypoint%time = temp_real*temp_real2
2626) waypoint%print_snap_output = PETSC_TRUE
2627) call WaypointInsertInList(waypoint,waypoint_list)
2628) case('INITIAL_TIMESTEP_SIZE')
2629) call InputReadDouble(input,option,temp_real)
2630) call InputErrorMsg(input,option,'Initial Timestep Size','TIME')
2631) call InputReadWord(input,option,word,PETSC_TRUE)
2632) call InputErrorMsg(input,option,'Initial Timestep Size Time &
2633) &Units','TIME')
2634) internal_units = 'sec'
2635) dt_init = temp_real*UnitsConvertToInternal(word, &
2636) internal_units,option)
2637) case('MINIMUM_TIMESTEP_SIZE')
2638) call InputReadDouble(input,option,temp_real)
2639) call InputErrorMsg(input,option,'Minimum Timestep Size','TIME')
2640) call InputReadWord(input,option,word,PETSC_TRUE)
2641) call InputErrorMsg(input,option,'Minimum Timestep Size Time &
2642) &Units','TIME')
2643) internal_units = 'sec'
2644) dt_min = temp_real*UnitsConvertToInternal(word, &
2645) internal_units,option)
2646) case('MAXIMUM_TIMESTEP_SIZE')
2647) call InputReadDouble(input,option,temp_real)
2648) call InputErrorMsg(input,option,'Maximum Timestep Size','TIME')
2649) call InputReadWord(input,option,word,PETSC_TRUE)
2650) call InputErrorMsg(input,option,'Maximum Timestep Size Time &
2651) &Units','TIME')
2652) waypoint => WaypointCreate()
2653) internal_units = 'sec'
2654) waypoint%dt_max = temp_real*UnitsConvertToInternal(word, &
2655) internal_units,option)
2656) call InputReadWord(input,option,word,PETSC_TRUE)
2657) if (input%ierr == 0) then
2658) call StringToUpper(word)
2659) if (StringCompare(word,'AT',TWO_INTEGER)) then
2660) call InputReadDouble(input,option,temp_real)
2661) call InputErrorMsg(input,option,'Maximum Timestep Size &
2662) &Update Time','TIME')
2663) call InputReadWord(input,option,word,PETSC_TRUE)
2664) call InputErrorMsg(input,option,'Maximum Timestep Size &
2665) &Update Time Units','TIME')
2666) internal_units = 'sec'
2667) waypoint%time = temp_real*UnitsConvertToInternal(word, &
2668) internal_units,option)
2669) else
2670) option%io_buffer = 'Keyword under "MAXIMUM_TIMESTEP_SIZE" &
2671) &after maximum timestep size should &
2672) &be "at".'
2673) call printErrMsg(option)
2674) endif
2675) else
2676) waypoint%time = 0.d0
2677) endif
2678) call WaypointInsertInList(waypoint,waypoint_list)
2679) case default
2680) call InputKeywordUnrecognized(word,'TIME',option)
2681) end select
2682) enddo
2683) if (Initialized(dt_init)) then
2684) if (associated(flow_timestepper)) then
2685) flow_timestepper%dt_init = dt_init
2686) option%flow_dt = dt_init
2687) endif
2688) if (associated(tran_timestepper)) then
2689) tran_timestepper%dt_init = dt_init
2690) option%tran_dt = dt_init
2691) endif
2692) endif
2693) if (Initialized(dt_min)) then
2694) if (associated(flow_timestepper)) then
2695) flow_timestepper%dt_min = dt_min
2696) endif
2697) if (associated(tran_timestepper)) then
2698) tran_timestepper%dt_min = dt_min
2699) endif
2700) endif
2701)
2702) !......................
2703) case ('HDF5_READ_GROUP_SIZE')
2704) call InputReadInt(input,option,option%hdf5_read_group_size)
2705) call InputErrorMsg(input,option,'HDF5_READ_GROUP_SIZE','Group size')
2706)
2707) !......................
2708) case ('HDF5_WRITE_GROUP_SIZE')
2709) call InputReadInt(input,option,option%hdf5_write_group_size)
2710) call InputErrorMsg(input,option,'HDF5_WRITE_GROUP_SIZE','Group size')
2711)
2712) !....................
2713) case('WIPP')
2714) wipp => WIPPGetPtr()
2715) call WIPPRead(input,option)
2716)
2717) !....................
2718) case('KLINKENBERG_EFFECT')
2719) wipp => WIPPGetPtr()
2720) call KlinkenbergInit()
2721) klinkenberg => KlinkenbergCreate()
2722) call Klinkenberg%Read(input,option)
2723)
2724) !....................
2725) case ('ONLY_VERTICAL_FLOW')
2726) option%flow%only_vertical_flow = PETSC_TRUE
2727) if (option%iflowmode /= TH_MODE .and. &
2728) option%iflowmode /= RICHARDS_MODE) then
2729) option%io_buffer = 'ONLY_VERTICAL_FLOW implemented in RICHARDS and TH mode.'
2730) call printErrMsg(option)
2731) endif
2732)
2733) !....................
2734) case ('QUASI_3D')
2735) option%flow%quasi_3d = PETSC_TRUE
2736) option%flow%only_vertical_flow = PETSC_TRUE
2737) if (option%iflowmode /= RICHARDS_MODE) then
2738) option%io_buffer = 'QUASI_3D implemented in RICHARDS mode.'
2739) call printErrMsg(option)
2740) endif
2741)
2742) !....................
2743) case ('ONLY_ENERGY_EQ')
2744) option%flow%only_energy_eq = PETSC_TRUE
2745) if (option%iflowmode /= TH_MODE) then
2746) option%io_buffer = 'ONLY_ENERGY_EQ applicable only in TH mode.'
2747) call printErrMsg(option)
2748) endif
2749)
2750) !....................
2751) case ('RELATIVE_PERMEABILITY_AVERAGE')
2752) call InputReadWord(input,option,word,PETSC_FALSE)
2753) call StringToUpper(word)
2754) select case (trim(word))
2755) case ('UPWIND')
2756) option%rel_perm_aveg = UPWIND
2757) case ('HARMONIC')
2758) option%rel_perm_aveg = HARMONIC
2759) case ('DYNAMIC_HARMONIC')
2760) option%rel_perm_aveg = DYNAMIC_HARMONIC
2761) case default
2762) option%io_buffer = 'Cannot identify the specificed ' // &
2763) 'RELATIVE_PERMEABILITY_AVERAGE.'
2764) call printErrMsg(option)
2765) end select
2766)
2767) !....................
2768) case ('DBASE_FILENAME')
2769)
2770) !....................
2771) case ('END_SUBSURFACE')
2772) exit
2773)
2774) !....................
2775) case ('MIN_ALLOWABLE_SCALE')
2776) call InputReadDouble(input,option,option%min_allowable_scale)
2777) call InputErrorMsg(input,option,'minimium allowable scaling factor', &
2778) 'InitSubsurface')
2779)
2780) !....................
2781) case default
2782) call InputKeywordUnrecognized(word,'SubsurfaceReadInput()',option)
2783) end select
2784)
2785) enddo
2786)
2787) if (associated(simulation%flow_process_model_coupler)) then
2788) flow_timestepper%name = 'FLOW'
2789) if (option%steady_state) call TimestepperSteadyCreateFromBE(flow_timestepper)
2790) simulation%flow_process_model_coupler%timestepper => flow_timestepper
2791) else
2792) call flow_timestepper%Destroy()
2793) deallocate(flow_timestepper)
2794) nullify(flow_timestepper)
2795) endif
2796) if (associated(simulation%rt_process_model_coupler)) then
2797) tran_timestepper%name = 'TRAN'
2798) if (option%steady_state) call TimestepperSteadyCreateFromBE(tran_timestepper)
2799) simulation%rt_process_model_coupler%timestepper => tran_timestepper
2800) else
2801) call tran_timestepper%Destroy()
2802) deallocate(tran_timestepper)
2803) nullify(tran_timestepper)
2804) endif
2805)
2806) end subroutine SubsurfaceReadInput
2807)
2808) end module Factory_Subsurface_module