MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
exg-gwfgwe.f90
Go to the documentation of this file.
2  use kindmodule, only: dp, i4b, lgp
6  use simmodule, only: store_error
7  use simvariablesmodule, only: errmsg
16  use gwfmodule, only: gwfmodeltype
17  use gwemodule, only: gwemodeltype
18  use bndmodule, only: bndtype, getbndfromlist
19 
20  implicit none
21  public :: gwfgweexchangetype
22  public :: gwfgwe_cr
23 
25 
26  integer(I4B), pointer :: m1_idx => null() !< index into the list of base exchanges for model 1
27  integer(I4B), pointer :: m2_idx => null() !< index into the list of base exchanges for model 2
28 
29  contains
30 
31  procedure :: exg_df
32  procedure :: exg_ar
33  procedure :: exg_da
34  procedure, private :: set_model_pointers
35  procedure, private :: allocate_scalars
36  procedure, private :: gwfbnd2gwefmi
37  procedure, private :: gwfconn2gweconn
38  procedure, private :: link_connections
39 
40  end type gwfgweexchangetype
41 
42 contains
43 
44  !> @brief Create a new GWF to GWE exchange object
45  !<
46  subroutine gwfgwe_cr(filename, id, m1_id, m2_id)
47  ! -- modules
49  ! -- dummy
50  character(len=*), intent(in) :: filename
51  integer(I4B), intent(in) :: id
52  integer(I4B), intent(in) :: m1_id
53  integer(I4B), intent(in) :: m2_id
54  ! -- local
55  class(baseexchangetype), pointer :: baseexchange => null()
56  type(gwfgweexchangetype), pointer :: exchange => null()
57  character(len=20) :: cint
58  !
59  ! -- Create a new exchange and add it to the baseexchangelist container
60  allocate (exchange)
61  baseexchange => exchange
62  call addbaseexchangetolist(baseexchangelist, baseexchange)
63  !
64  ! -- Assign id and name
65  exchange%id = id
66  write (cint, '(i0)') id
67  exchange%name = 'GWF-GWE_'//trim(adjustl(cint))
68  exchange%memoryPath = exchange%name
69  !
70  ! -- allocate scalars
71  call exchange%allocate_scalars()
72  !
73  ! -- NB: convert from id to local model index in base model list
74  exchange%m1_idx = model_loc_idx(m1_id)
75  exchange%m2_idx = model_loc_idx(m2_id)
76  !
77  ! -- set model pointers
78  call exchange%set_model_pointers()
79  end subroutine gwfgwe_cr
80 
81  !> @brief Allocate and read
82  !<
83  subroutine set_model_pointers(this)
84  ! -- dummy
85  class(gwfgweexchangetype) :: this
86  ! -- local
87  class(basemodeltype), pointer :: mb => null()
88  type(gwfmodeltype), pointer :: gwfmodel => null()
89  type(gwemodeltype), pointer :: gwemodel => null()
90  !
91  ! -- set gwfmodel
92  gwfmodel => null()
93  mb => getbasemodelfromlist(basemodellist, this%m1_idx)
94  select type (mb)
95  type is (gwfmodeltype)
96  gwfmodel => mb
97  end select
98  !
99  ! -- set gwemodel
100  gwemodel => null()
101  mb => getbasemodelfromlist(basemodellist, this%m2_idx)
102  select type (mb)
103  type is (gwemodeltype)
104  gwemodel => mb
105  end select
106  !
107  ! -- Verify that gwf model is of the correct type
108  if (.not. associated(gwfmodel)) then
109  write (errmsg, '(3a)') 'Problem with GWF-GWE exchange ', trim(this%name), &
110  '. Specified GWF Model does not appear to be of the correct type.'
111  call store_error(errmsg, terminate=.true.)
112  end if
113  !
114  ! -- Verify that gwe model is of the correct type
115  if (.not. associated(gwemodel)) then
116  write (errmsg, '(3a)') 'Problem with GWF-GWE exchange ', trim(this%name), &
117  '. Specified GWF Model does not appear to be of the correct type.'
118  call store_error(errmsg, terminate=.true.)
119  end if
120  !
121  ! -- Tell transport model fmi flows are not read from file
122  gwemodel%fmi%flows_from_file = .false.
123  !
124  ! -- Set a pointer to the GWF bndlist. This will allow the transport model
125  ! to look through the flow packages and establish a link to GWF flows
126  gwemodel%fmi%gwfbndlist => gwfmodel%bndlist
127  end subroutine set_model_pointers
128 
129  !> @brief Define the GwfGwe Exchange object
130  !<
131  subroutine exg_df(this)
132  ! -- modules
134  ! -- dummy
135  class(gwfgweexchangetype) :: this
136  ! -- local
137  class(basemodeltype), pointer :: mb => null()
138  type(gwfmodeltype), pointer :: gwfmodel => null()
139  type(gwemodeltype), pointer :: gwemodel => null()
140  !
141  ! -- set gwfmodel
142  mb => getbasemodelfromlist(basemodellist, this%m1_idx)
143  select type (mb)
144  type is (gwfmodeltype)
145  gwfmodel => mb
146  end select
147  !
148  ! -- set gwemodel
149  mb => getbasemodelfromlist(basemodellist, this%m2_idx)
150  select type (mb)
151  type is (gwemodeltype)
152  gwemodel => mb
153  end select
154  !
155  ! -- Check to make sure that flow is solved before transport and in a
156  ! different IMS solution
157  if (gwfmodel%idsoln >= gwemodel%idsoln) then
158  write (errmsg, '(3a)') 'Problem with GWF-GWE exchange ', trim(this%name), &
159  '. The GWF model must be solved by a different IMS than the GWE model. &
160  &Furthermore, the IMS specified for GWF must be listed in mfsim.nam &
161  &before the IMS for GWE.'
162  call store_error(errmsg, terminate=.true.)
163  end if
164  !
165  ! -- Set pointer to flowja
166  gwemodel%fmi%gwfflowja => gwfmodel%flowja
167  call mem_checkin(gwemodel%fmi%gwfflowja, &
168  'GWFFLOWJA', gwemodel%fmi%memoryPath, &
169  'FLOWJA', gwfmodel%memoryPath)
170 
171  !
172  ! -- Set the npf flag so that specific discharge is available for
173  ! transport calculations if dispersion is active
174  if (gwemodel%incnd > 0) then
175  gwfmodel%npf%icalcspdis = 1
176  end if
177  end subroutine exg_df
178 
179  !> @brief Allocate and read
180  !<
181  subroutine exg_ar(this)
182  ! -- modules
184  use dismodule, only: distype
185  use disvmodule, only: disvtype
186  use disumodule, only: disutype
187  ! -- dummy
188  class(gwfgweexchangetype) :: this
189  ! -- local
190  class(basemodeltype), pointer :: mb => null()
191  type(gwfmodeltype), pointer :: gwfmodel => null()
192  type(gwemodeltype), pointer :: gwemodel => null()
193  ! -- formats
194  character(len=*), parameter :: fmtdiserr = &
195  "('GWF and GWE Models do not have the same discretization for exchange&
196  & ',a,'.&
197  & GWF Model has ', i0, ' user nodes and ', i0, ' reduced nodes.&
198  & GWE Model has ', i0, ' user nodes and ', i0, ' reduced nodes.&
199  & Ensure discretization packages, including IDOMAIN, are identical.')"
200  character(len=*), parameter :: fmtidomerr = &
201  "('GWF and GWE Models do not have the same discretization for exchange&
202  & ',a,'.&
203  & GWF Model and GWE Model have different IDOMAIN arrays.&
204  & Ensure discretization packages, including IDOMAIN, are identical.')"
205  !
206  ! -- set gwfmodel
207  mb => getbasemodelfromlist(basemodellist, this%m1_idx)
208  select type (mb)
209  type is (gwfmodeltype)
210  gwfmodel => mb
211  end select
212  !
213  ! -- set gwemodel
214  mb => getbasemodelfromlist(basemodellist, this%m2_idx)
215  select type (mb)
216  type is (gwemodeltype)
217  gwemodel => mb
218  end select
219  !
220  ! -- Check to make sure sizes are identical
221  if (gwemodel%dis%nodes /= gwfmodel%dis%nodes .or. &
222  gwemodel%dis%nodesuser /= gwfmodel%dis%nodesuser) then
223  write (errmsg, fmtdiserr) trim(this%name), &
224  gwfmodel%dis%nodesuser, &
225  gwfmodel%dis%nodes, &
226  gwemodel%dis%nodesuser, &
227  gwemodel%dis%nodes
228  call store_error(errmsg, terminate=.true.)
229  end if
230  !
231  ! -- Make sure idomains are identical
232  select type (gwfdis => gwfmodel%dis)
233  type is (distype)
234  select type (gwedis => gwemodel%dis)
235  type is (distype)
236  if (.not. all(gwfdis%idomain == gwedis%idomain)) then
237  write (errmsg, fmtidomerr) trim(this%name)
238  call store_error(errmsg, terminate=.true.)
239  end if
240  end select
241  type is (disvtype)
242  select type (gwedis => gwemodel%dis)
243  type is (disvtype)
244  if (.not. all(gwfdis%idomain == gwedis%idomain)) then
245  write (errmsg, fmtidomerr) trim(this%name)
246  call store_error(errmsg, terminate=.true.)
247  end if
248  end select
249  type is (disutype)
250  select type (gwedis => gwemodel%dis)
251  type is (disutype)
252  if (.not. all(gwfdis%idomain == gwedis%idomain)) then
253  write (errmsg, fmtidomerr) trim(this%name)
254  call store_error(errmsg, terminate=.true.)
255  end if
256  end select
257  end select
258  !
259  ! -- setup pointers to gwf variables allocated in gwf_ar
260  gwemodel%fmi%gwfhead => gwfmodel%x
261  call mem_checkin(gwemodel%fmi%gwfhead, &
262  'GWFHEAD', gwemodel%fmi%memoryPath, &
263  'X', gwfmodel%memoryPath)
264  gwemodel%fmi%gwfsat => gwfmodel%npf%sat
265  call mem_checkin(gwemodel%fmi%gwfsat, &
266  'GWFSAT', gwemodel%fmi%memoryPath, &
267  'SAT', gwfmodel%npf%memoryPath)
268  gwemodel%fmi%gwfspdis => gwfmodel%npf%spdis
269  call mem_checkin(gwemodel%fmi%gwfspdis, &
270  'GWFSPDIS', gwemodel%fmi%memoryPath, &
271  'SPDIS', gwfmodel%npf%memoryPath)
272  !
273  ! -- setup pointers to the flow storage rates. GWF strg arrays are
274  ! available after the gwf_ar routine is called.
275  if (gwemodel%inest > 0) then
276  if (gwfmodel%insto > 0) then
277  gwemodel%fmi%gwfstrgss => gwfmodel%sto%strgss
278  gwemodel%fmi%igwfstrgss = 1
279  if (gwfmodel%sto%iusesy == 1) then
280  gwemodel%fmi%gwfstrgsy => gwfmodel%sto%strgsy
281  gwemodel%fmi%igwfstrgsy = 1
282  end if
283  end if
284  end if
285  !
286  ! -- Set a pointer to conc in buy
287  if (gwfmodel%inbuy > 0) then
288  call gwfmodel%buy%set_concentration_pointer(gwemodel%name, gwemodel%x, &
289  gwemodel%ibound)
290  end if
291  !
292  ! -- Set a pointer to conc (which could be a temperature) in vsc
293  if (gwfmodel%invsc > 0) then
294  call gwfmodel%vsc%set_concentration_pointer(gwemodel%name, gwemodel%x, &
295  gwemodel%ibound, 1)
296  end if
297  !
298  ! -- transfer the boundary package information from gwf to gwe
299  call this%gwfbnd2gwefmi()
300  !
301  ! -- if mover package is active, then set a pointer to it's budget object
302  if (gwfmodel%inmvr /= 0) then
303  gwemodel%fmi%mvrbudobj => gwfmodel%mvr%budobj
304  end if
305  !
306  ! -- connect Connections
307  call this%gwfconn2gweconn(gwfmodel, gwemodel)
308  end subroutine exg_ar
309 
310  !> @brief Link GWE connections to GWF connections or exchanges
311  !<
312  subroutine gwfconn2gweconn(this, gwfModel, gweModel)
313  ! -- modules
314  use simmodule, only: store_error
315  use simvariablesmodule, only: iout
317  ! -- dummy
318  class(gwfgweexchangetype) :: this !< this exchange
319  type(gwfmodeltype), pointer :: gwfModel !< the flow model
320  type(gwemodeltype), pointer :: gweModel !< the energy transport model
321  ! -- local
322  class(spatialmodelconnectiontype), pointer :: conn => null()
323  class(*), pointer :: objPtr => null()
324  class(gwegweconnectiontype), pointer :: gweConn => null()
325  class(gwfgwfconnectiontype), pointer :: gwfConn => null()
326  class(gwfexchangetype), pointer :: gwfEx => null()
327  integer(I4B) :: ic1, ic2, iex
328  integer(I4B) :: gwfConnIdx, gwfExIdx
329  logical(LGP) :: areEqual
330  !
331  ! loop over all connections
332  gweloop: do ic1 = 1, baseconnectionlist%Count()
333  !
335  if (.not. associated(conn%owner, gwemodel)) cycle gweloop
336  !
337  ! start with a GWE conn.
338  objptr => conn
339  gweconn => castasgwegweconnection(objptr)
340  gwfconnidx = -1
341  gwfexidx = -1
342  !
343  ! find matching GWF conn. in same list
344  gwfloop: do ic2 = 1, baseconnectionlist%Count()
346  !
347  if (associated(conn%owner, gwfmodel)) then
348  objptr => conn
349  gwfconn => castasgwfgwfconnection(objptr)
350  !
351  ! for now, connecting the same nodes nrs will be
352  ! sufficient evidence of equality
353  areequal = all(gwfconn%prim_exchange%nodem1 == &
354  gweconn%prim_exchange%nodem1)
355  areequal = areequal .and. all(gwfconn%prim_exchange%nodem2 == &
356  gweconn%prim_exchange%nodem2)
357  if (areequal) then
358  ! same DIS, same exchange: link and go to next GWE conn.
359  write (iout, '(/6a)') 'Linking exchange ', &
360  trim(gweconn%prim_exchange%name), &
361  ' to ', trim(gwfconn%prim_exchange%name), &
362  ' (using interface model) for GWE model ', &
363  trim(gwemodel%name)
364  gwfconnidx = ic2
365  call this%link_connections(gweconn, gwfconn)
366  exit gwfloop
367  end if
368  end if
369  end do gwfloop
370  !
371  ! fallback option: coupling to old gwfgwf exchange,
372  ! (this will go obsolete at some point)
373  if (gwfconnidx == -1) then
374  gwfloopexg: do iex = 1, baseexchangelist%Count()
376  !
377  ! -- There is no guarantee that iex is a gwfExg, in which case
378  ! it will return as null. cycle if so.
379  if (.not. associated(gwfex)) cycle gwfloopexg
380  !
381  if (associated(gwfex%model1, gwfmodel) .or. &
382  associated(gwfex%model2, gwfmodel)) then
383 
384  ! check exchanges have same node counts
385  areequal = size(gwfex%nodem1) == size(gweconn%prim_exchange%nodem1)
386  ! then, connecting the same nodes nrs will be
387  ! sufficient evidence of equality
388  if (areequal) &
389  areequal = all(gwfex%nodem1 == gweconn%prim_exchange%nodem1)
390  if (areequal) &
391  areequal = all(gwfex%nodem2 == gweconn%prim_exchange%nodem2)
392  if (areequal) then
393  ! link exchange to connection
394  write (iout, '(/6a)') 'Linking exchange ', &
395  trim(gweconn%prim_exchange%name), &
396  ' to ', trim(gwfex%name), ' for GWE model ', &
397  trim(gwemodel%name)
398  gwfexidx = iex
399  if (gweconn%owns_exchange) then
400  gweconn%gweExchange%gwfsimvals => gwfex%simvals
401  call mem_checkin(gweconn%gweExchange%gwfsimvals, &
402  'GWFSIMVALS', gweconn%gweExchange%memoryPath, &
403  'SIMVALS', gwfex%memoryPath)
404  end if
405  !
406  !cdl link up mvt to mvr
407  if (gwfex%inmvr > 0) then
408  if (gweconn%owns_exchange) then
409  !cdl todo: check and make sure gweEx has mvt active
410  call gweconn%gweExchange%mvt%set_pointer_mvrbudobj( &
411  gwfex%mvr%budobj)
412  end if
413  end if
414  !
415  if (associated(gwfex%model2, gwfmodel)) gweconn%exgflowSign = -1
416  gweconn%gweInterfaceModel%fmi%flows_from_file = .false.
417  !
418  exit gwfloopexg
419  end if
420  end if
421  !
422  end do gwfloopexg
423  end if
424  !
425  if (gwfconnidx == -1 .and. gwfexidx == -1) then
426  ! none found, report
427  write (errmsg, '(/6a)') 'Missing GWF-GWF exchange when connecting GWE'// &
428  ' model ', trim(gwemodel%name), ' with exchange ', &
429  trim(gweconn%prim_exchange%name), ' to GWF model ', &
430  trim(gwfmodel%name)
431  call store_error(errmsg, terminate=.true.)
432  end if
433  !
434  end do gweloop
435  end subroutine gwfconn2gweconn
436 
437  !> @brief Links a GWE connection to its GWF counterpart
438  !<
439  subroutine link_connections(this, gweConn, gwfConn)
440  ! -- modules
442  ! -- dummy
443  class(gwfgweexchangetype) :: this !< this exchange
444  class(gwegweconnectiontype), pointer :: gweConn !< GWE connection
445  class(gwfgwfconnectiontype), pointer :: gwfConn !< GWF connection
446  !
447  !gweConn%exgflowja => gwfConn%exgflowja
448  if (gweconn%owns_exchange) then
449  gweconn%gweExchange%gwfsimvals => gwfconn%gwfExchange%simvals
450  call mem_checkin(gweconn%gweExchange%gwfsimvals, &
451  'GWFSIMVALS', gweconn%gweExchange%memoryPath, &
452  'SIMVALS', gwfconn%gwfExchange%memoryPath)
453  end if
454  !
455  !cdl link up mvt to mvr
456  if (gwfconn%gwfExchange%inmvr > 0) then
457  if (gweconn%owns_exchange) then
458  !cdl todo: check and make sure gweEx has mvt active
459  call gweconn%gweExchange%mvt%set_pointer_mvrbudobj( &
460  gwfconn%gwfExchange%mvr%budobj)
461  end if
462  end if
463  !
464  if (associated(gwfconn%gwfExchange%model2, gwfconn%owner)) then
465  gweconn%exgflowSign = -1
466  end if
467  !
468  ! fmi flows are not read from file
469  gweconn%gweInterfaceModel%fmi%flows_from_file = .false.
470  !
471  ! set concentration pointer for buoyancy
472  !call gwfConn%gwfInterfaceModel%buy%set_concentration_pointer( &
473  ! gweConn%gweModel%name, &
474  ! gweConn%conc, &
475  ! gweConn%icbound)
476  end subroutine link_connections
477 
478  !> @brief Deallocate memory
479  !<
480  subroutine exg_da(this)
481  ! -- modules
483  ! -- dummy
484  class(gwfgweexchangetype) :: this
485  !
486  call mem_deallocate(this%m1_idx)
487  call mem_deallocate(this%m2_idx)
488  end subroutine exg_da
489 
490  !> @brief Allocate GwfGwe exchange scalars
491  !<
492  subroutine allocate_scalars(this)
493  ! -- modules
495  ! -- dummy
496  class(gwfgweexchangetype) :: this
497  !
498  call mem_allocate(this%m1_idx, 'M1ID', this%memoryPath)
499  call mem_allocate(this%m2_idx, 'M2ID', this%memoryPath)
500  this%m1_idx = 0
501  this%m2_idx = 0
502  end subroutine allocate_scalars
503 
504  !> @brief Call routines in FMI that will set pointers to the necessary flow
505  !! data (SIMVALS and SIMTOMVR) stored within each GWF flow package
506  !<
507  subroutine gwfbnd2gwefmi(this)
508  ! -- dummy
509  class(gwfgweexchangetype) :: this
510  ! -- local
511  integer(I4B) :: ngwfpack, ip, iterm, imover
512  class(basemodeltype), pointer :: mb => null()
513  type(gwfmodeltype), pointer :: gwfmodel => null()
514  type(gwemodeltype), pointer :: gwemodel => null()
515  class(bndtype), pointer :: packobj => null()
516  !
517  ! -- set gwfmodel
518  mb => getbasemodelfromlist(basemodellist, this%m1_idx)
519  select type (mb)
520  type is (gwfmodeltype)
521  gwfmodel => mb
522  end select
523  !
524  ! -- set gwemodel
525  mb => getbasemodelfromlist(basemodellist, this%m2_idx)
526  select type (mb)
527  type is (gwemodeltype)
528  gwemodel => mb
529  end select
530  !
531  ! -- Call routines in FMI that will set pointers to the necessary flow
532  ! data (SIMVALS and SIMTOMVR) stored within each GWF flow package
533  ngwfpack = gwfmodel%bndlist%Count()
534  iterm = 1
535  do ip = 1, ngwfpack
536  packobj => getbndfromlist(gwfmodel%bndlist, ip)
537  call gwemodel%fmi%gwfpackages(iterm)%set_pointers( &
538  'SIMVALS', &
539  packobj%memoryPath, packobj%input_mempath)
540  iterm = iterm + 1
541  !
542  ! -- If a mover is active for this package, then establish a separate
543  ! pointer link for the mover flows stored in SIMTOMVR
544  imover = packobj%imover
545  if (packobj%isadvpak /= 0) imover = 0
546  if (imover /= 0) then
547  call gwemodel%fmi%gwfpackages(iterm)%set_pointers( &
548  'SIMTOMVR', &
549  packobj%memoryPath, packobj%input_mempath)
550  iterm = iterm + 1
551  end if
552  end do
553  end subroutine gwfbnd2gwefmi
554 
555 end module gwfgweexchangemodule
subroutine, public addbaseexchangetolist(list, exchange)
Add the exchange object (BaseExchangeType) to a list.
class(basemodeltype) function, pointer, public getbasemodelfromlist(list, idx)
Definition: BaseModel.f90:172
This module contains the base boundary package.
class(bndtype) function, pointer, public getbndfromlist(list, idx)
Get boundary from package list.
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter lenpackagename
maximum length of the package name
Definition: Constants.f90:23
Definition: Dis.f90:1
class(gwegweconnectiontype) function, pointer, public castasgwegweconnection(obj)
Cast to GweGweConnectionType.
Definition: gwe.f90:3
subroutine gwfbnd2gwefmi(this)
Call routines in FMI that will set pointers to the necessary flow data (SIMVALS and SIMTOMVR) stored ...
Definition: exg-gwfgwe.f90:508
subroutine gwfconn2gweconn(this, gwfModel, gweModel)
Link GWE connections to GWF connections or exchanges.
Definition: exg-gwfgwe.f90:313
subroutine set_model_pointers(this)
Allocate and read.
Definition: exg-gwfgwe.f90:84
subroutine allocate_scalars(this)
Allocate GwfGwe exchange scalars.
Definition: exg-gwfgwe.f90:493
subroutine exg_da(this)
Deallocate memory.
Definition: exg-gwfgwe.f90:481
subroutine, public gwfgwe_cr(filename, id, m1_id, m2_id)
Create a new GWF to GWE exchange object.
Definition: exg-gwfgwe.f90:47
subroutine link_connections(this, gweConn, gwfConn)
Links a GWE connection to its GWF counterpart.
Definition: exg-gwfgwe.f90:440
class(gwfgwfconnectiontype) function, pointer, public castasgwfgwfconnection(obj)
Cast to GwfGwfConnectionType.
This module contains the GwfGwfExchangeModule Module.
Definition: exg-gwfgwf.f90:10
class(gwfexchangetype) function, pointer, public getgwfexchangefromlist(list, idx)
@ brief Get exchange from list
Definition: gwf.f90:1
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 baseconnectionlist
Definition: mf6lists.f90:28
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=maxcharlen) errmsg
error message string
integer(i4b), dimension(:), allocatable model_loc_idx
equals the local index into the basemodel list (-1 when not available)
integer(i4b) iout
file unit number for simulation output
class(spatialmodelconnectiontype) function, pointer, public get_smc_from_list(list, idx)
Get the connection from a list.
Highest level model type. All models extend this parent type.
Definition: BaseModel.f90:13
@ brief BndType
Structured grid discretization.
Definition: Dis.f90:23
Unstructured grid discretization.
Definition: Disu.f90:28
Vertex grid discretization.
Definition: Disv.f90:24
Connects a GWE model to other GWE models in space. Derives from NumericalExchangeType so the solution...
Connecting a GWF model to other models in space, implements NumericalExchangeType so the solution can...
Derived type for GwfExchangeType.
Definition: exg-gwfgwf.f90:47
Class to manage spatial connection of a model to one or more models of the same type....