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

Functions/Subroutines

subroutine, public simulation_cr ()
 Read the simulation name file and initialize the models, exchanges. More...
 
subroutine, public simulation_da ()
 Deallocate simulation variables. More...
 
subroutine source_simulation_nam ()
 Source the simulation name file. More...
 
subroutine options_create ()
 Set the simulation options. More...
 
subroutine timing_create ()
 Set the timing module to be used for the simulation. More...
 
subroutine models_create ()
 Set the models to be used for the simulation. More...
 
subroutine exchanges_create ()
 Set the exchanges to be used for the simulation. More...
 
subroutine solution_group_check (sgp, sgid, isgpsoln)
 Check a solution_group to be used for the simulation. More...
 
subroutine solution_groups_create ()
 Set the solution_groups to be used for the simulation. More...
 
subroutine check_model_assignment ()
 Check for dangling models, and break with error when found. More...
 
subroutine assign_exchanges ()
 Assign exchanges to solutions. More...
 
subroutine check_model_name (mtype, mname)
 Check that the model name is valid. More...
 

Function/Subroutine Documentation

◆ assign_exchanges()

subroutine simulationcreatemodule::assign_exchanges
private

This assigns NumericalExchanges to NumericalSolutions, based on the link between the models in the solution and those exchanges. The BaseExchangeconnects_model() function should be overridden to indicate if such a link exists.

Definition at line 800 of file SimulationCreate.f90.

801  ! -- local
802  class(BaseSolutionType), pointer :: sp
803  class(BaseExchangeType), pointer :: ep
804  class(BaseModelType), pointer :: mp
805  type(ListType), pointer :: models_in_solution
806  integer(I4B) :: is, ie, im
807 
808  do is = 1, basesolutionlist%Count()
809  sp => getbasesolutionfromlist(basesolutionlist, is)
810  !
811  ! -- now loop over exchanges
812  do ie = 1, baseexchangelist%Count()
813  ep => getbaseexchangefromlist(baseexchangelist, ie)
814  !
815  ! -- and add when it affects (any model in) the solution matrix
816  models_in_solution => sp%get_models()
817  do im = 1, models_in_solution%Count()
818  mp => getbasemodelfromlist(models_in_solution, im)
819  if (ep%connects_model(mp)) then
820  !
821  ! -- add to solution (and only once)
822  call sp%add_exchange(ep)
823  exit
824  end if
825  end do
826  end do
827  end do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ check_model_assignment()

subroutine simulationcreatemodule::check_model_assignment

Definition at line 774 of file SimulationCreate.f90.

775  character(len=LINELENGTH) :: errmsg
776  class(BaseModelType), pointer :: mp
777  integer(I4B) :: im
778 
779  do im = 1, basemodellist%Count()
780  mp => getbasemodelfromlist(basemodellist, im)
781  if (mp%idsoln == 0) then
782  write (errmsg, '(a,a)') &
783  'Model was not assigned to a solution: ', mp%name
784  call store_error(errmsg)
785  end if
786  end do
787  if (count_errors() > 0) then
788  call store_error_filename('mfsim.nam')
789  end if
790 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ check_model_name()

subroutine simulationcreatemodule::check_model_name ( character(len=*), intent(in)  mtype,
character(len=*), intent(inout)  mname 
)
private

Definition at line 832 of file SimulationCreate.f90.

833  ! -- dummy
834  character(len=*), intent(in) :: mtype
835  character(len=*), intent(inout) :: mname
836  ! -- local
837  integer :: ilen
838  integer :: i
839  character(len=LINELENGTH) :: errmsg
840  logical :: terminate = .true.
841 
842  ilen = len_trim(mname)
843  if (ilen > lenmodelname) then
844  write (errmsg, '(a,a)') 'Invalid model name: ', trim(mname)
845  call store_error(errmsg)
846  write (errmsg, '(a,i0,a,i0)') &
847  'Name length of ', ilen, ' exceeds maximum length of ', &
848  lenmodelname
849  call store_error(errmsg, terminate)
850  end if
851  do i = 1, ilen
852  if (mname(i:i) == ' ') then
853  write (errmsg, '(a,a)') 'Invalid model name: ', trim(mname)
854  call store_error(errmsg)
855  write (errmsg, '(a)') &
856  'Model name cannot have spaces within it.'
857  call store_error(errmsg, terminate)
858  end if
859  end do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ exchanges_create()

subroutine simulationcreatemodule::exchanges_create

Definition at line 360 of file SimulationCreate.f90.

361  ! -- modules
375  ! use VirtualPrtExchangeModule, only: add_virtual_prt_exchange
376  ! -- dummy
377  ! -- locals
378  character(len=LENMEMPATH) :: input_mempath
379  type(CharacterStringType), dimension(:), contiguous, &
380  pointer :: etypes !< exg types
381  type(CharacterStringType), dimension(:), contiguous, &
382  pointer :: efiles !< exg file names
383  type(CharacterStringType), dimension(:), contiguous, &
384  pointer :: emnames_a !< model a names
385  type(CharacterStringType), dimension(:), contiguous, &
386  pointer :: emnames_b !< model b names
387  type(CharacterStringType), dimension(:), contiguous, &
388  pointer :: emempaths
389  character(len=LINELENGTH) :: exgtype
390  integer(I4B) :: exg_id
391  integer(I4B) :: m1_id, m2_id
392  character(len=LINELENGTH) :: fname, name1, name2
393  character(len=LENEXCHANGENAME) :: exg_name
394  character(len=LENMEMPATH) :: exg_mempath
395  integer(I4B) :: n
396  character(len=LINELENGTH) :: errmsg
397  logical(LGP) :: terminate = .true.
398  logical(LGP) :: both_remote, both_local
399  ! -- formats
400  character(len=*), parameter :: fmtmerr = "('Error in simulation control ', &
401  &'file. Could not find model: ', a)"
402  !
403  ! -- set input memory path
404  input_mempath = create_mem_path('SIM', 'NAM', idm_context)
405  !
406  ! -- set pointers to input context exchange attribute arrays
407  call mem_setptr(etypes, 'EXGTYPE', input_mempath)
408  call mem_setptr(efiles, 'EXGFILE', input_mempath)
409  call mem_setptr(emnames_a, 'EXGMNAMEA', input_mempath)
410  call mem_setptr(emnames_b, 'EXGMNAMEB', input_mempath)
411  call mem_setptr(emempaths, 'EXGMEMPATHS', input_mempath)
412  !
413  ! -- open exchange logging block
414  write (iout, '(/1x,a)') 'READING SIMULATION EXCHANGES'
415  !
416  ! -- initialize
417  exg_id = 0
418  !
419  ! -- create exchanges
420  do n = 1, size(etypes)
421  !
422  ! -- attributes for this exchange
423  exgtype = etypes(n)
424  fname = efiles(n)
425  name1 = emnames_a(n)
426  name2 = emnames_b(n)
427  exg_mempath = emempaths(n)
428 
429  exg_id = exg_id + 1
430 
431  ! find model index in list
432  m1_id = ifind(model_names, name1)
433  if (m1_id < 0) then
434  write (errmsg, fmtmerr) trim(name1)
435  call store_error(errmsg, terminate)
436  end if
437  m2_id = ifind(model_names, name2)
438  if (m2_id < 0) then
439  write (errmsg, fmtmerr) trim(name2)
440  call store_error(errmsg, terminate)
441  end if
442 
443  ! both models on other process? then don't create it here...
444  both_remote = (model_loc_idx(m1_id) == -1 .and. &
445  model_loc_idx(m2_id) == -1)
446  both_local = (model_loc_idx(m1_id) > 0 .and. &
447  model_loc_idx(m2_id) > 0)
448  if (.not. both_remote) then
449  write (iout, '(4x,a,a,i0,a,i0,a,i0)') trim(exgtype), ' exchange ', &
450  exg_id, ' will be created to connect model ', m1_id, &
451  ' with model ', m2_id
452  end if
453 
454  select case (exgtype)
455  case ('CHF6-GWF6')
456  write (exg_name, '(a,i0)') 'CHF-GWF_', exg_id
457  if (both_local) then
458  call chfgwf_cr(fname, exg_name, exg_id, m1_id, m2_id, exg_mempath)
459  end if
460  case ('GWF6-GWF6')
461  write (exg_name, '(a,i0)') 'GWF-GWF_', exg_id
462  if (.not. both_remote) then
463  call gwfexchange_create(fname, exg_name, exg_id, m1_id, m2_id, &
464  exg_mempath)
465  end if
466  call add_virtual_gwf_exchange(exg_name, exg_id, m1_id, m2_id)
467  case ('GWF6-GWT6')
468  if (both_local) then
469  call gwfgwt_cr(fname, exg_id, m1_id, m2_id)
470  end if
471  case ('GWF6-GWE6')
472  if (both_local) then
473  call gwfgwe_cr(fname, exg_id, m1_id, m2_id)
474  end if
475  case ('GWF6-PRT6')
476  call gwfprt_cr(fname, exg_id, m1_id, m2_id)
477  case ('GWT6-GWT6')
478  write (exg_name, '(a,i0)') 'GWT-GWT_', exg_id
479  if (.not. both_remote) then
480  call gwtexchange_create(fname, exg_name, exg_id, m1_id, m2_id, &
481  exg_mempath)
482  end if
483  call add_virtual_tsp_exchange(exg_name, exg_id, m1_id, m2_id, &
484  "concentration")
485  case ('GWE6-GWE6')
486  write (exg_name, '(a,i0)') 'GWE-GWE_', exg_id
487  if (.not. both_remote) then
488  call gweexchange_create(fname, exg_name, exg_id, m1_id, m2_id, &
489  exg_mempath)
490  end if
491  call add_virtual_tsp_exchange(exg_name, exg_id, m1_id, m2_id, &
492  "temperature")
493  case ('OLF6-GWF6')
494  write (exg_name, '(a,i0)') 'OLF-GWF_', exg_id
495  if (both_local) then
496  call olfgwf_cr(fname, exg_name, exg_id, m1_id, m2_id, exg_mempath)
497  end if
498  case default
499  write (errmsg, '(a,a)') &
500  'Unknown simulation exchange type: ', trim(exgtype)
501  call store_error(errmsg, terminate)
502  end select
503  end do
504  !
505  ! -- close exchange logging block
506  write (iout, '(1x,a)') 'END OF SIMULATION EXCHANGES'
This module contains the ChfGwfExchangeModule Module.
Definition: exg-chfgwf.f90:6
subroutine, public chfgwf_cr(filename, name, id, m1_id, m2_id, input_mempath)
@ brief Create CHF GWF exchange
Definition: exg-chfgwf.f90:41
This module contains the GweGweExchangeModule Module.
Definition: exg-gwegwe.f90:10
subroutine, public gweexchange_create(filename, name, id, m1_id, m2_id, input_mempath)
@ brief Create GWT GWT exchange
Definition: exg-gwegwe.f90:114
subroutine, public gwfgwe_cr(filename, id, m1_id, m2_id)
Create a new GWF to GWE exchange object.
Definition: exg-gwfgwe.f90:47
This module contains the GwfGwfExchangeModule Module.
Definition: exg-gwfgwf.f90:10
subroutine, public gwfexchange_create(filename, name, id, m1_id, m2_id, input_mempath)
@ brief Create GWF GWF exchange
Definition: exg-gwfgwf.f90:122
subroutine, public gwfgwt_cr(filename, id, m1_id, m2_id)
Create a new GWF to GWT exchange object.
Definition: exg-gwfgwt.f90:49
subroutine, public gwfprt_cr(filename, id, m1id, m2id)
Create a new GWF to PRT exchange object.
Definition: exg-gwfprt.f90:40
This module contains the GwtGwtExchangeModule Module.
Definition: exg-gwtgwt.f90:10
subroutine, public gwtexchange_create(filename, name, id, m1_id, m2_id, input_mempath)
@ brief Create GWT GWT exchange
Definition: exg-gwtgwt.f90:112
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
This module contains the OlfGwfExchangeModule Module.
Definition: exg-olfgwf.f90:6
subroutine, public olfgwf_cr(filename, name, id, m1_id, m2_id, input_mempath)
@ brief Create OLF GWF exchange
Definition: exg-olfgwf.f90:41
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=linelength) idm_context
subroutine, public add_virtual_gwf_exchange(name, exchange_id, model1_id, model2_id)
Add a virtual GWF-GWF exchange to the simulation.
subroutine, public add_virtual_tsp_exchange(name, exchange_id, m1_id, m2_id, qtype)
Add a virtual GWT-GWT or GWE-GWE exchange to the simulation.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ models_create()

subroutine simulationcreatemodule::models_create

Definition at line 200 of file SimulationCreate.f90.

201  ! -- modules
206  use chfmodule, only: chf_cr
207  use gwfmodule, only: gwf_cr
208  use gwtmodule, only: gwt_cr
209  use gwemodule, only: gwe_cr
210  use olfmodule, only: olf_cr
211  use prtmodule, only: prt_cr
217  ! use VirtualPrtModelModule, only: add_virtual_prt_model
218  use constantsmodule, only: lenmodelname
219  ! -- dummy
220  ! -- locals
221  type(DistributedSimType), pointer :: ds
222  character(len=LENMEMPATH) :: input_mempath
223  type(CharacterStringType), dimension(:), contiguous, &
224  pointer :: mtypes !< model types
225  type(CharacterStringType), dimension(:), contiguous, &
226  pointer :: mfnames !< model file names
227  type(CharacterStringType), dimension(:), contiguous, &
228  pointer :: mnames !< model names
229  integer(I4B) :: im
230  class(NumericalModelType), pointer :: num_model
231  class(ExplicitModelType), pointer :: exp_model
232  character(len=LINELENGTH) :: model_type
233  character(len=LINELENGTH) :: fname, model_name
234  integer(I4B) :: n, nr_models_glob
235  integer(I4B), dimension(:), pointer :: model_ranks => null()
236  logical :: terminate = .true.
237  !
238  ! -- set input memory path
239  input_mempath = create_mem_path('SIM', 'NAM', idm_context)
240  !
241  ! -- set pointers to input context model attribute arrays
242  call mem_setptr(mtypes, 'MTYPE', input_mempath)
243  call mem_setptr(mfnames, 'MFNAME', input_mempath)
244  call mem_setptr(mnames, 'MNAME', input_mempath)
245  !
246  ! -- allocate global arrays
247  nr_models_glob = size(mnames)
248  allocate (model_names(nr_models_glob))
249  allocate (model_loc_idx(nr_models_glob))
250  !
251  ! -- get model-to-cpu assignment (in serial all to rank 0)
252  ds => get_dsim()
253  model_ranks => ds%get_load_balance()
254  !
255  ! -- open model logging block
256  write (iout, '(/1x,a)') 'READING SIMULATION MODELS'
257  !
258  ! -- create models
259  im = 0
260  do n = 1, size(mtypes)
261  !
262  ! -- attributes for this model
263  model_type = mtypes(n)
264  fname = mfnames(n)
265  model_name = mnames(n)
266  !
267  call check_model_name(model_type, model_name)
268  !
269  ! increment global model id
270  model_names(n) = model_name(1:lenmodelname)
271  model_loc_idx(n) = -1
272  num_model => null()
273  exp_model => null()
274  !
275  ! -- add a new (local or global) model
276  select case (model_type)
277  case ('GWF6')
278  if (model_ranks(n) == proc_id) then
279  im = im + 1
280  write (iout, '(4x,2a,i0,a)') trim(model_type), ' model ', &
281  n, ' will be created'
282  call gwf_cr(fname, n, model_names(n))
283  num_model => getnumericalmodelfromlist(basemodellist, im)
284  model_loc_idx(n) = im
285  end if
286  call add_virtual_gwf_model(n, model_names(n), num_model)
287  case ('GWT6')
288  if (model_ranks(n) == proc_id) then
289  im = im + 1
290  write (iout, '(4x,2a,i0,a)') trim(model_type), ' model ', &
291  n, ' will be created'
292  call gwt_cr(fname, n, model_names(n))
293  num_model => getnumericalmodelfromlist(basemodellist, im)
294  model_loc_idx(n) = im
295  end if
296  call add_virtual_gwt_model(n, model_names(n), num_model)
297  case ('GWE6')
298  if (model_ranks(n) == proc_id) then
299  im = im + 1
300  write (iout, '(4x,2a,i0,a)') trim(model_type), ' model ', &
301  n, ' will be created'
302  call gwe_cr(fname, n, model_names(n))
303  num_model => getnumericalmodelfromlist(basemodellist, im)
304  model_loc_idx(n) = im
305  end if
306  call add_virtual_gwe_model(n, model_names(n), num_model)
307  case ('CHF6')
308  if (model_ranks(n) == proc_id) then
309  im = im + 1
310  write (iout, '(4x,2a,i0,a)') trim(model_type), " model ", &
311  n, " will be created"
312  call chf_cr(fname, n, model_names(n))
313  call developmode('CHF is still under development, install the &
314  &nightly build or compile from source with IDEVELOPMODE = 1.')
315  num_model => getnumericalmodelfromlist(basemodellist, im)
316  model_loc_idx(n) = im
317  end if
318  case ('OLF6')
319  if (model_ranks(n) == proc_id) then
320  im = im + 1
321  write (iout, '(4x,2a,i0,a)') trim(model_type), " model ", &
322  n, " will be created"
323  call olf_cr(fname, n, model_names(n))
324  call developmode('OLF is still under development, install the &
325  &nightly build or compile from source with IDEVELOPMODE = 1.')
326  num_model => getnumericalmodelfromlist(basemodellist, im)
327  model_loc_idx(n) = im
328  end if
329  case ('PRT6')
330  if (model_ranks(n) == proc_id) then
331  im = im + 1
332  write (iout, '(4x,2a,i0,a)') trim(model_type), ' model ', &
333  n, ' will be created'
334  call prt_cr(fname, n, model_names(n))
335  exp_model => getexplicitmodelfromlist(basemodellist, im)
336  model_loc_idx(n) = im
337  end if
338  ! When virtual PRT is implemented, uncomment:
339  ! call add_virtual_prt_model(n, model_names(n), exp_model)
340  case default
341  write (errmsg, '(a,a)') &
342  'Unknown simulation model type: ', trim(model_type)
343  call store_error(errmsg, terminate)
344  end select
345  end do
346  !
347  ! -- close model logging block
348  write (iout, '(1x,a)') 'END OF SIMULATION MODELS'
349  !
350  ! -- sanity check
351  if (simulation_mode == 'PARALLEL' .and. im == 0) then
352  write (errmsg, '(a, i0)') &
353  'No MODELS assigned to process ', proc_id
354  call store_error(errmsg, terminate)
355  end if
Channel Flow (CHF) Module.
Definition: chf.f90:3
subroutine, public chf_cr(filename, id, modelname)
Create a new surface water flow model object.
Definition: chf.f90:56
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter lenmodelname
maximum length of the model name
Definition: Constants.f90:22
class(distributedsimtype) function, pointer, public get_dsim()
Get pointer to the distributed simulation object.
Models that solve themselves.
class(explicitmodeltype) function, pointer, public getexplicitmodelfromlist(list, idx)
@ brief Get generic object from list and return as explicit model
Definition: gwe.f90:3
subroutine, public gwe_cr(filename, id, modelname)
Create a new groundwater energy transport model object.
Definition: gwe.f90:98
Definition: gwf.f90:1
subroutine, public gwf_cr(filename, id, modelname)
Create a new groundwater flow model object.
Definition: gwf.f90:138
Definition: gwt.f90:8
subroutine, public gwt_cr(filename, id, modelname)
Create a new groundwater transport model object.
Definition: gwt.f90:101
class(numericalmodeltype) function, pointer, public getnumericalmodelfromlist(list, idx)
Channel Flow (OLF) Module.
Definition: olf.f90:3
subroutine, public olf_cr(filename, id, modelname)
Create a new overland flow model object.
Definition: olf.f90:56
Definition: prt.f90:1
subroutine, public prt_cr(filename, id, modelname)
Create a new particle tracking model object.
Definition: prt.f90:133
character(len=maxcharlen) errmsg
error message string
subroutine, public add_virtual_gwe_model(model_id, model_name, model)
subroutine, public add_virtual_gwf_model(model_id, model_name, model)
Add virtual GWF model.
subroutine, public add_virtual_gwt_model(model_id, model_name, model)
Base type for models that solve themselves.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ options_create()

subroutine simulationcreatemodule::options_create
private

Definition at line 95 of file SimulationCreate.f90.

96  ! -- modules
102  ! -- dummy
103  ! -- locals
104  character(len=LENMEMPATH) :: input_mempath
105  integer(I4B), pointer :: simcontinue, nocheck, maxerror
106  character(len=:), pointer :: prmem
107  character(len=LINELENGTH) :: errmsg
108  !
109  ! -- set input memory path
110  input_mempath = create_mem_path('SIM', 'NAM', idm_context)
111  !
112  ! -- set pointers to input context option params
113  call mem_setptr(simcontinue, 'CONTINUE', input_mempath)
114  call mem_setptr(nocheck, 'NOCHECK', input_mempath)
115  call mem_setptr(prmem, 'PRMEM', input_mempath)
116  call mem_setptr(maxerror, 'MAXERRORS', input_mempath)
117  !
118  ! -- update sim options
119  isimcontinue = simcontinue
120  if (nocheck == 1) then
121  isimcheck = 0
122  else
123  isimcheck = 1
124  end if
125 
126  call maxerrors(maxerror)
127 
128  if (prmem /= '') then
129  errmsg = ''
130  call mem_set_print_option(iout, prmem, errmsg)
131  if (errmsg /= '') then
132  call store_error(errmsg, .true.)
133  end if
134  end if
135 
136  !
137  ! -- log values to list file
138  if (iout > 0) then
139  write (iout, '(/1x,a)') 'READING SIMULATION OPTIONS'
140  !
141  if (isimcontinue == 1) then
142  write (iout, '(4x, a)') &
143  'SIMULATION WILL CONTINUE EVEN IF THERE IS NONCONVERGENCE.'
144  end if
145  !
146  if (isimcheck == 0) then
147  write (iout, '(4x, a)') &
148  'MODEL DATA WILL NOT BE CHECKED FOR ERRORS.'
149  end if
150  !
151  write (iout, '(4x, a, i0)') &
152  'MAXIMUM NUMBER OF ERRORS THAT WILL BE STORED IS ', maxerror
153  !
154  if (prmem /= '') then
155  write (iout, '(4x, a, a, a)') &
156  'MEMORY_PRINT_OPTION SET TO "', trim(prmem), '".'
157  end if
158  !
159  write (iout, '(1x,a)') 'END OF SIMULATION OPTIONS'
160  end if
subroutine, public mem_set_print_option(iout, keyword, error_msg)
Set the memory print option.
integer(i4b) isimcheck
simulation input check flag (1) to check input, (0) to ignore checks
integer(i4b) isimcontinue
simulation continue flag (1) to continue if isimcnvg = 0, (0) to terminate
Here is the call graph for this function:
Here is the caller graph for this function:

◆ simulation_cr()

subroutine, public simulationcreatemodule::simulation_cr

Definition at line 34 of file SimulationCreate.f90.

35  ! -- modules
36  ! -- local
37  !
38  ! -- Source simulation nam input context and create objects
39  call source_simulation_nam()
Here is the call graph for this function:
Here is the caller graph for this function:

◆ simulation_da()

subroutine, public simulationcreatemodule::simulation_da

Definition at line 44 of file SimulationCreate.f90.

45  ! -- modules
48  ! -- local
49  type(DistributedSimType), pointer :: ds
50 
51  ! -- variables
52  ds => get_dsim()
53  call ds%destroy()
54  !
55  deallocate (model_names)
56  deallocate (model_loc_idx)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ solution_group_check()

subroutine simulationcreatemodule::solution_group_check ( type(solutiongrouptype), intent(inout), pointer  sgp,
integer(i4b), intent(in)  sgid,
integer(i4b), intent(in)  isgpsoln 
)

Definition at line 511 of file SimulationCreate.f90.

512  ! -- modules
513  ! -- dummy
514  type(SolutionGroupType), pointer, intent(inout) :: sgp
515  integer(I4B), intent(in) :: sgid
516  integer(I4B), intent(in) :: isgpsoln
517  ! -- local
518  character(len=LINELENGTH) :: errmsg
519  logical :: terminate = .true.
520  ! -- formats
521  character(len=*), parameter :: fmterrmxiter = &
522  "('MXITER is set to ', i0, ' but there is only one solution', &
523  &' in SOLUTION GROUP ', i0, '. Set MXITER to 1 in simulation control', &
524  &' file.')"
525  !
526  ! -- error check completed group
527  if (sgid > 0) then
528  !
529  ! -- Make sure there is a solution in this solution group
530  if (isgpsoln == 0) then
531  write (errmsg, '(a,i0)') &
532  'There are no solutions for solution group ', sgid
533  call store_error(errmsg, terminate)
534  end if
535  !
536  ! -- If there is only one solution then mxiter should be 1.
537  if (isgpsoln == 1 .and. sgp%mxiter > 1) then
538  write (errmsg, fmterrmxiter) sgp%mxiter, isgpsoln
539  call store_error(errmsg, terminate)
540  end if
541  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ solution_groups_create()

subroutine simulationcreatemodule::solution_groups_create
private

Definition at line 546 of file SimulationCreate.f90.

547  ! -- modules
555  use basemodelmodule, only: basemodeltype
560  ! -- dummy
561  ! -- local
562  character(len=LENMEMPATH) :: input_mempath
563  type(CharacterStringType), dimension(:), contiguous, &
564  pointer :: slntype
565  type(CharacterStringType), dimension(:), contiguous, &
566  pointer :: slnfname
567  type(CharacterStringType), dimension(:), contiguous, &
568  pointer :: slnmnames
569  integer(I4B), dimension(:), contiguous, pointer :: blocknum
570  character(len=LINELENGTH) :: stype, fname
571  character(len=:), allocatable :: mnames
572  type(SolutionGroupType), pointer :: sgp
573  class(BaseSolutionType), pointer :: sp
574  class(BaseModelType), pointer :: mp
575  integer(I4B) :: isoln
576  integer(I4B) :: isgpsoln
577  integer(I4B) :: sgid
578  integer(I4B) :: glo_mid
579  integer(I4B) :: loc_idx
580  integer(I4B) :: i, j, istat, mxiter
581  integer(I4B) :: nwords
582  character(len=LENMODELNAME), dimension(:), allocatable :: words
583  character(len=:), allocatable :: parse_str
584  character(len=LINELENGTH) :: errmsg
585  logical :: terminate = .true.
586  !
587  ! -- set memory path
588  input_mempath = create_mem_path('SIM', 'NAM', idm_context)
589  !
590  ! -- set pointers to input context solution attribute arrays
591  call mem_setptr(slntype, 'SLNTYPE', input_mempath)
592  call mem_setptr(slnfname, 'SLNFNAME', input_mempath)
593  call mem_setptr(slnmnames, 'SLNMNAMES', input_mempath)
594  call mem_setptr(blocknum, 'SOLUTIONGROUPNUM', input_mempath)
595  !
596  ! -- open solution group logging block
597  write (iout, '(/1x,a)') 'READING SOLUTIONGROUP'
598  !
599  ! -- initialize
600  sgid = 0 ! integer id of soln group, tracks with blocknum
601  isoln = 0 ! cumulative solution number
602  !
603  ! -- create solution groups
604  do i = 1, size(blocknum)
605  !
606  ! -- allocate slnmnames string
607  allocate (character(slnmnames(i)%strlen()) :: mnames)
608  !
609  ! -- attributes for this solution
610  stype = slntype(i)
611  fname = slnfname(i)
612  mnames = slnmnames(i)
613 
614  if (blocknum(i) /= sgid) then
615  !
616  ! -- check for new soln group
617  if (blocknum(i) == sgid + 1) then
618  !
619  ! -- error check completed group
620  call solution_group_check(sgp, sgid, isgpsoln)
621  !
622  ! -- reinitialize
623  nullify (sgp)
624  isgpsoln = 0 ! solution counter for this solution group
625  !
626  ! -- set sgid
627  sgid = blocknum(i)
628  !
629  ! -- create new soln group and add to global list
630  call solutiongroup_create(sgp, sgid)
631  call addsolutiongrouptolist(solutiongrouplist, sgp)
632  else
633  write (errmsg, '(a,i0,a,i0,a)') &
634  'Solution group blocks are not listed consecutively. Found ', &
635  blocknum(i), ' when looking for ', sgid + 1, '.'
636  call store_error(errmsg, terminate)
637  end if
638  end if
639  !
640  ! --
641  select case (stype)
642  !
643  case ('MXITER')
644  read (fname, *, iostat=istat) mxiter
645  if (istat == 0) then
646  sgp%mxiter = mxiter
647  end if
648  case ('IMS6')
649  !
650  ! -- increment solution counters
651  isoln = isoln + 1
652  isgpsoln = isgpsoln + 1
653  !
654  ! -- create soln and add to group
655  sp => create_ims_solution(simulation_mode, fname, isoln)
656  call sgp%add_solution(isoln, sp)
657  !
658  ! -- parse model names
659  parse_str = trim(mnames)//' '
660  call parseline(parse_str, nwords, words)
661  !
662  ! -- Find each model id and get model
663  do j = 1, nwords
664  call upcase(words(j))
665  glo_mid = ifind(model_names, words(j))
666  if (glo_mid == -1) then
667  write (errmsg, '(a,a)') 'Invalid model name: ', trim(words(j))
668  call store_error(errmsg, terminate)
669  end if
670  !
671  loc_idx = model_loc_idx(glo_mid)
672  if (loc_idx == -1) then
673  if (simulation_mode == 'PARALLEL') then
674  ! this is still ok
675  cycle
676  end if
677  end if
678  !
679  mp => getbasemodelfromlist(basemodellist, loc_idx)
680  !
681  ! -- Check for solver-model type mismatch
682  select type (mp)
683  class is (explicitmodeltype)
684  write (errmsg, '(4a)') &
685  'Model "', trim(words(j)), &
686  '" is an explicit model and cannot be added to an IMS6 ', &
687  'solution. Explicit models require EMS6.'
688  call store_error(errmsg)
689  end select
690  !
691  ! -- Add the model to the solution
692  call sp%add_model(mp)
693  mp%idsoln = isoln
694  end do
695  !
696  ! -- Check for any errors
697  if (count_errors() > 0) then
698  call store_error_filename('mfsim.nam')
699  end if
700  case ('EMS6')
701  !
702  ! -- increment solution counters
703  isoln = isoln + 1
704  isgpsoln = isgpsoln + 1
705  !
706  ! -- create soln and add to group
707  sp => create_ems_solution(simulation_mode, fname, isoln)
708  call sgp%add_solution(isoln, sp)
709  !
710  ! -- parse model names
711  parse_str = trim(mnames)//' '
712  call parseline(parse_str, nwords, words)
713  !
714  ! -- Find each model id and get model
715  do j = 1, nwords
716  call upcase(words(j))
717  glo_mid = ifind(model_names, words(j))
718  if (glo_mid == -1) then
719  write (errmsg, '(a,a)') 'Invalid model name: ', trim(words(j))
720  call store_error(errmsg, terminate)
721  end if
722  !
723  loc_idx = model_loc_idx(glo_mid)
724  if (loc_idx == -1) then
725  if (simulation_mode == 'PARALLEL') then
726  ! this is still ok
727  cycle
728  end if
729  end if
730  !
731  mp => getbasemodelfromlist(basemodellist, loc_idx)
732  !
733  ! -- Check for solver-model type mismatch
734  select type (mp)
735  class is (numericalmodeltype)
736  write (errmsg, '(4a)') &
737  'Model "', trim(words(j)), &
738  '" is a numerical model and cannot be added to an EMS6 ', &
739  'solution. Numerical models require IMS6.'
740  call store_error(errmsg)
741  end select
742  !
743  ! -- Add the model to the solution
744  call sp%add_model(mp)
745  mp%idsoln = isoln
746  end do
747  !
748  ! -- Check for any errors
749  if (count_errors() > 0) then
750  call store_error_filename('mfsim.nam')
751  end if
752  case default
753  end select
754  !
755  ! -- clean up
756  deallocate (mnames)
757  end do
758  !
759  ! -- error check final group
760  call solution_group_check(sgp, sgid, isgpsoln)
761  !
762  ! -- close exchange logging block
763  write (iout, '(1x,a)') 'END OF SOLUTIONGROUP'
764  !
765  ! -- Check and make sure at least one solution group was found
766  if (solutiongrouplist%Count() == 0) then
767  call store_error('There are no solution groups.', terminate)
768  end if
subroutine, public parseline(line, nwords, words, inunit, filename)
Parse a line into words.
subroutine, public upcase(word)
Convert to upper case.
character(len=linelength) simulation_mode
class(basesolutiontype) function, pointer, public create_ims_solution(sim_mode, filename, sol_id)
Create an IMS solution of type NumericalSolution for serial runs or its sub-type ParallelSolution for...
class(basesolutiontype) function, pointer, public create_ems_solution(sim_mode, filename, sol_id)
Create an EMS solution of type ExplicitSolution for serial runs or its sub-type ParallelSolution for.
subroutine, public solutiongroup_create(sgp, id)
Create a new solution group.
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:

◆ source_simulation_nam()

subroutine simulationcreatemodule::source_simulation_nam

Source from the simulation nam input context to initialize the models, exchanges, solutions, solutions groups. Then add the exchanges to the appropriate solutions.

Definition at line 66 of file SimulationCreate.f90.

67  ! -- dummy
68  ! -- local
69  !
70  ! -- Process OPTIONS block in namfile
71  call options_create()
72  !
73  ! -- Process TIMING block in namfile
74  call timing_create()
75  !
76  ! -- Process MODELS block in namfile
77  call models_create()
78  !
79  ! -- Process EXCHANGES block in namfile
80  call exchanges_create()
81  !
82  ! -- Process SOLUTION_GROUPS blocks in namfile
83  call solution_groups_create()
84  !
85  ! -- Go through each model and make sure that it has been assigned to
86  ! a solution.
87  call check_model_assignment()
88  !
89  ! -- Go through each solution and assign exchanges accordingly
90  call assign_exchanges()
Here is the call graph for this function:
Here is the caller graph for this function:

◆ timing_create()

subroutine simulationcreatemodule::timing_create

Definition at line 165 of file SimulationCreate.f90.

166  ! -- modules
170  use tdismodule, only: tdis_cr
171  ! -- dummy
172  ! -- locals
173  character(len=LENMEMPATH) :: input_mempath
174  character(len=LENMEMPATH) :: tdis_input_mempath
175  character(len=:), pointer :: tdis6
176  logical :: terminate = .true.
177  !
178  ! -- set input memory path
179  input_mempath = create_mem_path('SIM', 'NAM', idm_context)
180  tdis_input_mempath = create_mem_path('SIM', 'TDIS', idm_context)
181  !
182  write (iout, '(/1x,a)') 'READING SIMULATION TIMING'
183  !
184  ! -- set pointers to input context timing params
185  call mem_setptr(tdis6, 'TDIS6', input_mempath)
186  !
187  ! -- create timing
188  if (tdis6 /= '') then
189  call tdis_cr(tdis6, tdis_input_mempath)
190  else
191  call store_error('TIMING block variable TDIS6 is unset'// &
192  ' in simulation control input.', terminate)
193  end if
194  !
195  write (iout, '(1x,a)') 'END OF SIMULATION TIMING'
subroutine, public tdis_cr(fname, inmempath)
Create temporal discretization.
Definition: tdis.f90:55
Here is the call graph for this function:
Here is the caller graph for this function: