MODFLOW 6  version 6.6.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  ! -- dummy
185  class(gwfgweexchangetype) :: this
186  ! -- local
187  class(basemodeltype), pointer :: mb => null()
188  type(gwfmodeltype), pointer :: gwfmodel => null()
189  type(gwemodeltype), pointer :: gwemodel => null()
190  ! -- formats
191  character(len=*), parameter :: fmtdiserr = &
192  "('GWF and GWE Models do not have the same discretization for exchange&
193  & ',a,'.&
194  & GWF Model has ', i0, ' user nodes and ', i0, ' reduced nodes.&
195  & GWE Model has ', i0, ' user nodes and ', i0, ' reduced nodes.&
196  & Ensure discretization packages, including IDOMAIN, are identical.')"
197  !
198  ! -- set gwfmodel
199  mb => getbasemodelfromlist(basemodellist, this%m1_idx)
200  select type (mb)
201  type is (gwfmodeltype)
202  gwfmodel => mb
203  end select
204  !
205  ! -- set gwemodel
206  mb => getbasemodelfromlist(basemodellist, this%m2_idx)
207  select type (mb)
208  type is (gwemodeltype)
209  gwemodel => mb
210  end select
211  !
212  ! -- Check to make sure sizes are identical
213  if (gwemodel%dis%nodes /= gwfmodel%dis%nodes .or. &
214  gwemodel%dis%nodesuser /= gwfmodel%dis%nodesuser) then
215  write (errmsg, fmtdiserr) trim(this%name), &
216  gwfmodel%dis%nodesuser, &
217  gwfmodel%dis%nodes, &
218  gwemodel%dis%nodesuser, &
219  gwemodel%dis%nodes
220  call store_error(errmsg, terminate=.true.)
221  end if
222  !
223  ! -- setup pointers to gwf variables allocated in gwf_ar
224  gwemodel%fmi%gwfhead => gwfmodel%x
225  call mem_checkin(gwemodel%fmi%gwfhead, &
226  'GWFHEAD', gwemodel%fmi%memoryPath, &
227  'X', gwfmodel%memoryPath)
228  gwemodel%fmi%gwfsat => gwfmodel%npf%sat
229  call mem_checkin(gwemodel%fmi%gwfsat, &
230  'GWFSAT', gwemodel%fmi%memoryPath, &
231  'SAT', gwfmodel%npf%memoryPath)
232  gwemodel%fmi%gwfspdis => gwfmodel%npf%spdis
233  call mem_checkin(gwemodel%fmi%gwfspdis, &
234  'GWFSPDIS', gwemodel%fmi%memoryPath, &
235  'SPDIS', gwfmodel%npf%memoryPath)
236  !
237  ! -- setup pointers to the flow storage rates. GWF strg arrays are
238  ! available after the gwf_ar routine is called.
239  if (gwemodel%inest > 0) then
240  if (gwfmodel%insto > 0) then
241  gwemodel%fmi%gwfstrgss => gwfmodel%sto%strgss
242  gwemodel%fmi%igwfstrgss = 1
243  if (gwfmodel%sto%iusesy == 1) then
244  gwemodel%fmi%gwfstrgsy => gwfmodel%sto%strgsy
245  gwemodel%fmi%igwfstrgsy = 1
246  end if
247  end if
248  end if
249  !
250  ! -- Set a pointer to conc in buy
251  if (gwfmodel%inbuy > 0) then
252  call gwfmodel%buy%set_concentration_pointer(gwemodel%name, gwemodel%x, &
253  gwemodel%ibound)
254  end if
255  !
256  ! -- Set a pointer to conc (which could be a temperature) in vsc
257  if (gwfmodel%invsc > 0) then
258  call gwfmodel%vsc%set_concentration_pointer(gwemodel%name, gwemodel%x, &
259  gwemodel%ibound, 1)
260  end if
261  !
262  ! -- transfer the boundary package information from gwf to gwe
263  call this%gwfbnd2gwefmi()
264  !
265  ! -- if mover package is active, then set a pointer to it's budget object
266  if (gwfmodel%inmvr /= 0) then
267  gwemodel%fmi%mvrbudobj => gwfmodel%mvr%budobj
268  end if
269  !
270  ! -- connect Connections
271  call this%gwfconn2gweconn(gwfmodel, gwemodel)
272  end subroutine exg_ar
273 
274  !> @brief Link GWE connections to GWF connections or exchanges
275  !<
276  subroutine gwfconn2gweconn(this, gwfModel, gweModel)
277  ! -- modules
278  use simmodule, only: store_error
279  use simvariablesmodule, only: iout
281  ! -- dummy
282  class(gwfgweexchangetype) :: this !< this exchange
283  type(gwfmodeltype), pointer :: gwfModel !< the flow model
284  type(gwemodeltype), pointer :: gweModel !< the energy transport model
285  ! -- local
286  class(spatialmodelconnectiontype), pointer :: conn => null()
287  class(*), pointer :: objPtr => null()
288  class(gwegweconnectiontype), pointer :: gweConn => null()
289  class(gwfgwfconnectiontype), pointer :: gwfConn => null()
290  class(gwfexchangetype), pointer :: gwfEx => null()
291  integer(I4B) :: ic1, ic2, iex
292  integer(I4B) :: gwfConnIdx, gwfExIdx
293  logical(LGP) :: areEqual
294  !
295  ! loop over all connections
296  gweloop: do ic1 = 1, baseconnectionlist%Count()
297  !
299  if (.not. associated(conn%owner, gwemodel)) cycle gweloop
300  !
301  ! start with a GWE conn.
302  objptr => conn
303  gweconn => castasgwegweconnection(objptr)
304  gwfconnidx = -1
305  gwfexidx = -1
306  !
307  ! find matching GWF conn. in same list
308  gwfloop: do ic2 = 1, baseconnectionlist%Count()
310  !
311  if (associated(conn%owner, gwfmodel)) then
312  objptr => conn
313  gwfconn => castasgwfgwfconnection(objptr)
314  !
315  ! for now, connecting the same nodes nrs will be
316  ! sufficient evidence of equality
317  areequal = all(gwfconn%prim_exchange%nodem1 == &
318  gweconn%prim_exchange%nodem1)
319  areequal = areequal .and. all(gwfconn%prim_exchange%nodem2 == &
320  gweconn%prim_exchange%nodem2)
321  if (areequal) then
322  ! same DIS, same exchange: link and go to next GWE conn.
323  write (iout, '(/6a)') 'Linking exchange ', &
324  trim(gweconn%prim_exchange%name), &
325  ' to ', trim(gwfconn%prim_exchange%name), &
326  ' (using interface model) for GWE model ', &
327  trim(gwemodel%name)
328  gwfconnidx = ic2
329  call this%link_connections(gweconn, gwfconn)
330  exit gwfloop
331  end if
332  end if
333  end do gwfloop
334  !
335  ! fallback option: coupling to old gwfgwf exchange,
336  ! (this will go obsolete at some point)
337  if (gwfconnidx == -1) then
338  gwfloopexg: do iex = 1, baseexchangelist%Count()
340  !
341  ! -- There is no guarantee that iex is a gwfExg, in which case
342  ! it will return as null. cycle if so.
343  if (.not. associated(gwfex)) cycle gwfloopexg
344  !
345  if (associated(gwfex%model1, gwfmodel) .or. &
346  associated(gwfex%model2, gwfmodel)) then
347 
348  ! check exchanges have same node counts
349  areequal = size(gwfex%nodem1) == size(gweconn%prim_exchange%nodem1)
350  ! then, connecting the same nodes nrs will be
351  ! sufficient evidence of equality
352  if (areequal) &
353  areequal = all(gwfex%nodem1 == gweconn%prim_exchange%nodem1)
354  if (areequal) &
355  areequal = all(gwfex%nodem2 == gweconn%prim_exchange%nodem2)
356  if (areequal) then
357  ! link exchange to connection
358  write (iout, '(/6a)') 'Linking exchange ', &
359  trim(gweconn%prim_exchange%name), &
360  ' to ', trim(gwfex%name), ' for GWE model ', &
361  trim(gwemodel%name)
362  gwfexidx = iex
363  if (gweconn%owns_exchange) then
364  gweconn%gweExchange%gwfsimvals => gwfex%simvals
365  call mem_checkin(gweconn%gweExchange%gwfsimvals, &
366  'GWFSIMVALS', gweconn%gweExchange%memoryPath, &
367  'SIMVALS', gwfex%memoryPath)
368  end if
369  !
370  !cdl link up mvt to mvr
371  if (gwfex%inmvr > 0) then
372  if (gweconn%owns_exchange) then
373  !cdl todo: check and make sure gweEx has mvt active
374  call gweconn%gweExchange%mvt%set_pointer_mvrbudobj( &
375  gwfex%mvr%budobj)
376  end if
377  end if
378  !
379  if (associated(gwfex%model2, gwfmodel)) gweconn%exgflowSign = -1
380  gweconn%gweInterfaceModel%fmi%flows_from_file = .false.
381  !
382  exit gwfloopexg
383  end if
384  end if
385  !
386  end do gwfloopexg
387  end if
388  !
389  if (gwfconnidx == -1 .and. gwfexidx == -1) then
390  ! none found, report
391  write (errmsg, '(/6a)') 'Missing GWF-GWF exchange when connecting GWE'// &
392  ' model ', trim(gwemodel%name), ' with exchange ', &
393  trim(gweconn%prim_exchange%name), ' to GWF model ', &
394  trim(gwfmodel%name)
395  call store_error(errmsg, terminate=.true.)
396  end if
397  !
398  end do gweloop
399  end subroutine gwfconn2gweconn
400 
401  !> @brief Links a GWE connection to its GWF counterpart
402  !<
403  subroutine link_connections(this, gweConn, gwfConn)
404  ! -- modules
406  ! -- dummy
407  class(gwfgweexchangetype) :: this !< this exchange
408  class(gwegweconnectiontype), pointer :: gweConn !< GWE connection
409  class(gwfgwfconnectiontype), pointer :: gwfConn !< GWF connection
410  !
411  !gweConn%exgflowja => gwfConn%exgflowja
412  if (gweconn%owns_exchange) then
413  gweconn%gweExchange%gwfsimvals => gwfconn%gwfExchange%simvals
414  call mem_checkin(gweconn%gweExchange%gwfsimvals, &
415  'GWFSIMVALS', gweconn%gweExchange%memoryPath, &
416  'SIMVALS', gwfconn%gwfExchange%memoryPath)
417  end if
418  !
419  !cdl link up mvt to mvr
420  if (gwfconn%gwfExchange%inmvr > 0) then
421  if (gweconn%owns_exchange) then
422  !cdl todo: check and make sure gweEx has mvt active
423  call gweconn%gweExchange%mvt%set_pointer_mvrbudobj( &
424  gwfconn%gwfExchange%mvr%budobj)
425  end if
426  end if
427  !
428  if (associated(gwfconn%gwfExchange%model2, gwfconn%owner)) then
429  gweconn%exgflowSign = -1
430  end if
431  !
432  ! fmi flows are not read from file
433  gweconn%gweInterfaceModel%fmi%flows_from_file = .false.
434  !
435  ! set concentration pointer for buoyancy
436  !call gwfConn%gwfInterfaceModel%buy%set_concentration_pointer( &
437  ! gweConn%gweModel%name, &
438  ! gweConn%conc, &
439  ! gweConn%icbound)
440  end subroutine link_connections
441 
442  !> @brief Deallocate memory
443  !<
444  subroutine exg_da(this)
445  ! -- modules
447  ! -- dummy
448  class(gwfgweexchangetype) :: this
449  !
450  call mem_deallocate(this%m1_idx)
451  call mem_deallocate(this%m2_idx)
452  end subroutine exg_da
453 
454  !> @brief Allocate GwfGwe exchange scalars
455  !<
456  subroutine allocate_scalars(this)
457  ! -- modules
459  ! -- dummy
460  class(gwfgweexchangetype) :: this
461  !
462  call mem_allocate(this%m1_idx, 'M1ID', this%memoryPath)
463  call mem_allocate(this%m2_idx, 'M2ID', this%memoryPath)
464  this%m1_idx = 0
465  this%m2_idx = 0
466  end subroutine allocate_scalars
467 
468  !> @brief Call routines in FMI that will set pointers to the necessary flow
469  !! data (SIMVALS and SIMTOMVR) stored within each GWF flow package
470  !<
471  subroutine gwfbnd2gwefmi(this)
472  ! -- dummy
473  class(gwfgweexchangetype) :: this
474  ! -- local
475  integer(I4B) :: ngwfpack, ip, iterm, imover
476  class(basemodeltype), pointer :: mb => null()
477  type(gwfmodeltype), pointer :: gwfmodel => null()
478  type(gwemodeltype), pointer :: gwemodel => null()
479  class(bndtype), pointer :: packobj => null()
480  !
481  ! -- set gwfmodel
482  mb => getbasemodelfromlist(basemodellist, this%m1_idx)
483  select type (mb)
484  type is (gwfmodeltype)
485  gwfmodel => mb
486  end select
487  !
488  ! -- set gwemodel
489  mb => getbasemodelfromlist(basemodellist, this%m2_idx)
490  select type (mb)
491  type is (gwemodeltype)
492  gwemodel => mb
493  end select
494  !
495  ! -- Call routines in FMI that will set pointers to the necessary flow
496  ! data (SIMVALS and SIMTOMVR) stored within each GWF flow package
497  ngwfpack = gwfmodel%bndlist%Count()
498  iterm = 1
499  do ip = 1, ngwfpack
500  packobj => getbndfromlist(gwfmodel%bndlist, ip)
501  call gwemodel%fmi%gwfpackages(iterm)%set_pointers( &
502  'SIMVALS', &
503  packobj%memoryPath, packobj%input_mempath)
504  iterm = iterm + 1
505  !
506  ! -- If a mover is active for this package, then establish a separate
507  ! pointer link for the mover flows stored in SIMTOMVR
508  imover = packobj%imover
509  if (packobj%isadvpak /= 0) imover = 0
510  if (imover /= 0) then
511  call gwemodel%fmi%gwfpackages(iterm)%set_pointers( &
512  'SIMTOMVR', &
513  packobj%memoryPath, packobj%input_mempath)
514  iterm = iterm + 1
515  end if
516  end do
517  end subroutine gwfbnd2gwefmi
518 
519 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
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:472
subroutine gwfconn2gweconn(this, gwfModel, gweModel)
Link GWE connections to GWF connections or exchanges.
Definition: exg-gwfgwe.f90:277
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:457
subroutine exg_da(this)
Deallocate memory.
Definition: exg-gwfgwe.f90:445
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:404
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
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....