MODFLOW 6  version 6.8.0.dev0
USGS Modular Hydrologic Model
mf6coremodule Module Reference

Core MODFLOW 6 module. More...

Functions/Subroutines

subroutine mf6run
 Main controller. More...
 
subroutine mf6initialize ()
 Initialize a simulation. More...
 
logical(lgp) function mf6update ()
 Run a time step. More...
 
subroutine mf6finalize ()
 Finalize the simulation. More...
 
subroutine print_info ()
 print initial message More...
 
subroutine create_lstfile ()
 Set up mfsim list file output logging. More...
 
subroutine static_input_load ()
 Create simulation input context. More...
 
subroutine simulation_df ()
 Define the simulation. More...
 
subroutine simulation_ar ()
 Simulation allocate and read. More...
 
subroutine connections_cr ()
 Create the model connections from the exchanges. More...
 
subroutine mf6preparetimestep ()
 Read and prepare time step. More...
 
subroutine mf6dotimestep ()
 Run time step. More...
 
subroutine sim_step_retry (finishedTrying)
 Rerun time step. More...
 
logical(lgp) function mf6finalizetimestep ()
 Finalize time step. More...
 

Variables

class(runcontroltype), pointer run_ctrl => null()
 the run controller for this simulation More...
 

Detailed Description

This module contains the core components for MODFLOW 6. This module is used by the stand-alone executable and the share object versions of MODFLOW 6.

Function/Subroutine Documentation

◆ connections_cr()

subroutine mf6coremodule::connections_cr

This will upgrade the numerical exchanges in the solution, whenever the configuration requires this, to Connection objects. Currently we anticipate:

GWF-GWF => GwfGwfConnection GWT-GWT => GwtGwtConecction

Definition at line 451 of file mf6core.f90.

453  use simvariablesmodule, only: iout
454  use versionmodule, only: idevelopmode
455  integer(I4B) :: isol
456  type(ConnectionBuilderType) :: connectionBuilder
457  class(BaseSolutionType), pointer :: sol => null()
458  integer(I4B) :: status
459  character(len=16) :: envvar
460 
461  write (iout, '(/a)') 'PROCESSING MODEL CONNECTIONS'
462 
463  if (baseexchangelist%Count() == 0) then
464  ! if this is not a coupled simulation in any way,
465  ! then we will not need model connections
466  return
467  end if
468 
469  if (idevelopmode == 1) then
470  call get_environment_variable('DEV_ALWAYS_USE_IFMOD', &
471  value=envvar, status=status)
472  if (status == 0 .and. envvar == '1') then
473  connectionbuilder%dev_always_ifmod = .true.
474  write (iout, '(/a)') "Development option: forcing interface model"
475  end if
476  end if
477 
478  do isol = 1, basesolutionlist%Count()
479  sol => getbasesolutionfromlist(basesolutionlist, isol)
480  call connectionbuilder%processSolution(sol)
481  end do
482 
483  write (iout, '(a)') 'END OF MODEL CONNECTIONS'
This module contains simulation variables.
Definition: SimVariables.f90:9
integer(i4b) iout
file unit number for simulation output
This module contains version information.
Definition: version.f90:7
integer(i4b), parameter idevelopmode
Definition: version.f90:19
Here is the call graph for this function:
Here is the caller graph for this function:

◆ create_lstfile()

subroutine mf6coremodule::create_lstfile

This subroutine creates the mfsim list file and writes the header.

Definition at line 267 of file mf6core.f90.

268  use constantsmodule, only: linelength
271  use messagemodule, only: write_message
273  character(len=LINELENGTH) :: line
274  !
275  ! -- Open simulation list file
276  iout = getunit()
277  !
278  if (nr_procs > 1) then
280  end if
281  !
282  call openfile(iout, 0, simlstfile, 'LIST', filstat_opt='REPLACE')
283  !
284  ! -- write simlstfile to stdout
285  write (line, '(2(1x,A))') 'Writing simulation list file:', &
286  trim(adjustl(simlstfile))
287  !
288  call write_message(line)
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
integer(i4b) function, public getunit()
Get a free unit number.
subroutine, public append_processor_id(name, proc_id)
Append processor id to a string.
subroutine, public openfile(iu, iout, fname, ftype, fmtarg_opt, accarg_opt, filstat_opt, mode_opt)
Open a file.
Definition: InputOutput.f90:30
Store and issue logging messages to output units.
Definition: Message.f90:2
subroutine, public write_message(text, iunit, fmt, skipbefore, skipafter, advance)
Write a message to an output unit.
Definition: Message.f90:210
integer(i4b) nr_procs
character(len=linelength) simlstfile
simulation listing file name
integer(i4b) proc_id
subroutine write_listfile_header(iout, cmodel_type, write_sys_command, write_kind_info)
@ brief Write program header
Definition: version.f90:98
Here is the call graph for this function:
Here is the caller graph for this function:

◆ mf6dotimestep()

subroutine mf6coremodule::mf6dotimestep

This subroutine runs a single time step for the simulation. Steps include:

  • formulate the system of equations for each model and exchange
  • solve each solution

Definition at line 624 of file mf6core.f90.

625  ! -- modules
626  use kindmodule, only: i4b
627  use listsmodule, only: solutiongrouplist
630  use idmloadmodule, only: idm_ad
631  ! -- local variables
632  class(SolutionGroupType), pointer :: sgp => null()
633  integer(I4B) :: isg
634  logical :: finishedTrying
635 
636  ! start timer
637  call g_prof%start("Do time step", g_prof%tmr_do_tstp)
638 
639  ! -- By default, the solution groups will be solved once, and
640  ! may fail. But if adaptive stepping is active, then
641  ! the solution groups may be solved over and over with
642  ! progressively smaller time steps to see if convergence
643  ! can be obtained.
644  ifailedstepretry = 0
645  retryloop: do
646 
647  if (ifailedstepretry > 0) then
648  ! advance IDM
649  call idm_ad()
650  end if
651 
652  do isg = 1, solutiongrouplist%Count()
654  call sgp%sgp_ca()
655  end do
656 
657  call sim_step_retry(finishedtrying)
658  if (finishedtrying) exit retryloop
660 
661  end do retryloop
662 
663  ! stop timer
664  call g_prof%stop(g_prof%tmr_do_tstp)
665 
This module contains the IdmLoadModule.
Definition: IdmLoad.f90:7
subroutine, public idm_ad()
advance package dynamic data for period steps
Definition: IdmLoad.f90:62
This module defines variable data types.
Definition: kind.f90:8
type(listtype), public solutiongrouplist
Definition: mf6lists.f90:22
integer(i4b) ifailedstepretry
current retry for this time step
class(solutiongrouptype) function, pointer, public getsolutiongroupfromlist(list, idx)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ mf6finalize()

subroutine mf6coremodule::mf6finalize

This subroutine finalizes a simulation. Steps include:

  • final processing
  • deallocate memory

Definition at line 145 of file mf6core.f90.

146  ! -- modules
147  use, intrinsic :: iso_fortran_env, only: output_unit
148  use listsmodule, only: lists_da
150  use tdismodule, only: tdis_da
151  use idmloadmodule, only: idm_da
152  use sourceloadmodule, only: export_da
153  use simvariablesmodule, only: iout
154  ! -- local variables
155  integer(I4B) :: im
156  integer(I4B) :: ic
157  integer(I4B) :: is
158  integer(I4B) :: isg
159  class(SolutionGroupType), pointer :: sgp => null()
160  class(BaseSolutionType), pointer :: sp => null()
161  class(BaseModelType), pointer :: mp => null()
162  class(BaseExchangeType), pointer :: ep => null()
163  class(SpatialModelConnectionType), pointer :: mc => null()
164  integer(I4B) :: tmr_dealloc
165 
166  ! start timer
167  call g_prof%start("Finalize", g_prof%tmr_finalize)
168 
169  !
170  ! -- FINAL PROCESSING (FP)
171  ! -- Final processing for each model
172  do im = 1, basemodellist%Count()
173  mp => getbasemodelfromlist(basemodellist, im)
174  call mp%model_fp()
175  end do
176  !
177  ! -- Final processing for each exchange
178  do ic = 1, baseexchangelist%Count()
179  ep => getbaseexchangefromlist(baseexchangelist, ic)
180  call ep%exg_fp()
181  end do
182  !
183  ! -- Final processing for each solution
184  do is = 1, basesolutionlist%Count()
185  sp => getbasesolutionfromlist(basesolutionlist, is)
186  call sp%sln_fp()
187  end do
188 
189  ! start timer for deallocation
190  tmr_dealloc = -1
191  call g_prof%start("Deallocate", tmr_dealloc)
192 
193  !
194  ! -- DEALLOCATE (DA)
195  ! -- Deallocate tdis
196  call tdis_da()
197  !
198  ! -- Deallocate for each model
199  do im = 1, basemodellist%Count()
200  mp => getbasemodelfromlist(basemodellist, im)
201  call mp%model_da()
202  deallocate (mp)
203  end do
204  !
205  ! -- Deallocate for each exchange
206  do ic = 1, baseexchangelist%Count()
207  ep => getbaseexchangefromlist(baseexchangelist, ic)
208  call ep%exg_da()
209  deallocate (ep)
210  end do
211  !
212  ! -- Deallocate for each connection
213  do ic = 1, baseconnectionlist%Count()
214  mc => get_smc_from_list(baseconnectionlist, ic)
215  call mc%exg_da()
216  deallocate (mc)
217  end do
218  !
219  ! -- Deallocate for each solution
220  do is = 1, basesolutionlist%Count()
221  sp => getbasesolutionfromlist(basesolutionlist, is)
222  call sp%sln_da()
223  deallocate (sp)
224  end do
225  !
226  ! -- Deallocate solution group and simulation variables
227  do isg = 1, solutiongrouplist%Count()
228  sgp => getsolutiongroupfromlist(solutiongrouplist, isg)
229  call sgp%sgp_da()
230  deallocate (sgp)
231  end do
232  !
233  call idm_da(iout)
234  call export_da()
235  call simulation_da()
236  call lists_da()
237 
238  ! stop timer
239  call g_prof%stop(tmr_dealloc)
240 
241  ! finish gently (No calls after this)
242  ! timer is stopped inside because this call does not return
243  call run_ctrl%finish()
244 
subroutine, public idm_da(iout)
idm deallocate routine
Definition: IdmLoad.f90:74
subroutine, public lists_da()
Definition: mf6lists.f90:34
subroutine, public simulation_da()
Deallocate simulation variables.
This module contains the SourceLoadModule.
Definition: SourceLoad.F90:8
subroutine, public export_da()
deallocate model export objects and list
Definition: SourceLoad.F90:340
subroutine, public tdis_da()
Deallocate memory.
Definition: tdis.f90:345
Here is the call graph for this function:
Here is the caller graph for this function:

◆ mf6finalizetimestep()

logical(lgp) function mf6coremodule::mf6finalizetimestep

This function finalizes a single time step for the simulation and writes output for the time step. Steps include:

  • write output for each model
  • write output for each exchange
  • write output for each solutions
  • perform a final convergence check and whether the simulation can continue if convergence was not achieved
Returns
hasConverged boolean indicating if convergence was achieved for the time step

Definition at line 714 of file mf6core.f90.

715  ! -- modules
716  use kindmodule, only: i4b
722  use simmodule, only: converge_check
723  use simvariablesmodule, only: isim_mode
725  ! -- return variable
726  logical(LGP) :: hasConverged
727  ! -- local variables
728  class(BaseSolutionType), pointer :: sp => null()
729  class(BaseModelType), pointer :: mp => null()
730  class(BaseExchangeType), pointer :: ep => null()
731  class(SpatialModelConnectionType), pointer :: mc => null()
732  character(len=LINELENGTH) :: line
733  character(len=LINELENGTH) :: fmt
734  integer(I4B) :: im
735  integer(I4B) :: ix
736  integer(I4B) :: ic
737  integer(I4B) :: is
738  !
739  ! -- initialize format and line
740  fmt = "(/,a,/)"
741  line = 'end timestep'
742 
743  ! start timer
744  call g_prof%start("Finalize time step", g_prof%tmr_final_tstp)
745 
746  !
747  ! -- evaluate simulation mode
748  select case (isim_mode)
749  case (mvalidate)
750  !
751  ! -- Write final message for timestep for each model
752  do im = 1, basemodellist%Count()
754  call mp%model_message(line, fmt=fmt)
755  end do
756  case (mnormal)
757 
758  call g_prof%start("Write output", g_prof%tmr_output)
759  !
760  ! -- Write output and final message for timestep for each model
761  do im = 1, basemodellist%Count()
763  call mp%model_ot()
764  call mp%model_message(line, fmt=fmt)
765  end do
766  !
767  ! -- Write output for each exchange
768  do ix = 1, baseexchangelist%Count()
770  call ep%exg_ot()
771  end do
772  !
773  ! -- Write output for each connection
774  do ic = 1, baseconnectionlist%Count()
775  mc => get_smc_from_list(baseconnectionlist, ic)
776  call mc%exg_ot()
777  end do
778  !
779  ! -- Write output for each solution
780  do is = 1, basesolutionlist%Count()
782  call sp%sln_ot()
783  end do
784  !
785  ! -- update exports
786  call g_prof%start("NetCDF export", g_prof%tmr_nc_export)
787  call export_post_step()
788  call g_prof%stop(g_prof%tmr_nc_export)
789 
790  call g_prof%stop(g_prof%tmr_output)
791  end select
792  !
793  ! -- Check if we're done
794  call converge_check(hasconverged)
795 
796  ! stop timer
797  call g_prof%stop(g_prof%tmr_final_tstp)
798 
class(baseexchangetype) function, pointer, public getbaseexchangefromlist(list, idx)
Retrieve a specific BaseExchangeType object from a list.
class(basemodeltype) function, pointer, public getbasemodelfromlist(list, idx)
Definition: BaseModel.f90:171
class(basesolutiontype) function, pointer, public getbasesolutionfromlist(list, idx)
@ mvalidate
validation mode - do not run time steps
Definition: Constants.f90:205
@ mnormal
normal output mode
Definition: Constants.f90:206
type(listtype), public basemodellist
Definition: mf6lists.f90:16
type(listtype), public baseexchangelist
Definition: mf6lists.f90:25
type(listtype), public basesolutionlist
Definition: mf6lists.f90:19
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public converge_check(hasConverged)
Simulation convergence check.
Definition: Sim.f90:402
integer(i4b) isim_mode
simulation mode
subroutine, public export_post_step()
model exports post step actions
Definition: SourceLoad.F90:333
Highest level model type. All models extend this parent type.
Definition: BaseModel.f90:16
Here is the call graph for this function:
Here is the caller graph for this function:

◆ mf6initialize()

subroutine mf6coremodule::mf6initialize

This subroutine initializes a MODFLOW 6 simulation. The subroutine:

  • creates the simulation
  • defines
  • allocates and reads static data

Definition at line 70 of file mf6core.f90.

71  ! -- modules
74  use sourceloadmodule, only: export_cr
75 
76  call g_prof%pre_init()
77 
78  ! -- get the run controller for sequential or parallel builds
79  run_ctrl => create_run_control()
80  call run_ctrl%start()
81 
82  ! -- print info and start timer
83  call print_info()
84 
85  ! -- create mfsim.lst
86  call create_lstfile()
87 
88  ! -- load input context
89  call static_input_load()
90 
91  ! init timer and start
92  call g_prof%initialize()
93  call g_prof%start("Initialize", g_prof%tmr_init)
94 
95  ! -- create
96  call simulation_cr()
97 
98  ! -- define
99  call simulation_df()
100 
101  ! -- allocate and read
102  call simulation_ar()
103 
104  ! -- create model exports
105  call export_cr()
106 
107  ! -- stop the timer
108  call g_prof%stop(g_prof%tmr_init)
109 
class(runcontroltype) function, pointer, public create_run_control()
subroutine, public simulation_cr()
Read the simulation name file and initialize the models, exchanges.
subroutine, public export_cr()
create model exports list
Definition: SourceLoad.F90:305
Here is the call graph for this function:
Here is the caller graph for this function:

◆ mf6preparetimestep()

subroutine mf6coremodule::mf6preparetimestep

This subroutine reads and prepares period data for the simulation. Steps include:

  • read and prepare for each model
  • read and prepare for each exchange
  • reset convergence flag
  • calculate maximum time step for each model
  • calculate maximum time step for each exchange
  • calculate maximum time step for each solution
  • set time discretization timestep using smallest maximum timestep

Definition at line 499 of file mf6core.f90.

500  ! -- modules
501  use kindmodule, only: i4b
504  kstp, kper
509  use simmodule, only: converge_reset
510  use simvariablesmodule, only: isim_mode
511  use idmloadmodule, only: idm_rp, idm_ad
513  ! -- local variables
514  class(BaseModelType), pointer :: mp => null()
515  class(BaseExchangeType), pointer :: ep => null()
516  class(SpatialModelConnectionType), pointer :: mc => null()
517  class(BaseSolutionType), pointer :: sp => null()
518  character(len=LINELENGTH) :: line
519  character(len=LINELENGTH) :: fmt
520  integer(I4B) :: im
521  integer(I4B) :: ie
522  integer(I4B) :: ic
523  integer(I4B) :: is
524 
525  ! start timer
526  call g_prof%start("Prepare time step", g_prof%tmr_prep_tstp)
527 
528  !
529  ! -- initialize fmt
530  fmt = "(/,a,/)"
531  !
532  ! -- period update
533  call tdis_set_counters()
534  !
535  ! -- set base line
536  write (line, '(a,i0,a,i0,a)') &
537  'start timestep kper="', kper, '" kstp="', kstp, '" mode="'
538  !
539  ! -- evaluate simulation mode
540  select case (isim_mode)
541  case (mvalidate)
542  line = trim(line)//'validate"'
543  case (mnormal)
544  line = trim(line)//'normal"'
545  end select
546 
547  ! -- load dynamic input
548  call idm_rp()
549 
550  ! -- Read and prepare each model
551  do im = 1, basemodellist%Count()
553  call mp%model_message(line, fmt=fmt)
554  call mp%model_rp()
555  end do
556  !
557  ! -- Synchronize
558  call run_ctrl%at_stage(stg_bfr_exg_rp)
559  !
560  ! -- Read and prepare each exchange
561  do ie = 1, baseexchangelist%Count()
563  call ep%exg_rp()
564  end do
565  !
566  ! -- Read and prepare each connection
567  do ic = 1, baseconnectionlist%Count()
568  mc => get_smc_from_list(baseconnectionlist, ic)
569  call mc%exg_rp()
570  end do
571  !
572  ! -- Synchronize
573  call run_ctrl%at_stage(stg_aft_con_rp)
574  !
575  ! -- reset simulation convergence flag
576  call converge_reset()
577  !
578  ! -- time update for each model
579  do im = 1, basemodellist%Count()
581  call mp%model_dt()
582  end do
583  !
584  ! -- time update for each exchange
585  do ie = 1, baseexchangelist%Count()
587  call ep%exg_dt()
588  end do
589  !
590  ! -- time update for each connection
591  do ic = 1, baseconnectionlist%Count()
592  mc => get_smc_from_list(baseconnectionlist, ic)
593  call mc%exg_dt()
594  end do
595  !
596  ! -- time update for each solution
597  do is = 1, basesolutionlist%Count()
598  sp => getbasesolutionfromlist(basesolutionlist, is)
599  call sp%sln_dt()
600  end do
601  !
602  ! -- update exports
603  call export_post_prepare()
604  !
605  ! -- set time step
606  call tdis_set_timestep()
607 
608  ! advance IDM
609  call idm_ad()
610 
611  ! stop timer
612  call g_prof%stop(g_prof%tmr_prep_tstp)
613 
subroutine, public idm_rp()
load package dynamic data for period
Definition: IdmLoad.f90:50
subroutine, public converge_reset()
Reset the simulation convergence flag.
Definition: Sim.f90:389
subroutine, public export_post_prepare()
model exports post prepare step actions
Definition: SourceLoad.F90:326
subroutine, public tdis_set_timestep()
Set time step length.
Definition: tdis.f90:153
subroutine, public tdis_set_counters()
Set kstp and kper.
Definition: tdis.f90:91
integer(i4b), pointer, public kstp
current time step number
Definition: tdis.f90:24
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
Here is the call graph for this function:
Here is the caller graph for this function:

◆ mf6run()

subroutine mf6coremodule::mf6run

This subroutine is the main controller for MODFLOW 6.

Definition at line 33 of file mf6core.f90.

34  ! -- modules
36  use tdismodule, only: endofsimulation
37  ! -- local
38  logical(LGP) :: hasConverged
39  !
40  ! -- parse any command line arguments
42  !
43  ! initialize simulation
44  call mf6initialize()
45  !
46  ! -- time loop
47  do while (.not. endofsimulation)
48 
49  ! perform a time step
50  hasconverged = mf6update()
51 
52  ! if not converged, break
53  if (.not. hasconverged) exit
54 
55  end do
56  !
57  ! -- finalize simulation
58  call mf6finalize()
59 
subroutine, public getcommandlinearguments()
Get command line arguments.
Definition: comarg.f90:29
logical(lgp), pointer, public endofsimulation
flag indicating end of simulation
Definition: tdis.f90:28
Here is the call graph for this function:
Here is the caller graph for this function:

◆ mf6update()

logical(lgp) function mf6coremodule::mf6update

This function runs a single time step to completion.

Returns
hasConverged boolean indicating if convergence was achieved for the time step

Definition at line 119 of file mf6core.f90.

120  logical(LGP) :: hasConverged
121  ! start timer
122  call g_prof%start("Update", g_prof%tmr_update)
123  !
124  ! -- prepare timestep
125  call mf6preparetimestep()
126  !
127  ! -- do timestep
128  call mf6dotimestep()
129  !
130  ! -- after timestep
131  hasconverged = mf6finalizetimestep()
132 
133  ! stop timer
134  call g_prof%stop(g_prof%tmr_update)
135 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ print_info()

subroutine mf6coremodule::print_info

Definition at line 249 of file mf6core.f90.

250  use simmodule, only: initial_message
251  use timermodule, only: print_start_time
252 
253  ! print initial message
254  call initial_message()
255 
256  ! get start time
257  call print_start_time()
258 
subroutine, public initial_message()
Print the header and initializes messaging.
Definition: Sim.f90:442
subroutine, public print_start_time()
Start simulation timer.
Definition: Timer.f90:19
Here is the call graph for this function:
Here is the caller graph for this function:

◆ sim_step_retry()

subroutine mf6coremodule::sim_step_retry ( logical, intent(out)  finishedTrying)

This subroutine reruns a single time step for the simulation when the adaptive time step option is used.

Parameters
[out]finishedtryingboolean that indicates if no

Definition at line 674 of file mf6core.f90.

675  ! -- modules
676  use kindmodule, only: dp
678  use simmodule, only: converge_reset
679  use tdismodule, only: kstp, kper, delt, tdis_delt_reset
681  ! -- dummy variables
682  logical, intent(out) :: finishedTrying !< boolean that indicates if no
683  ! additional reruns of the time step are required
684  !
685  ! -- Check with ats to reset delt and keep trying
686  finishedtrying = .true.
687  call ats_reset_delt(kstp, kper, laststepfailed, delt, finishedtrying)
688  !
689  if (.not. finishedtrying) then
690  !
691  ! -- Reset delt, which requires updating pertim, totim
692  ! and end of period and simulation indicators
693  call tdis_delt_reset(delt)
694  !
695  ! -- Reset state of the simulation convergence flag
696  call converge_reset()
697 
698  end if
subroutine, public ats_reset_delt(kstp, kper, lastStepFailed, delt, finishedTrying)
@ brief Reset time step because failure has occurred
Definition: ats.f90:606
integer(i4b) laststepfailed
flag indicating if the last step failed (1) if last step failed; (0) otherwise (set in converge_check...
subroutine, public tdis_delt_reset(deltnew)
Reset delt and update timing variables and indicators.
Definition: tdis.f90:213
real(dp), pointer, public delt
length of the current time step
Definition: tdis.f90:29
Here is the call graph for this function:
Here is the caller graph for this function:

◆ simulation_ar()

subroutine mf6coremodule::simulation_ar

This subroutine allocates and reads static data for the simulation. Steps include:

  • allocate and read for each model
  • allocate and read for each exchange
  • allocate and read for each solution

Definition at line 399 of file mf6core.f90.

401  ! -- local variables
402  integer(I4B) :: im
403  integer(I4B) :: ic
404  integer(I4B) :: is
405  class(BaseSolutionType), pointer :: sp => null()
406  class(BaseModelType), pointer :: mp => null()
407  class(BaseExchangeType), pointer :: ep => null()
408  class(SpatialModelConnectionType), pointer :: mc => null()
409 
410  ! -- Allocate and read each model
411  do im = 1, basemodellist%Count()
412  mp => getbasemodelfromlist(basemodellist, im)
413  call mp%model_ar()
414  end do
415  !
416  ! -- Allocate and read each exchange
417  do ic = 1, baseexchangelist%Count()
418  ep => getbaseexchangefromlist(baseexchangelist, ic)
419  call ep%exg_ar()
420  end do
421  !
422  ! -- Synchronize
423  call run_ctrl%at_stage(stg_bfr_con_ar)
424  !
425  ! -- Allocate and read all model connections
426  do ic = 1, baseconnectionlist%Count()
427  mc => get_smc_from_list(baseconnectionlist, ic)
428  call mc%exg_ar()
429  end do
430  !
431  ! -- Synchronize
432  call run_ctrl%at_stage(stg_aft_con_ar)
433  !
434  ! -- Allocate and read each solution
435  do is = 1, basesolutionlist%Count()
436  sp => getbasesolutionfromlist(basesolutionlist, is)
437  call sp%sln_ar()
438  end do
439  !
Here is the call graph for this function:
Here is the caller graph for this function:

◆ simulation_df()

subroutine mf6coremodule::simulation_df

This subroutine defined the simulation. Steps include:

  • define each model
  • define each solution

Definition at line 327 of file mf6core.f90.

328  ! -- modules
329  use idmloadmodule, only: idm_df
330  ! -- local variables
331  integer(I4B) :: im
332  integer(I4B) :: ic
333  integer(I4B) :: is
334  class(BaseSolutionType), pointer :: sp => null()
335  class(BaseModelType), pointer :: mp => null()
336  class(BaseExchangeType), pointer :: ep => null()
337  class(SpatialModelConnectionType), pointer :: mc => null()
338 
339  ! -- init virtual data environment
340  call run_ctrl%at_stage(stg_bfr_mdl_df)
341 
342  ! -- Define each model
343  do im = 1, basemodellist%Count()
344  mp => getbasemodelfromlist(basemodellist, im)
345  call mp%model_df()
346  end do
347  !
348  ! -- synchronize
349  call run_ctrl%at_stage(stg_aft_mdl_df)
350  !
351  ! -- Define each exchange
352  do ic = 1, baseexchangelist%Count()
353  ep => getbaseexchangefromlist(baseexchangelist, ic)
354  call ep%exg_df()
355  end do
356  !
357  ! -- synchronize
358  call run_ctrl%at_stage(stg_aft_exg_df)
359  !
360  ! -- when needed, this is where the interface models are
361  ! created and added to the numerical solutions
362  call connections_cr()
363  !
364  ! -- synchronize
365  call run_ctrl%at_stage(stg_aft_con_cr)
366  !
367  ! -- synchronize
368  call run_ctrl%at_stage(stg_bfr_con_df)
369  !
370  ! -- Define each connection
371  do ic = 1, baseconnectionlist%Count()
372  mc => get_smc_from_list(baseconnectionlist, ic)
373  call mc%exg_df()
374  end do
375  !
376  ! -- synchronize
377  call run_ctrl%at_stage(stg_aft_con_df)
378  !
379  ! -- Define each solution
380  do is = 1, basesolutionlist%Count()
381  sp => getbasesolutionfromlist(basesolutionlist, is)
382  call sp%sln_df()
383  end do
384 
385  ! idm df
386  call idm_df()
387 
subroutine, public idm_df()
advance package dynamic data for period steps
Definition: IdmLoad.f90:38
Here is the call graph for this function:
Here is the caller graph for this function:

◆ static_input_load()

subroutine mf6coremodule::static_input_load

This subroutine creates the simulation input context

Definition at line 297 of file mf6core.f90.

298  ! -- modules
299  use constantsmodule, only: lenmempath
300  use simvariablesmodule, only: iout
301  use idmloadmodule, only: simnam_load, simtdis_load, &
305  use simvariablesmodule, only: iparamlog
306  !
307  ! -- load simnam input context
308  call simnam_load(iparamlog)
309  !
310  ! -- load tdis to input context
311  call simtdis_load()
312  !
313  ! -- load in scope models
314  call load_models(iout)
315  !
316  ! -- load in scope exchanges
317  call load_exchanges(iout)
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:27
subroutine, public simnam_load(paramlog)
MODFLOW 6 mfsim.nam input load routine.
Definition: IdmLoad.f90:442
subroutine, public simtdis_load()
MODFLOW 6 tdis input load routine.
Definition: IdmLoad.f90:457
subroutine, public load_models(iout)
load model namfiles and model package files
Definition: IdmLoad.f90:218
subroutine, public load_exchanges(iout)
load exchange files
Definition: IdmLoad.f90:285
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
integer(i4b) iparamlog
input (idm) parameter logging to simulation listing file
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ run_ctrl

class(runcontroltype), pointer mf6coremodule::run_ctrl => null()

Definition at line 24 of file mf6core.f90.

24  class(RunControlType), pointer :: run_ctrl => null() !< the run controller for this simulation