MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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 441 of file mf6core.f90.

443  use simvariablesmodule, only: iout
444  use versionmodule, only: idevelopmode
445  integer(I4B) :: isol
446  type(ConnectionBuilderType) :: connectionBuilder
447  class(BaseSolutionType), pointer :: sol => null()
448  integer(I4B) :: status
449  character(len=16) :: envvar
450 
451  write (iout, '(/a)') 'PROCESSING MODEL CONNECTIONS'
452 
453  if (baseexchangelist%Count() == 0) then
454  ! if this is not a coupled simulation in any way,
455  ! then we will not need model connections
456  return
457  end if
458 
459  if (idevelopmode == 1) then
460  call get_environment_variable('DEV_ALWAYS_USE_IFMOD', &
461  value=envvar, status=status)
462  if (status == 0 .and. envvar == '1') then
463  connectionbuilder%dev_always_ifmod = .true.
464  write (iout, '(/a)') "Development option: forcing interface model"
465  end if
466  end if
467 
468  do isol = 1, basesolutionlist%Count()
469  sol => getbasesolutionfromlist(basesolutionlist, isol)
470  call connectionbuilder%processSolution(sol)
471  end do
472 
473  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 234 of file mf6core.f90.

235  use constantsmodule, only: linelength
238  use messagemodule, only: write_message
240  character(len=LINELENGTH) :: line
241  !
242  ! -- Open simulation list file
243  iout = getunit()
244  !
245  if (nr_procs > 1) then
247  end if
248  !
249  call openfile(iout, 0, simlstfile, 'LIST', filstat_opt='REPLACE')
250  !
251  ! -- write simlstfile to stdout
252  write (line, '(2(1x,A))') 'Writing simulation list file:', &
253  trim(adjustl(simlstfile))
254  !
255  call write_message(line)
257  !
258  ! -- return
259  return
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:44
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 600 of file mf6core.f90.

601  ! -- modules
602  use kindmodule, only: i4b
603  use listsmodule, only: solutiongrouplist
606  use idmloadmodule, only: idm_ad
607  ! -- local variables
608  class(SolutionGroupType), pointer :: sgp => null()
609  integer(I4B) :: isg
610  logical :: finishedTrying
611 
612  ! -- By default, the solution groups will be solved once, and
613  ! may fail. But if adaptive stepping is active, then
614  ! the solution groups may be solved over and over with
615  ! progressively smaller time steps to see if convergence
616  ! can be obtained.
617  ifailedstepretry = 0
618  retryloop: do
619 
620  ! -- idm advance
621  call idm_ad()
622 
623  do isg = 1, solutiongrouplist%Count()
625  call sgp%sgp_ca()
626  end do
627 
628  call sim_step_retry(finishedtrying)
629  if (finishedtrying) exit retryloop
631 
632  end do retryloop
633 
This module contains the IdmLoadModule.
Definition: IdmLoad.f90:7
subroutine, public idm_ad()
advance package dynamic data for period steps
Definition: IdmLoad.f90:71
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 127 of file mf6core.f90.

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

686  ! -- modules
687  use kindmodule, only: i4b
693  use simmodule, only: converge_check
694  use simvariablesmodule, only: isim_mode
695  ! -- return variable
696  logical(LGP) :: hasConverged
697  ! -- local variables
698  class(BaseSolutionType), pointer :: sp => null()
699  class(BaseModelType), pointer :: mp => null()
700  class(BaseExchangeType), pointer :: ep => null()
701  class(SpatialModelConnectionType), pointer :: mc => null()
702  character(len=LINELENGTH) :: line
703  character(len=LINELENGTH) :: fmt
704  integer(I4B) :: im
705  integer(I4B) :: ix
706  integer(I4B) :: ic
707  integer(I4B) :: is
708  !
709  ! -- initialize format and line
710  fmt = "(/,a,/)"
711  line = 'end timestep'
712  !
713  ! -- evaluate simulation mode
714  select case (isim_mode)
715  case (mvalidate)
716  !
717  ! -- Write final message for timestep for each model
718  do im = 1, basemodellist%Count()
720  call mp%model_message(line, fmt=fmt)
721  end do
722  case (mnormal)
723  !
724  ! -- Write output and final message for timestep for each model
725  do im = 1, basemodellist%Count()
727  call mp%model_ot()
728  call mp%model_message(line, fmt=fmt)
729  end do
730  !
731  ! -- Write output for each exchange
732  do ix = 1, baseexchangelist%Count()
734  call ep%exg_ot()
735  end do
736  !
737  ! -- Write output for each connection
738  do ic = 1, baseconnectionlist%Count()
739  mc => get_smc_from_list(baseconnectionlist, ic)
740  call mc%exg_ot()
741  end do
742  !
743  ! -- Write output for each solution
744  do is = 1, basesolutionlist%Count()
746  call sp%sln_ot()
747  end do
748  end select
749  !
750  ! -- Check if we're done
751  call converge_check(hasconverged)
752  !
753  ! -- return
754  return
755 
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:204
@ mnormal
normal output mode
Definition: Constants.f90:205
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:405
integer(i4b) isim_mode
simulation mode
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 
74  ! -- get the run controller for sequential or parallel builds
75  run_ctrl => create_run_control()
76  call run_ctrl%start()
77 
78  ! -- print info and start timer
79  call print_info()
80 
81  ! -- create mfsim.lst
82  call create_lstfile()
83 
84  ! -- load input context
85  call static_input_load()
86 
87  ! -- create
88  call simulation_cr()
89 
90  ! -- define
91  call simulation_df()
92 
93  ! -- allocate and read
94  call simulation_ar()
95 
class(runcontroltype) function, pointer, public create_run_control()
subroutine, public simulation_cr()
Read the simulation name file and initialize the models, exchanges.
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 489 of file mf6core.f90.

490  ! -- modules
491  use kindmodule, only: i4b
494  kstp, kper
499  use simmodule, only: converge_reset
500  use simvariablesmodule, only: isim_mode
501  use idmloadmodule, only: idm_rp
502  ! -- local variables
503  class(BaseModelType), pointer :: mp => null()
504  class(BaseExchangeType), pointer :: ep => null()
505  class(SpatialModelConnectionType), pointer :: mc => null()
506  class(BaseSolutionType), pointer :: sp => null()
507  character(len=LINELENGTH) :: line
508  character(len=LINELENGTH) :: fmt
509  integer(I4B) :: im
510  integer(I4B) :: ie
511  integer(I4B) :: ic
512  integer(I4B) :: is
513  !
514  ! -- initialize fmt
515  fmt = "(/,a,/)"
516  !
517  ! -- period update
518  call tdis_set_counters()
519  !
520  ! -- set base line
521  write (line, '(a,i0,a,i0,a)') &
522  'start timestep kper="', kper, '" kstp="', kstp, '" mode="'
523  !
524  ! -- evaluate simulation mode
525  select case (isim_mode)
526  case (mvalidate)
527  line = trim(line)//'validate"'
528  case (mnormal)
529  line = trim(line)//'normal"'
530  end select
531 
532  ! -- load dynamic input
533  call idm_rp()
534 
535  ! -- Read and prepare each model
536  do im = 1, basemodellist%Count()
538  call mp%model_message(line, fmt=fmt)
539  call mp%model_rp()
540  end do
541  !
542  ! -- Synchronize
543  call run_ctrl%at_stage(stg_bfr_exg_rp)
544  !
545  ! -- Read and prepare each exchange
546  do ie = 1, baseexchangelist%Count()
548  call ep%exg_rp()
549  end do
550  !
551  ! -- Read and prepare each connection
552  do ic = 1, baseconnectionlist%Count()
553  mc => get_smc_from_list(baseconnectionlist, ic)
554  call mc%exg_rp()
555  end do
556  !
557  ! -- Synchronize
558  call run_ctrl%at_stage(stg_aft_con_rp)
559  !
560  ! -- reset simulation convergence flag
561  call converge_reset()
562  !
563  ! -- time update for each model
564  do im = 1, basemodellist%Count()
566  call mp%model_calculate_delt()
567  end do
568  !
569  ! -- time update for each exchange
570  do ie = 1, baseexchangelist%Count()
572  call ep%exg_calculate_delt()
573  end do
574  !
575  ! -- time update for each connection
576  do ic = 1, baseconnectionlist%Count()
577  mc => get_smc_from_list(baseconnectionlist, ic)
578  call mc%exg_calculate_delt()
579  end do
580  !
581  ! -- time update for each solution
582  do is = 1, basesolutionlist%Count()
583  sp => getbasesolutionfromlist(basesolutionlist, is)
584  call sp%sln_calculate_delt()
585  end do
586  !
587  ! -- set time step
588  call tdis_set_timestep()
589 
subroutine, public idm_rp()
load package dynamic data for period
Definition: IdmLoad.f90:55
subroutine, public converge_reset()
Reset the simulation convergence flag.
Definition: Sim.f90:392
subroutine, public tdis_set_timestep()
Set time step length.
Definition: tdis.f90:159
subroutine, public tdis_set_counters()
Set kstp and kper.
Definition: tdis.f90:94
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 105 of file mf6core.f90.

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

217  use simmodule, only: initial_message
218  use timermodule, only: print_start_time
219 
220  ! print initial message
221  call initial_message()
222 
223  ! get start time
224  call print_start_time()
225 
subroutine, public initial_message()
Print the header and initializes messaging.
Definition: Sim.f90:445
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 642 of file mf6core.f90.

643  ! -- modules
644  use kindmodule, only: dp
646  use simmodule, only: converge_reset
647  use tdismodule, only: kstp, kper, delt, tdis_delt_reset
649  ! -- dummy variables
650  logical, intent(out) :: finishedTrying !< boolean that indicates if no
651  ! additional reruns of the time step are required
652  !
653  ! -- Check with ats to reset delt and keep trying
654  finishedtrying = .true.
655  call ats_reset_delt(kstp, kper, laststepfailed, delt, finishedtrying)
656  !
657  if (.not. finishedtrying) then
658  !
659  ! -- Reset delt, which requires updating pertim, totim
660  ! and end of period and simulation indicators
661  call tdis_delt_reset(delt)
662  !
663  ! -- Reset state of the simulation convergence flag
664  call converge_reset()
665 
666  end if
667  !
668  ! -- return
669  return
subroutine, public ats_reset_delt(kstp, kper, lastStepFailed, delt, finishedTrying)
@ brief Reset time step because failure has occurred
Definition: ats.f90:633
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:222
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 read 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 389 of file mf6core.f90.

391  ! -- local variables
392  integer(I4B) :: im
393  integer(I4B) :: ic
394  integer(I4B) :: is
395  class(BaseSolutionType), pointer :: sp => null()
396  class(BaseModelType), pointer :: mp => null()
397  class(BaseExchangeType), pointer :: ep => null()
398  class(SpatialModelConnectionType), pointer :: mc => null()
399 
400  ! -- Allocate and read each model
401  do im = 1, basemodellist%Count()
402  mp => getbasemodelfromlist(basemodellist, im)
403  call mp%model_ar()
404  end do
405  !
406  ! -- Allocate and read each exchange
407  do ic = 1, baseexchangelist%Count()
408  ep => getbaseexchangefromlist(baseexchangelist, ic)
409  call ep%exg_ar()
410  end do
411  !
412  ! -- Synchronize
413  call run_ctrl%at_stage(stg_bfr_con_ar)
414  !
415  ! -- Allocate and read all model connections
416  do ic = 1, baseconnectionlist%Count()
417  mc => get_smc_from_list(baseconnectionlist, ic)
418  call mc%exg_ar()
419  end do
420  !
421  ! -- Synchronize
422  call run_ctrl%at_stage(stg_aft_con_ar)
423  !
424  ! -- Allocate and read each solution
425  do is = 1, basesolutionlist%Count()
426  sp => getbasesolutionfromlist(basesolutionlist, is)
427  call sp%sln_ar()
428  end do
429  !
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 317 of file mf6core.f90.

318  ! -- modules
319  use idmloadmodule, only: idm_df
320  ! -- local variables
321  integer(I4B) :: im
322  integer(I4B) :: ic
323  integer(I4B) :: is
324  class(BaseSolutionType), pointer :: sp => null()
325  class(BaseModelType), pointer :: mp => null()
326  class(BaseExchangeType), pointer :: ep => null()
327  class(SpatialModelConnectionType), pointer :: mc => null()
328 
329  ! -- init virtual data environment
330  call run_ctrl%at_stage(stg_bfr_mdl_df)
331 
332  ! -- Define each model
333  do im = 1, basemodellist%Count()
334  mp => getbasemodelfromlist(basemodellist, im)
335  call mp%model_df()
336  end do
337  !
338  ! -- synchronize
339  call run_ctrl%at_stage(stg_aft_mdl_df)
340  !
341  ! -- Define each exchange
342  do ic = 1, baseexchangelist%Count()
343  ep => getbaseexchangefromlist(baseexchangelist, ic)
344  call ep%exg_df()
345  end do
346  !
347  ! -- synchronize
348  call run_ctrl%at_stage(stg_aft_exg_df)
349  !
350  ! -- when needed, this is were the interface models are
351  ! created and added to the numerical solutions
352  call connections_cr()
353  !
354  ! -- synchronize
355  call run_ctrl%at_stage(stg_aft_con_cr)
356  !
357  ! -- synchronize
358  call run_ctrl%at_stage(stg_bfr_con_df)
359  !
360  ! -- Define each connection
361  do ic = 1, baseconnectionlist%Count()
362  mc => get_smc_from_list(baseconnectionlist, ic)
363  call mc%exg_df()
364  end do
365  !
366  ! -- synchronize
367  call run_ctrl%at_stage(stg_aft_con_df)
368  !
369  ! -- Define each solution
370  do is = 1, basesolutionlist%Count()
371  sp => getbasesolutionfromlist(basesolutionlist, is)
372  call sp%sln_df()
373  end do
374 
375  ! idm df
376  call idm_df()
377 
subroutine, public idm_df()
advance package dynamic data for period steps
Definition: IdmLoad.f90:39
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 267 of file mf6core.f90.

268  ! -- modules
269  use constantsmodule, only: lenmempath
270  use simvariablesmodule, only: iout
271  use idmloadmodule, only: simnam_load, simtdis_load, &
277  ! -- dummy
278  ! -- locals
279  character(len=LENMEMPATH) :: input_mempath
280  integer(I4B), dimension(:), pointer, contiguous :: model_loadmask
281  integer(I4B), pointer :: nummodels => null()
282  !
283  ! -- load simnam input context
284  call simnam_load(iparamlog)
285  !
286  ! -- load tdis to input context
287  call simtdis_load()
288  !
289  ! -- allocate model load mask
290  input_mempath = create_mem_path(component='SIM', context=idm_context)
291  call mem_setptr(nummodels, 'NUMMODELS', input_mempath)
292  allocate (model_loadmask(nummodels))
293  !
294  ! -- initialize mask
295  call create_load_mask(model_loadmask)
296  !
297  ! -- load in scope models
298  call load_models(model_loadmask, iout)
299  !
300  ! -- load in scope exchanges
301  call load_exchanges(model_loadmask, iout)
302  !
303  ! -- cleanup
304  deallocate (model_loadmask)
305  !
306  ! -- return
307  return
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:26
subroutine, public simnam_load(paramlog)
MODFLOW 6 mfsim.nam input load routine.
Definition: IdmLoad.f90:435
subroutine, public simtdis_load()
MODFLOW 6 tdis input load routine.
Definition: IdmLoad.f90:457
subroutine, public load_models(model_loadmask, iout)
load model namfiles and model package files
Definition: IdmLoad.f90:211
subroutine, public load_exchanges(model_loadmask, iout)
load exchange files
Definition: IdmLoad.f90:279
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
subroutine, public create_load_mask(mask_array)
Create a load mask to determine which models should be loaded by idm on this process....
character(len=linelength) idm_context
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