MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
mf6core.f90
Go to the documentation of this file.
1 !> @brief Core MODFLOW 6 module
2 !!
3 !! This module contains the core components for MODFLOW 6. This module
4 !! is used by the stand-alone executable and the share object versions
5 !! of MODFLOW 6.
6 !!
7 !<
9  use kindmodule, only: i4b, lgp
20  use simstagesmodule
21  implicit none
22 
23  class(runcontroltype), pointer :: run_ctrl => null() !< the run controller for this simulation
24 
25 contains
26 
27  !> @brief Main controller
28  !!
29  !! This subroutine is the main controller for MODFLOW 6.
30  !!
31  !<
32  subroutine mf6run
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 
59  end subroutine mf6run
60 
61  !> @brief Initialize a simulation
62  !!
63  !! This subroutine initializes a MODFLOW 6 simulation. The subroutine:
64  !! - creates the simulation
65  !! - defines
66  !! - allocates and reads static data
67  !!
68  !<
69  subroutine mf6initialize()
70  ! -- modules
73 
74  ! -- get the run controller for sequential or parallel builds
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 
96  end subroutine mf6initialize
97 
98  !> @brief Run a time step
99  !!
100  !! This function runs a single time step to completion.
101  !!
102  !! @return hasConverged boolean indicating if convergence was achieved for the time step
103  !!
104  !<
105  function mf6update() result(hasConverged)
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  !
118  end function mf6update
119 
120  !> @brief Finalize the simulation
121  !!
122  !! This subroutine finalizes a simulation. Steps include:
123  !! - final processing
124  !! - deallocate memory
125  !!
126  !<
127  subroutine mf6finalize()
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  !
212  end subroutine mf6finalize
213 
214  !> @brief print initial message
215  !<
216  subroutine print_info()
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 
226  end subroutine print_info
227 
228  !> @brief Set up mfsim list file output logging
229  !!
230  !! This subroutine creates the mfsim list file
231  !! and writes the header.
232  !!
233  !<
234  subroutine create_lstfile()
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
260  end subroutine create_lstfile
261 
262  !> @brief Create simulation input context
263  !!
264  !! This subroutine creates the simulation input context
265  !!
266  !<
267  subroutine static_input_load()
268  ! -- modules
269  use constantsmodule, only: lenmempath
270  use simvariablesmodule, only: iout
271  use idmloadmodule, only: simnam_load, simtdis_load, &
275  use simvariablesmodule, only: iparamlog
276  !
277  ! -- load simnam input context
278  call simnam_load(iparamlog)
279  !
280  ! -- load tdis to input context
281  call simtdis_load()
282  !
283  ! -- load in scope models
284  call load_models(iout)
285  !
286  ! -- load in scope exchanges
287  call load_exchanges(iout)
288  !
289  ! -- return
290  return
291  end subroutine static_input_load
292 
293  !> @brief Define the simulation
294  !!
295  !! This subroutine defined the simulation. Steps include:
296  !! - define each model
297  !! - define each solution
298  !!
299  !<
300  subroutine simulation_df()
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()
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()
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 were 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()
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()
355  call sp%sln_df()
356  end do
357 
358  ! idm df
359  call idm_df()
360 
361  end subroutine simulation_df
362 
363  !> @brief Simulation allocate and read
364  !!
365  !! This subroutine allocates and read static data for the simulation.
366  !! Steps include:
367  !! - allocate and read for each model
368  !! - allocate and read for each exchange
369  !! - allocate and read for each solution
370  !!
371  !<
372  subroutine simulation_ar()
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()
386  call mp%model_ar()
387  end do
388  !
389  ! -- Allocate and read each exchange
390  do ic = 1, baseexchangelist%Count()
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()
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()
410  call sp%sln_ar()
411  end do
412  !
413  end subroutine simulation_ar
414 
415  !> @brief Create the model connections from the exchanges
416  !!
417  !! This will upgrade the numerical exchanges in the solution,
418  !! whenever the configuration requires this, to Connection
419  !! objects. Currently we anticipate:
420  !!
421  !! GWF-GWF => GwfGwfConnection
422  !! GWT-GWT => GwtGwtConecction
423  !<
424  subroutine connections_cr()
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()
453  call connectionbuilder%processSolution(sol)
454  end do
455 
456  write (iout, '(a)') 'END OF MODEL CONNECTIONS'
457  end subroutine connections_cr
458 
459  !> @brief Read and prepare time step
460  !!
461  !! This subroutine reads and prepares period data for the simulation.
462  !! Steps include:
463  !! - read and prepare for each model
464  !! - read and prepare for each exchange
465  !! - reset convergence flag
466  !! - calculate maximum time step for each model
467  !! - calculate maximum time step for each exchange
468  !! - calculate maximum time step for each solution
469  !! - set time discretization timestep using smallest maximum timestep
470  !!
471  !<
472  subroutine mf6preparetimestep()
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
485  ! -- local variables
486  class(basemodeltype), pointer :: mp => null()
487  class(baseexchangetype), pointer :: ep => null()
488  class(spatialmodelconnectiontype), pointer :: mc => null()
489  class(basesolutiontype), pointer :: sp => null()
490  character(len=LINELENGTH) :: line
491  character(len=LINELENGTH) :: fmt
492  integer(I4B) :: im
493  integer(I4B) :: ie
494  integer(I4B) :: ic
495  integer(I4B) :: is
496  !
497  ! -- initialize fmt
498  fmt = "(/,a,/)"
499  !
500  ! -- period update
501  call tdis_set_counters()
502  !
503  ! -- set base line
504  write (line, '(a,i0,a,i0,a)') &
505  'start timestep kper="', kper, '" kstp="', kstp, '" mode="'
506  !
507  ! -- evaluate simulation mode
508  select case (isim_mode)
509  case (mvalidate)
510  line = trim(line)//'validate"'
511  case (mnormal)
512  line = trim(line)//'normal"'
513  end select
514 
515  ! -- load dynamic input
516  call idm_rp()
517 
518  ! -- Read and prepare each model
519  do im = 1, basemodellist%Count()
521  call mp%model_message(line, fmt=fmt)
522  call mp%model_rp()
523  end do
524  !
525  ! -- Synchronize
526  call run_ctrl%at_stage(stg_bfr_exg_rp)
527  !
528  ! -- Read and prepare each exchange
529  do ie = 1, baseexchangelist%Count()
531  call ep%exg_rp()
532  end do
533  !
534  ! -- Read and prepare each connection
535  do ic = 1, baseconnectionlist%Count()
536  mc => get_smc_from_list(baseconnectionlist, ic)
537  call mc%exg_rp()
538  end do
539  !
540  ! -- Synchronize
541  call run_ctrl%at_stage(stg_aft_con_rp)
542  !
543  ! -- reset simulation convergence flag
544  call converge_reset()
545  !
546  ! -- time update for each model
547  do im = 1, basemodellist%Count()
549  call mp%model_calculate_delt()
550  end do
551  !
552  ! -- time update for each exchange
553  do ie = 1, baseexchangelist%Count()
555  call ep%exg_calculate_delt()
556  end do
557  !
558  ! -- time update for each connection
559  do ic = 1, baseconnectionlist%Count()
560  mc => get_smc_from_list(baseconnectionlist, ic)
561  call mc%exg_calculate_delt()
562  end do
563  !
564  ! -- time update for each solution
565  do is = 1, basesolutionlist%Count()
566  sp => getbasesolutionfromlist(basesolutionlist, is)
567  call sp%sln_calculate_delt()
568  end do
569  !
570  ! -- set time step
571  call tdis_set_timestep()
572 
573  end subroutine mf6preparetimestep
574 
575  !> @brief Run time step
576  !!
577  !! This subroutine runs a single time step for the simulation.
578  !! Steps include:
579  !! - formulate the system of equations for each model and exchange
580  !! - solve each solution
581  !!
582  !<
583  subroutine mf6dotimestep()
584  ! -- modules
585  use kindmodule, only: i4b
586  use listsmodule, only: solutiongrouplist
589  use idmloadmodule, only: idm_ad
590  ! -- local variables
591  class(solutiongrouptype), pointer :: sgp => null()
592  integer(I4B) :: isg
593  logical :: finishedTrying
594 
595  ! -- By default, the solution groups will be solved once, and
596  ! may fail. But if adaptive stepping is active, then
597  ! the solution groups may be solved over and over with
598  ! progressively smaller time steps to see if convergence
599  ! can be obtained.
600  ifailedstepretry = 0
601  retryloop: do
602 
603  ! -- idm advance
604  call idm_ad()
605 
606  do isg = 1, solutiongrouplist%Count()
608  call sgp%sgp_ca()
609  end do
610 
611  call sim_step_retry(finishedtrying)
612  if (finishedtrying) exit retryloop
614 
615  end do retryloop
616 
617  end subroutine mf6dotimestep
618 
619  !> @brief Rerun time step
620  !!
621  !! This subroutine reruns a single time step for the simulation when
622  !! the adaptive time step option is used.
623  !!
624  !<
625  subroutine sim_step_retry(finishedTrying)
626  ! -- modules
627  use kindmodule, only: dp
629  use simmodule, only: converge_reset
630  use tdismodule, only: kstp, kper, delt, tdis_delt_reset
632  ! -- dummy variables
633  logical, intent(out) :: finishedTrying !< boolean that indicates if no
634  ! additional reruns of the time step are required
635  !
636  ! -- Check with ats to reset delt and keep trying
637  finishedtrying = .true.
638  call ats_reset_delt(kstp, kper, laststepfailed, delt, finishedtrying)
639  !
640  if (.not. finishedtrying) then
641  !
642  ! -- Reset delt, which requires updating pertim, totim
643  ! and end of period and simulation indicators
644  call tdis_delt_reset(delt)
645  !
646  ! -- Reset state of the simulation convergence flag
647  call converge_reset()
648 
649  end if
650  !
651  ! -- return
652  return
653  end subroutine sim_step_retry
654 
655  !> @brief Finalize time step
656  !!
657  !! This function finalizes a single time step for the simulation
658  !! and writes output for the time step. Steps include:
659  !! - write output for each model
660  !! - write output for each exchange
661  !! - write output for each solutions
662  !! - perform a final convergence check and whether the simulation
663  !! can continue if convergence was not achieved
664  !!
665  !! @return hasConverged boolean indicating if convergence was achieved for the time step
666  !!
667  !<
668  function mf6finalizetimestep() result(hasConverged)
669  ! -- modules
670  use kindmodule, only: i4b
676  use simmodule, only: converge_check
677  use simvariablesmodule, only: isim_mode
678  ! -- return variable
679  logical(LGP) :: hasconverged
680  ! -- local variables
681  class(basesolutiontype), pointer :: sp => null()
682  class(basemodeltype), pointer :: mp => null()
683  class(baseexchangetype), pointer :: ep => null()
684  class(spatialmodelconnectiontype), pointer :: mc => null()
685  character(len=LINELENGTH) :: line
686  character(len=LINELENGTH) :: fmt
687  integer(I4B) :: im
688  integer(I4B) :: ix
689  integer(I4B) :: ic
690  integer(I4B) :: is
691  !
692  ! -- initialize format and line
693  fmt = "(/,a,/)"
694  line = 'end timestep'
695  !
696  ! -- evaluate simulation mode
697  select case (isim_mode)
698  case (mvalidate)
699  !
700  ! -- Write final message for timestep for each model
701  do im = 1, basemodellist%Count()
703  call mp%model_message(line, fmt=fmt)
704  end do
705  case (mnormal)
706  !
707  ! -- Write output and final message for timestep for each model
708  do im = 1, basemodellist%Count()
710  call mp%model_ot()
711  call mp%model_message(line, fmt=fmt)
712  end do
713  !
714  ! -- Write output for each exchange
715  do ix = 1, baseexchangelist%Count()
717  call ep%exg_ot()
718  end do
719  !
720  ! -- Write output for each connection
721  do ic = 1, baseconnectionlist%Count()
722  mc => get_smc_from_list(baseconnectionlist, ic)
723  call mc%exg_ot()
724  end do
725  !
726  ! -- Write output for each solution
727  do is = 1, basesolutionlist%Count()
729  call sp%sln_ot()
730  end do
731  end select
732  !
733  ! -- Check if we're done
734  call converge_check(hasconverged)
735  !
736  ! -- return
737  return
738 
739  end function mf6finalizetimestep
740 
741 end module mf6coremodule
subroutine, public ats_reset_delt(kstp, kper, lastStepFailed, delt, finishedTrying)
@ brief Reset time step because failure has occurred
Definition: ats.f90:633
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)
subroutine, public getcommandlinearguments()
Get command line arguments.
Definition: comarg.f90:29
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:44
@ mvalidate
validation mode - do not run time steps
Definition: Constants.f90:204
@ mnormal
normal output mode
Definition: Constants.f90:205
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:26
This module contains the IdmLoadModule.
Definition: IdmLoad.f90:7
subroutine, public simnam_load(paramlog)
MODFLOW 6 mfsim.nam input load routine.
Definition: IdmLoad.f90:448
subroutine, public idm_da(iout)
idm deallocate routine
Definition: IdmLoad.f90:87
subroutine, public idm_df()
advance package dynamic data for period steps
Definition: IdmLoad.f90:39
subroutine, public simtdis_load()
MODFLOW 6 tdis input load routine.
Definition: IdmLoad.f90:470
subroutine, public idm_rp()
load package dynamic data for period
Definition: IdmLoad.f90:55
subroutine, public idm_ad()
advance package dynamic data for period steps
Definition: IdmLoad.f90:71
subroutine, public load_models(iout)
load model namfiles and model package files
Definition: IdmLoad.f90:212
subroutine, public load_exchanges(iout)
load exchange files
Definition: IdmLoad.f90:286
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
This module defines variable data types.
Definition: kind.f90:8
type(listtype), public basemodellist
Definition: mf6lists.f90:16
type(listtype), public baseexchangelist
Definition: mf6lists.f90:25
type(listtype), public solutiongrouplist
Definition: mf6lists.f90:22
type(listtype), public baseconnectionlist
Definition: mf6lists.f90:28
type(listtype), public basesolutionlist
Definition: mf6lists.f90:19
subroutine, public lists_da()
Definition: mf6lists.f90:33
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
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
Core MODFLOW 6 module.
Definition: mf6core.f90:8
subroutine mf6dotimestep()
Run time step.
Definition: mf6core.f90:584
subroutine mf6run
Main controller.
Definition: mf6core.f90:33
subroutine static_input_load()
Create simulation input context.
Definition: mf6core.f90:268
subroutine mf6preparetimestep()
Read and prepare time step.
Definition: mf6core.f90:473
subroutine simulation_df()
Define the simulation.
Definition: mf6core.f90:301
subroutine simulation_ar()
Simulation allocate and read.
Definition: mf6core.f90:373
class(runcontroltype), pointer run_ctrl
the run controller for this simulation
Definition: mf6core.f90:23
logical(lgp) function mf6update()
Run a time step.
Definition: mf6core.f90:106
subroutine sim_step_retry(finishedTrying)
Rerun time step.
Definition: mf6core.f90:626
subroutine mf6initialize()
Initialize a simulation.
Definition: mf6core.f90:70
logical(lgp) function mf6finalizetimestep()
Finalize time step.
Definition: mf6core.f90:669
subroutine create_lstfile()
Set up mfsim list file output logging.
Definition: mf6core.f90:235
subroutine connections_cr()
Create the model connections from the exchanges.
Definition: mf6core.f90:425
subroutine mf6finalize()
Finalize the simulation.
Definition: mf6core.f90:128
subroutine print_info()
print initial message
Definition: mf6core.f90:217
class(runcontroltype) function, pointer, public create_run_control()
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public converge_reset()
Reset the simulation convergence flag.
Definition: Sim.f90:392
subroutine, public initial_message()
Print the header and initializes messaging.
Definition: Sim.f90:445
subroutine, public converge_check(hasConverged)
Simulation convergence check.
Definition: Sim.f90:405
integer(i4b), parameter, public stg_aft_con_ar
afterr connection allocate read
Definition: SimStages.f90:18
integer(i4b), parameter, public stg_aft_con_df
after connection define
Definition: SimStages.f90:15
integer(i4b), parameter, public stg_aft_mdl_df
after model define
Definition: SimStages.f90:11
integer(i4b), parameter, public stg_bfr_mdl_df
before model define
Definition: SimStages.f90:10
integer(i4b), parameter, public stg_aft_exg_df
after exchange define
Definition: SimStages.f90:12
integer(i4b), parameter, public stg_aft_con_cr
after connection create
Definition: SimStages.f90:13
integer(i4b), parameter, public stg_bfr_exg_rp
before exchange read prepare
Definition: SimStages.f90:19
integer(i4b), parameter, public stg_bfr_con_ar
before connection allocate read
Definition: SimStages.f90:17
integer(i4b), parameter, public stg_aft_con_rp
after connection read prepare
Definition: SimStages.f90:20
integer(i4b), parameter, public stg_bfr_con_df
before connection define
Definition: SimStages.f90:14
subroutine, public simulation_da()
Deallocate simulation variables.
subroutine, public simulation_cr()
Read the simulation name file and initialize the models, exchanges.
This module contains simulation variables.
Definition: SimVariables.f90:9
integer(i4b) ifailedstepretry
current retry for this time step
integer(i4b) nr_procs
integer(i4b) laststepfailed
flag indicating if the last step failed (1) if last step failed; (0) otherwise (set in converge_check...
integer(i4b) iout
file unit number for simulation output
integer(i4b) iparamlog
input (idm) parameter logging to simulation listing file
character(len=linelength) simlstfile
simulation listing file name
integer(i4b) proc_id
integer(i4b) isim_mode
simulation mode
class(solutiongrouptype) function, pointer, public getsolutiongroupfromlist(list, idx)
class(spatialmodelconnectiontype) function, pointer, public get_smc_from_list(list, idx)
Get the connection from a list.
logical(lgp), pointer, public endofsimulation
flag indicating end of simulation
Definition: tdis.f90:28
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
subroutine, public tdis_da()
Deallocate memory.
Definition: tdis.f90:423
subroutine, public tdis_delt_reset(deltnew)
Reset delt and update timing variables and indicators.
Definition: tdis.f90:222
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
real(dp), pointer, public delt
length of the current time step
Definition: tdis.f90:29
subroutine, public print_start_time()
Start simulation timer.
Definition: Timer.f90:19
This module contains version information.
Definition: version.f90:7
subroutine write_listfile_header(iout, cmodel_type, write_sys_command, write_kind_info)
@ brief Write program header
Definition: version.f90:98
integer(i4b), parameter idevelopmode
Definition: version.f90:19
Highest level model type. All models extend this parent type.
Definition: BaseModel.f90:13
Class to manage spatial connection of a model to one or more models of the same type....