MODFLOW 6  version 6.6.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 424 of file mf6core.f90.

426  use simvariablesmodule, only: iout
427  use versionmodule, only: idevelopmode
428  integer(I4B) :: isol
429  type(ConnectionBuilderType) :: connectionBuilder
430  class(BaseSolutionType), pointer :: sol => null()
431  integer(I4B) :: status
432  character(len=16) :: envvar
433 
434  write (iout, '(/a)') 'PROCESSING MODEL CONNECTIONS'
435 
436  if (baseexchangelist%Count() == 0) then
437  ! if this is not a coupled simulation in any way,
438  ! then we will not need model connections
439  return
440  end if
441 
442  if (idevelopmode == 1) then
443  call get_environment_variable('DEV_ALWAYS_USE_IFMOD', &
444  value=envvar, status=status)
445  if (status == 0 .and. envvar == '1') then
446  connectionbuilder%dev_always_ifmod = .true.
447  write (iout, '(/a)') "Development option: forcing interface model"
448  end if
449  end if
450 
451  do isol = 1, basesolutionlist%Count()
452  sol => getbasesolutionfromlist(basesolutionlist, isol)
453  call connectionbuilder%processSolution(sol)
454  end do
455 
456  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 240 of file mf6core.f90.

241  use constantsmodule, only: linelength
244  use messagemodule, only: write_message
246  character(len=LINELENGTH) :: line
247  !
248  ! -- Open simulation list file
249  iout = getunit()
250  !
251  if (nr_procs > 1) then
253  end if
254  !
255  call openfile(iout, 0, simlstfile, 'LIST', filstat_opt='REPLACE')
256  !
257  ! -- write simlstfile to stdout
258  write (line, '(2(1x,A))') 'Writing simulation list file:', &
259  trim(adjustl(simlstfile))
260  !
261  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 587 of file mf6core.f90.

588  ! -- modules
589  use kindmodule, only: i4b
590  use listsmodule, only: solutiongrouplist
593  use idmloadmodule, only: idm_ad
594  ! -- local variables
595  class(SolutionGroupType), pointer :: sgp => null()
596  integer(I4B) :: isg
597  logical :: finishedTrying
598 
599  ! -- By default, the solution groups will be solved once, and
600  ! may fail. But if adaptive stepping is active, then
601  ! the solution groups may be solved over and over with
602  ! progressively smaller time steps to see if convergence
603  ! can be obtained.
604  ifailedstepretry = 0
605  retryloop: do
606 
607  ! -- idm advance
608  call idm_ad()
609 
610  do isg = 1, solutiongrouplist%Count()
612  call sgp%sgp_ca()
613  end do
614 
615  call sim_step_retry(finishedtrying)
616  if (finishedtrying) exit retryloop
618 
619  end do retryloop
620 
This module contains the IdmLoadModule.
Definition: IdmLoad.f90:7
subroutine, public idm_ad()
advance package dynamic data for period steps
Definition: IdmLoad.f90:70
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 131 of file mf6core.f90.

132  ! -- modules
133  use, intrinsic :: iso_fortran_env, only: output_unit
134  use listsmodule, only: lists_da
136  use tdismodule, only: tdis_da
137  use idmloadmodule, only: idm_da
138  use sourceloadmodule, only: export_da
139  use simvariablesmodule, only: iout
140  ! -- local variables
141  integer(I4B) :: im
142  integer(I4B) :: ic
143  integer(I4B) :: is
144  integer(I4B) :: isg
145  class(SolutionGroupType), pointer :: sgp => null()
146  class(BaseSolutionType), pointer :: sp => null()
147  class(BaseModelType), pointer :: mp => null()
148  class(BaseExchangeType), pointer :: ep => null()
149  class(SpatialModelConnectionType), pointer :: mc => null()
150  !
151  !
152  ! -- FINAL PROCESSING (FP)
153  ! -- Final processing for each model
154  do im = 1, basemodellist%Count()
155  mp => getbasemodelfromlist(basemodellist, im)
156  call mp%model_fp()
157  end do
158  !
159  ! -- Final processing for each exchange
160  do ic = 1, baseexchangelist%Count()
161  ep => getbaseexchangefromlist(baseexchangelist, ic)
162  call ep%exg_fp()
163  end do
164  !
165  ! -- Final processing for each solution
166  do is = 1, basesolutionlist%Count()
167  sp => getbasesolutionfromlist(basesolutionlist, is)
168  call sp%sln_fp()
169  end do
170  !
171  ! -- DEALLOCATE (DA)
172  ! -- Deallocate tdis
173  call tdis_da()
174  !
175  ! -- Deallocate for each model
176  do im = 1, basemodellist%Count()
177  mp => getbasemodelfromlist(basemodellist, im)
178  call mp%model_da()
179  deallocate (mp)
180  end do
181  !
182  ! -- Deallocate for each exchange
183  do ic = 1, baseexchangelist%Count()
184  ep => getbaseexchangefromlist(baseexchangelist, ic)
185  call ep%exg_da()
186  deallocate (ep)
187  end do
188  !
189  ! -- Deallocate for each connection
190  do ic = 1, baseconnectionlist%Count()
191  mc => get_smc_from_list(baseconnectionlist, ic)
192  call mc%exg_da()
193  deallocate (mc)
194  end do
195  !
196  ! -- Deallocate for each solution
197  do is = 1, basesolutionlist%Count()
198  sp => getbasesolutionfromlist(basesolutionlist, is)
199  call sp%sln_da()
200  deallocate (sp)
201  end do
202  !
203  ! -- Deallocate solution group and simulation variables
204  do isg = 1, solutiongrouplist%Count()
205  sgp => getsolutiongroupfromlist(solutiongrouplist, isg)
206  call sgp%sgp_da()
207  deallocate (sgp)
208  end do
209  !
210  call idm_da(iout)
211  call export_da()
212  call simulation_da()
213  call lists_da()
214  !
215  ! -- finish gently (No calls after this)
216  call run_ctrl%finish()
217  !
subroutine, public idm_da(iout)
idm deallocate routine
Definition: IdmLoad.f90:86
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:385
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 669 of file mf6core.f90.

670  ! -- modules
671  use kindmodule, only: i4b
677  use simmodule, only: converge_check
678  use simvariablesmodule, only: isim_mode
680  ! -- return variable
681  logical(LGP) :: hasConverged
682  ! -- local variables
683  class(BaseSolutionType), pointer :: sp => null()
684  class(BaseModelType), pointer :: mp => null()
685  class(BaseExchangeType), pointer :: ep => null()
686  class(SpatialModelConnectionType), pointer :: mc => null()
687  character(len=LINELENGTH) :: line
688  character(len=LINELENGTH) :: fmt
689  integer(I4B) :: im
690  integer(I4B) :: ix
691  integer(I4B) :: ic
692  integer(I4B) :: is
693  !
694  ! -- initialize format and line
695  fmt = "(/,a,/)"
696  line = 'end timestep'
697  !
698  ! -- evaluate simulation mode
699  select case (isim_mode)
700  case (mvalidate)
701  !
702  ! -- Write final message for timestep for each model
703  do im = 1, basemodellist%Count()
705  call mp%model_message(line, fmt=fmt)
706  end do
707  case (mnormal)
708  !
709  ! -- Write output and final message for timestep for each model
710  do im = 1, basemodellist%Count()
712  call mp%model_ot()
713  call mp%model_message(line, fmt=fmt)
714  end do
715  !
716  ! -- Write output for each exchange
717  do ix = 1, baseexchangelist%Count()
719  call ep%exg_ot()
720  end do
721  !
722  ! -- Write output for each connection
723  do ic = 1, baseconnectionlist%Count()
724  mc => get_smc_from_list(baseconnectionlist, ic)
725  call mc%exg_ot()
726  end do
727  !
728  ! -- Write output for each solution
729  do is = 1, basesolutionlist%Count()
731  call sp%sln_ot()
732  end do
733  !
734  ! -- update exports
735  call export_post_step()
736  end select
737  !
738  ! -- Check if we're done
739  call converge_check(hasconverged)
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:172
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:377
Highest level model type. All models extend this parent type.
Definition: BaseModel.f90:13
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 69 of file mf6core.f90.

70  ! -- modules
73  use sourceloadmodule, only: export_cr
74 
75  ! -- get the run controller for sequential or parallel builds
76  run_ctrl => create_run_control()
77  call run_ctrl%start()
78 
79  ! -- print info and start timer
80  call print_info()
81 
82  ! -- create mfsim.lst
83  call create_lstfile()
84 
85  ! -- load input context
86  call static_input_load()
87 
88  ! -- create
89  call simulation_cr()
90 
91  ! -- define
92  call simulation_df()
93 
94  ! -- allocate and read
95  call simulation_ar()
96 
97  ! -- create model exports
98  call export_cr()
99 
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:346
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 472 of file mf6core.f90.

473  ! -- modules
474  use kindmodule, only: i4b
477  kstp, kper
482  use simmodule, only: converge_reset
483  use simvariablesmodule, only: isim_mode
484  use idmloadmodule, only: idm_rp
486  ! -- local variables
487  class(BaseModelType), pointer :: mp => null()
488  class(BaseExchangeType), pointer :: ep => null()
489  class(SpatialModelConnectionType), pointer :: mc => null()
490  class(BaseSolutionType), pointer :: sp => null()
491  character(len=LINELENGTH) :: line
492  character(len=LINELENGTH) :: fmt
493  integer(I4B) :: im
494  integer(I4B) :: ie
495  integer(I4B) :: ic
496  integer(I4B) :: is
497  !
498  ! -- initialize fmt
499  fmt = "(/,a,/)"
500  !
501  ! -- period update
502  call tdis_set_counters()
503  !
504  ! -- set base line
505  write (line, '(a,i0,a,i0,a)') &
506  'start timestep kper="', kper, '" kstp="', kstp, '" mode="'
507  !
508  ! -- evaluate simulation mode
509  select case (isim_mode)
510  case (mvalidate)
511  line = trim(line)//'validate"'
512  case (mnormal)
513  line = trim(line)//'normal"'
514  end select
515 
516  ! -- load dynamic input
517  call idm_rp()
518 
519  ! -- Read and prepare each model
520  do im = 1, basemodellist%Count()
522  call mp%model_message(line, fmt=fmt)
523  call mp%model_rp()
524  end do
525  !
526  ! -- Synchronize
527  call run_ctrl%at_stage(stg_bfr_exg_rp)
528  !
529  ! -- Read and prepare each exchange
530  do ie = 1, baseexchangelist%Count()
532  call ep%exg_rp()
533  end do
534  !
535  ! -- Read and prepare each connection
536  do ic = 1, baseconnectionlist%Count()
537  mc => get_smc_from_list(baseconnectionlist, ic)
538  call mc%exg_rp()
539  end do
540  !
541  ! -- Synchronize
542  call run_ctrl%at_stage(stg_aft_con_rp)
543  !
544  ! -- reset simulation convergence flag
545  call converge_reset()
546  !
547  ! -- time update for each model
548  do im = 1, basemodellist%Count()
550  call mp%model_dt()
551  end do
552  !
553  ! -- time update for each exchange
554  do ie = 1, baseexchangelist%Count()
556  call ep%exg_dt()
557  end do
558  !
559  ! -- time update for each connection
560  do ic = 1, baseconnectionlist%Count()
561  mc => get_smc_from_list(baseconnectionlist, ic)
562  call mc%exg_dt()
563  end do
564  !
565  ! -- time update for each solution
566  do is = 1, basesolutionlist%Count()
567  sp => getbasesolutionfromlist(basesolutionlist, is)
568  call sp%sln_dt()
569  end do
570  !
571  ! -- update exports
572  call export_post_prepare()
573  !
574  ! -- set time step
575  call tdis_set_timestep()
576 
subroutine, public idm_rp()
load package dynamic data for period
Definition: IdmLoad.f90:54
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:369
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 32 of file mf6core.f90.

33  ! -- modules
35  use tdismodule, only: endofsimulation
36  ! -- local
37  logical(LGP) :: hasConverged
38  !
39  ! -- parse any command line arguments
41  !
42  ! initialize simulation
43  call mf6initialize()
44  !
45  ! -- time loop
46  do while (.not. endofsimulation)
47 
48  ! perform a time step
49  hasconverged = mf6update()
50 
51  ! if not converged, break
52  if (.not. hasconverged) exit
53 
54  end do
55  !
56  ! -- finalize simulation
57  call mf6finalize()
58 
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 109 of file mf6core.f90.

110  ! -- return variable
111  logical(LGP) :: hasConverged
112  !
113  ! -- prepare timestep
114  call mf6preparetimestep()
115  !
116  ! -- do timestep
117  call mf6dotimestep()
118  !
119  ! -- after timestep
120  hasconverged = mf6finalizetimestep()
121  !
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 222 of file mf6core.f90.

223  use simmodule, only: initial_message
224  use timermodule, only: print_start_time
225 
226  ! print initial message
227  call initial_message()
228 
229  ! get start time
230  call print_start_time()
231 
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 629 of file mf6core.f90.

630  ! -- modules
631  use kindmodule, only: dp
633  use simmodule, only: converge_reset
634  use tdismodule, only: kstp, kper, delt, tdis_delt_reset
636  ! -- dummy variables
637  logical, intent(out) :: finishedTrying !< boolean that indicates if no
638  ! additional reruns of the time step are required
639  !
640  ! -- Check with ats to reset delt and keep trying
641  finishedtrying = .true.
642  call ats_reset_delt(kstp, kper, laststepfailed, delt, finishedtrying)
643  !
644  if (.not. finishedtrying) then
645  !
646  ! -- Reset delt, which requires updating pertim, totim
647  ! and end of period and simulation indicators
648  call tdis_delt_reset(delt)
649  !
650  ! -- Reset state of the simulation convergence flag
651  call converge_reset()
652 
653  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 372 of file mf6core.f90.

374  ! -- local variables
375  integer(I4B) :: im
376  integer(I4B) :: ic
377  integer(I4B) :: is
378  class(BaseSolutionType), pointer :: sp => null()
379  class(BaseModelType), pointer :: mp => null()
380  class(BaseExchangeType), pointer :: ep => null()
381  class(SpatialModelConnectionType), pointer :: mc => null()
382 
383  ! -- Allocate and read each model
384  do im = 1, basemodellist%Count()
385  mp => getbasemodelfromlist(basemodellist, im)
386  call mp%model_ar()
387  end do
388  !
389  ! -- Allocate and read each exchange
390  do ic = 1, baseexchangelist%Count()
391  ep => getbaseexchangefromlist(baseexchangelist, ic)
392  call ep%exg_ar()
393  end do
394  !
395  ! -- Synchronize
396  call run_ctrl%at_stage(stg_bfr_con_ar)
397  !
398  ! -- Allocate and read all model connections
399  do ic = 1, baseconnectionlist%Count()
400  mc => get_smc_from_list(baseconnectionlist, ic)
401  call mc%exg_ar()
402  end do
403  !
404  ! -- Synchronize
405  call run_ctrl%at_stage(stg_aft_con_ar)
406  !
407  ! -- Allocate and read each solution
408  do is = 1, basesolutionlist%Count()
409  sp => getbasesolutionfromlist(basesolutionlist, is)
410  call sp%sln_ar()
411  end do
412  !
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 300 of file mf6core.f90.

301  ! -- modules
302  use idmloadmodule, only: idm_df
303  ! -- local variables
304  integer(I4B) :: im
305  integer(I4B) :: ic
306  integer(I4B) :: is
307  class(BaseSolutionType), pointer :: sp => null()
308  class(BaseModelType), pointer :: mp => null()
309  class(BaseExchangeType), pointer :: ep => null()
310  class(SpatialModelConnectionType), pointer :: mc => null()
311 
312  ! -- init virtual data environment
313  call run_ctrl%at_stage(stg_bfr_mdl_df)
314 
315  ! -- Define each model
316  do im = 1, basemodellist%Count()
317  mp => getbasemodelfromlist(basemodellist, im)
318  call mp%model_df()
319  end do
320  !
321  ! -- synchronize
322  call run_ctrl%at_stage(stg_aft_mdl_df)
323  !
324  ! -- Define each exchange
325  do ic = 1, baseexchangelist%Count()
326  ep => getbaseexchangefromlist(baseexchangelist, ic)
327  call ep%exg_df()
328  end do
329  !
330  ! -- synchronize
331  call run_ctrl%at_stage(stg_aft_exg_df)
332  !
333  ! -- when needed, this is where the interface models are
334  ! created and added to the numerical solutions
335  call connections_cr()
336  !
337  ! -- synchronize
338  call run_ctrl%at_stage(stg_aft_con_cr)
339  !
340  ! -- synchronize
341  call run_ctrl%at_stage(stg_bfr_con_df)
342  !
343  ! -- Define each connection
344  do ic = 1, baseconnectionlist%Count()
345  mc => get_smc_from_list(baseconnectionlist, ic)
346  call mc%exg_df()
347  end do
348  !
349  ! -- synchronize
350  call run_ctrl%at_stage(stg_aft_con_df)
351  !
352  ! -- Define each solution
353  do is = 1, basesolutionlist%Count()
354  sp => getbasesolutionfromlist(basesolutionlist, is)
355  call sp%sln_df()
356  end do
357 
358  ! idm df
359  call idm_df()
360 
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 270 of file mf6core.f90.

271  ! -- modules
272  use constantsmodule, only: lenmempath
273  use simvariablesmodule, only: iout
274  use idmloadmodule, only: simnam_load, simtdis_load, &
278  use simvariablesmodule, only: iparamlog
279  !
280  ! -- load simnam input context
281  call simnam_load(iparamlog)
282  !
283  ! -- load tdis to input context
284  call simtdis_load()
285  !
286  ! -- load in scope models
287  call load_models(iout)
288  !
289  ! -- load in scope exchanges
290  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:493
subroutine, public simtdis_load()
MODFLOW 6 tdis input load routine.
Definition: IdmLoad.f90:515
subroutine, public load_models(iout)
load model namfiles and model package files
Definition: IdmLoad.f90:247
subroutine, public load_exchanges(iout)
load exchange files
Definition: IdmLoad.f90:326
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 23 of file mf6core.f90.

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