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 799 of file SimulationCreate.f90.

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

774  character(len=LINELENGTH) :: errmsg
775  class(BaseModelType), pointer :: mp
776  integer(I4B) :: im
777 
778  do im = 1, basemodellist%Count()
779  mp => getbasemodelfromlist(basemodellist, im)
780  if (mp%idsoln == 0) then
781  write (errmsg, '(a,a)') &
782  'Model was not assigned to a solution: ', mp%name
783  call store_error(errmsg)
784  end if
785  end do
786  if (count_errors() > 0) then
787  call store_error_filename('mfsim.nam')
788  end if
789 
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 831 of file SimulationCreate.f90.

832  ! -- dummy
833  character(len=*), intent(in) :: mtype
834  character(len=*), intent(inout) :: mname
835  ! -- local
836  integer :: ilen
837  integer :: i
838  character(len=LINELENGTH) :: errmsg
839  logical :: terminate = .true.
840 
841  ilen = len_trim(mname)
842  if (ilen > lenmodelname) then
843  write (errmsg, '(a,a)') 'Invalid model name: ', trim(mname)
844  call store_error(errmsg)
845  write (errmsg, '(a,i0,a,i0)') &
846  'Name length of ', ilen, ' exceeds maximum length of ', &
847  lenmodelname
848  call store_error(errmsg, terminate)
849  end if
850  do i = 1, ilen
851  if (mname(i:i) == ' ') then
852  write (errmsg, '(a,a)') 'Invalid model name: ', trim(mname)
853  call store_error(errmsg)
854  write (errmsg, '(a)') &
855  'Model name cannot have spaces within it.'
856  call store_error(errmsg, terminate)
857  end if
858  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
376  ! use VirtualPrtExchangeModule, only: add_virtual_prt_exchange
377  ! -- dummy
378  ! -- locals
379  character(len=LENMEMPATH) :: input_mempath
380  type(CharacterStringType), dimension(:), contiguous, &
381  pointer :: etypes !< exg types
382  type(CharacterStringType), dimension(:), contiguous, &
383  pointer :: efiles !< exg file names
384  type(CharacterStringType), dimension(:), contiguous, &
385  pointer :: emnames_a !< model a names
386  type(CharacterStringType), dimension(:), contiguous, &
387  pointer :: emnames_b !< model b names
388  type(CharacterStringType), dimension(:), contiguous, &
389  pointer :: emempaths
390  character(len=LINELENGTH) :: exgtype
391  integer(I4B) :: exg_id
392  integer(I4B) :: m1_id, m2_id
393  character(len=LINELENGTH) :: fname, name1, name2
394  character(len=LENEXCHANGENAME) :: exg_name
395  character(len=LENMEMPATH) :: exg_mempath
396  integer(I4B) :: n
397  character(len=LINELENGTH) :: errmsg
398  logical(LGP) :: terminate = .true.
399  logical(LGP) :: both_remote, both_local
400  ! -- formats
401  character(len=*), parameter :: fmtmerr = "('Error in simulation control ', &
402  &'file. Could not find model: ', a)"
403  !
404  ! -- set input memory path
405  input_mempath = create_mem_path('SIM', 'NAM', idm_context)
406  !
407  ! -- set pointers to input context exchange attribute arrays
408  call mem_setptr(etypes, 'EXGTYPE', input_mempath)
409  call mem_setptr(efiles, 'EXGFILE', input_mempath)
410  call mem_setptr(emnames_a, 'EXGMNAMEA', input_mempath)
411  call mem_setptr(emnames_b, 'EXGMNAMEB', input_mempath)
412  call mem_setptr(emempaths, 'EXGMEMPATHS', input_mempath)
413  !
414  ! -- open exchange logging block
415  write (iout, '(/1x,a)') 'READING SIMULATION EXCHANGES'
416  !
417  ! -- initialize
418  exg_id = 0
419  !
420  ! -- create exchanges
421  do n = 1, size(etypes)
422  !
423  ! -- attributes for this exchange
424  exgtype = etypes(n)
425  fname = efiles(n)
426  name1 = emnames_a(n)
427  name2 = emnames_b(n)
428  exg_mempath = emempaths(n)
429 
430  exg_id = exg_id + 1
431 
432  ! find model index in list
433  m1_id = ifind(model_names, name1)
434  if (m1_id < 0) then
435  write (errmsg, fmtmerr) trim(name1)
436  call store_error(errmsg, terminate)
437  end if
438  m2_id = ifind(model_names, name2)
439  if (m2_id < 0) then
440  write (errmsg, fmtmerr) trim(name2)
441  call store_error(errmsg, terminate)
442  end if
443 
444  ! both models on other process? then don't create it here...
445  both_remote = (model_loc_idx(m1_id) == -1 .and. &
446  model_loc_idx(m2_id) == -1)
447  both_local = (model_loc_idx(m1_id) > 0 .and. &
448  model_loc_idx(m2_id) > 0)
449  if (.not. both_remote) then
450  write (iout, '(4x,a,a,i0,a,i0,a,i0)') trim(exgtype), ' exchange ', &
451  exg_id, ' will be created to connect model ', m1_id, &
452  ' with model ', m2_id
453  end if
454 
455  select case (exgtype)
456  case ('CHF6-GWF6')
457  write (exg_name, '(a,i0)') 'CHF-GWF_', exg_id
458  if (both_local) then
459  call chfgwf_cr(fname, exg_name, exg_id, m1_id, m2_id, exg_mempath)
460  end if
461  case ('GWF6-GWF6')
462  write (exg_name, '(a,i0)') 'GWF-GWF_', exg_id
463  if (.not. both_remote) then
464  call gwfexchange_create(fname, exg_name, exg_id, m1_id, m2_id, &
465  exg_mempath)
466  end if
467  call add_virtual_gwf_exchange(exg_name, exg_id, m1_id, m2_id)
468  case ('GWF6-GWT6')
469  if (both_local) then
470  call gwfgwt_cr(fname, exg_id, m1_id, m2_id)
471  end if
472  case ('GWF6-GWE6')
473  if (both_local) then
474  call gwfgwe_cr(fname, exg_id, m1_id, m2_id)
475  end if
476  case ('GWF6-PRT6')
477  call gwfprt_cr(fname, exg_id, m1_id, m2_id)
478  case ('GWT6-GWT6')
479  write (exg_name, '(a,i0)') 'GWT-GWT_', exg_id
480  if (.not. both_remote) then
481  call gwtexchange_create(fname, exg_name, exg_id, m1_id, m2_id, &
482  exg_mempath)
483  end if
484  call add_virtual_gwt_exchange(exg_name, exg_id, m1_id, m2_id)
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_gwe_exchange(exg_name, exg_id, m1_id, m2_id)
492  case ('OLF6-GWF6')
493  write (exg_name, '(a,i0)') 'OLF-GWF_', exg_id
494  if (both_local) then
495  call olfgwf_cr(fname, exg_name, exg_id, m1_id, m2_id, exg_mempath)
496  end if
497  case default
498  write (errmsg, '(a,a)') &
499  'Unknown simulation exchange type: ', trim(exgtype)
500  call store_error(errmsg, terminate)
501  end select
502  end do
503  !
504  ! -- close exchange logging block
505  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:112
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:110
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_gwe_exchange(name, exchange_id, model1_id, model2_id)
Add a virtual GWE-GWE exchange to the simulation.
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_gwt_exchange(name, exchange_id, model1_id, model2_id)
Add a virtual GWT-GWT 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:128
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 510 of file SimulationCreate.f90.

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

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