MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
gwfgweexchangemodule Module Reference

Data Types

type  gwfgweexchangetype
 

Functions/Subroutines

subroutine, public gwfgwe_cr (filename, id, m1_id, m2_id)
 Create a new GWF to GWE exchange object. More...
 
subroutine set_model_pointers (this)
 Allocate and read. More...
 
subroutine exg_df (this)
 Define the GwfGwe Exchange object. More...
 
subroutine exg_ar (this)
 Allocate and read. More...
 
subroutine gwfconn2gweconn (this, gwfModel, gweModel)
 Link GWE connections to GWF connections or exchanges. More...
 
subroutine link_connections (this, gweConn, gwfConn)
 Links a GWE connection to its GWF counterpart. More...
 
subroutine exg_da (this)
 Deallocate memory. More...
 
subroutine allocate_scalars (this)
 Allocate GwfGwe exchange scalars. More...
 
subroutine gwfbnd2gwefmi (this)
 Call routines in FMI that will set pointers to the necessary flow data (SIMVALS and SIMTOMVR) stored within each GWF flow package. More...
 

Function/Subroutine Documentation

◆ allocate_scalars()

subroutine gwfgweexchangemodule::allocate_scalars ( class(gwfgweexchangetype this)

Definition at line 477 of file exg-gwfgwe.f90.

478  ! -- modules
480  ! -- dummy
481  class(GwfGweExchangeType) :: this
482  !
483  call mem_allocate(this%m1_idx, 'M1ID', this%memoryPath)
484  call mem_allocate(this%m2_idx, 'M2ID', this%memoryPath)
485  this%m1_idx = 0
486  this%m2_idx = 0
487  !
488  ! -- Return
489  return

◆ exg_ar()

subroutine gwfgweexchangemodule::exg_ar ( class(gwfgweexchangetype this)

Definition at line 190 of file exg-gwfgwe.f90.

191  ! -- modules
193  ! -- dummy
194  class(GwfGweExchangeType) :: this
195  ! -- local
196  class(BaseModelType), pointer :: mb => null()
197  type(GwfModelType), pointer :: gwfmodel => null()
198  type(GweModelType), pointer :: gwemodel => null()
199  ! -- formats
200  character(len=*), parameter :: fmtdiserr = &
201  "('GWF and GWE Models do not have the same discretization for exchange&
202  & ',a,'.&
203  & GWF Model has ', i0, ' user nodes and ', i0, ' reduced nodes.&
204  & GWE Model has ', i0, ' user nodes and ', i0, ' reduced nodes.&
205  & Ensure discretization packages, including IDOMAIN, are identical.')"
206  !
207  ! -- set gwfmodel
208  mb => getbasemodelfromlist(basemodellist, this%m1_idx)
209  select type (mb)
210  type is (gwfmodeltype)
211  gwfmodel => mb
212  end select
213  !
214  ! -- set gwemodel
215  mb => getbasemodelfromlist(basemodellist, this%m2_idx)
216  select type (mb)
217  type is (gwemodeltype)
218  gwemodel => mb
219  end select
220  !
221  ! -- Check to make sure sizes are identical
222  if (gwemodel%dis%nodes /= gwfmodel%dis%nodes .or. &
223  gwemodel%dis%nodesuser /= gwfmodel%dis%nodesuser) then
224  write (errmsg, fmtdiserr) trim(this%name), &
225  gwfmodel%dis%nodesuser, &
226  gwfmodel%dis%nodes, &
227  gwemodel%dis%nodesuser, &
228  gwemodel%dis%nodes
229  call store_error(errmsg, terminate=.true.)
230  end if
231  !
232  ! -- setup pointers to gwf variables allocated in gwf_ar
233  gwemodel%fmi%gwfhead => gwfmodel%x
234  call mem_checkin(gwemodel%fmi%gwfhead, &
235  'GWFHEAD', gwemodel%fmi%memoryPath, &
236  'X', gwfmodel%memoryPath)
237  gwemodel%fmi%gwfsat => gwfmodel%npf%sat
238  call mem_checkin(gwemodel%fmi%gwfsat, &
239  'GWFSAT', gwemodel%fmi%memoryPath, &
240  'SAT', gwfmodel%npf%memoryPath)
241  gwemodel%fmi%gwfspdis => gwfmodel%npf%spdis
242  call mem_checkin(gwemodel%fmi%gwfspdis, &
243  'GWFSPDIS', gwemodel%fmi%memoryPath, &
244  'SPDIS', gwfmodel%npf%memoryPath)
245  !
246  ! -- setup pointers to the flow storage rates. GWF strg arrays are
247  ! available after the gwf_ar routine is called.
248  if (gwemodel%inest > 0) then
249  if (gwfmodel%insto > 0) then
250  gwemodel%fmi%gwfstrgss => gwfmodel%sto%strgss
251  gwemodel%fmi%igwfstrgss = 1
252  if (gwfmodel%sto%iusesy == 1) then
253  gwemodel%fmi%gwfstrgsy => gwfmodel%sto%strgsy
254  gwemodel%fmi%igwfstrgsy = 1
255  end if
256  end if
257  end if
258  !
259  ! -- Set a pointer to conc in buy
260  if (gwfmodel%inbuy > 0) then
261  call gwfmodel%buy%set_concentration_pointer(gwemodel%name, gwemodel%x, &
262  gwemodel%ibound)
263  end if
264  !
265  ! -- Set a pointer to conc (which could be a temperature) in vsc
266  if (gwfmodel%invsc > 0) then
267  call gwfmodel%vsc%set_concentration_pointer(gwemodel%name, gwemodel%x, &
268  gwemodel%ibound, 1)
269  end if
270  !
271  ! -- transfer the boundary package information from gwf to gwe
272  call this%gwfbnd2gwefmi()
273  !
274  ! -- if mover package is active, then set a pointer to it's budget object
275  if (gwfmodel%inmvr /= 0) then
276  gwemodel%fmi%mvrbudobj => gwfmodel%mvr%budobj
277  end if
278  !
279  ! -- connect Connections
280  call this%gwfconn2gweconn(gwfmodel, gwemodel)
281  !
282  ! -- Return
283  return
Here is the call graph for this function:

◆ exg_da()

subroutine gwfgweexchangemodule::exg_da ( class(gwfgweexchangetype this)

Definition at line 462 of file exg-gwfgwe.f90.

463  ! -- modules
465  ! -- dummy
466  class(GwfGweExchangeType) :: this
467  !
468  call mem_deallocate(this%m1_idx)
469  call mem_deallocate(this%m2_idx)
470  !
471  ! -- Return
472  return

◆ exg_df()

subroutine gwfgweexchangemodule::exg_df ( class(gwfgweexchangetype this)

Definition at line 137 of file exg-gwfgwe.f90.

138  ! -- modules
140  ! -- dummy
141  class(GwfGweExchangeType) :: this
142  ! -- local
143  class(BaseModelType), pointer :: mb => null()
144  type(GwfModelType), pointer :: gwfmodel => null()
145  type(GweModelType), pointer :: gwemodel => null()
146  !
147  ! -- set gwfmodel
148  mb => getbasemodelfromlist(basemodellist, this%m1_idx)
149  select type (mb)
150  type is (gwfmodeltype)
151  gwfmodel => mb
152  end select
153  !
154  ! -- set gwemodel
155  mb => getbasemodelfromlist(basemodellist, this%m2_idx)
156  select type (mb)
157  type is (gwemodeltype)
158  gwemodel => mb
159  end select
160  !
161  ! -- Check to make sure that flow is solved before transport and in a
162  ! different IMS solution
163  if (gwfmodel%idsoln >= gwemodel%idsoln) then
164  write (errmsg, '(3a)') 'Problem with GWF-GWE exchange ', trim(this%name), &
165  '. The GWF model must be solved by a different IMS than the GWE model. &
166  &Furthermore, the IMS specified for GWF must be listed in mfsim.nam &
167  &before the IMS for GWE.'
168  call store_error(errmsg, terminate=.true.)
169  end if
170  !
171  ! -- Set pointer to flowja
172  gwemodel%fmi%gwfflowja => gwfmodel%flowja
173  call mem_checkin(gwemodel%fmi%gwfflowja, &
174  'GWFFLOWJA', gwemodel%fmi%memoryPath, &
175  'FLOWJA', gwfmodel%memoryPath)
176 
177  !
178  ! -- Set the npf flag so that specific discharge is available for
179  ! transport calculations if dispersion is active
180  if (gwemodel%incnd > 0) then
181  gwfmodel%npf%icalcspdis = 1
182  end if
183  !
184  ! -- Return
185  return
Here is the call graph for this function:

◆ gwfbnd2gwefmi()

subroutine gwfgweexchangemodule::gwfbnd2gwefmi ( class(gwfgweexchangetype this)

Definition at line 495 of file exg-gwfgwe.f90.

496  ! -- dummy
497  class(GwfGweExchangeType) :: this
498  ! -- local
499  integer(I4B) :: ngwfpack, ip, iterm, imover
500  class(BaseModelType), pointer :: mb => null()
501  type(GwfModelType), pointer :: gwfmodel => null()
502  type(GweModelType), pointer :: gwemodel => null()
503  class(BndType), pointer :: packobj => null()
504  !
505  ! -- set gwfmodel
506  mb => getbasemodelfromlist(basemodellist, this%m1_idx)
507  select type (mb)
508  type is (gwfmodeltype)
509  gwfmodel => mb
510  end select
511  !
512  ! -- set gwemodel
513  mb => getbasemodelfromlist(basemodellist, this%m2_idx)
514  select type (mb)
515  type is (gwemodeltype)
516  gwemodel => mb
517  end select
518  !
519  ! -- Call routines in FMI that will set pointers to the necessary flow
520  ! data (SIMVALS and SIMTOMVR) stored within each GWF flow package
521  ngwfpack = gwfmodel%bndlist%Count()
522  iterm = 1
523  do ip = 1, ngwfpack
524  packobj => getbndfromlist(gwfmodel%bndlist, ip)
525  call gwemodel%fmi%gwfpackages(iterm)%set_pointers( &
526  'SIMVALS', &
527  packobj%memoryPath, packobj%input_mempath)
528  iterm = iterm + 1
529  !
530  ! -- If a mover is active for this package, then establish a separate
531  ! pointer link for the mover flows stored in SIMTOMVR
532  imover = packobj%imover
533  if (packobj%isadvpak /= 0) imover = 0
534  if (imover /= 0) then
535  call gwemodel%fmi%gwfpackages(iterm)%set_pointers( &
536  'SIMTOMVR', &
537  packobj%memoryPath, packobj%input_mempath)
538  iterm = iterm + 1
539  end if
540  end do
541  !
542  ! -- Return
543  return
Here is the call graph for this function:

◆ gwfconn2gweconn()

subroutine gwfgweexchangemodule::gwfconn2gweconn ( class(gwfgweexchangetype this,
type(gwfmodeltype), pointer  gwfModel,
type(gwemodeltype), pointer  gweModel 
)
Parameters
thisthis exchange
gwfmodelthe flow model
gwemodelthe energy transport model

Definition at line 288 of file exg-gwfgwe.f90.

289  ! -- modules
290  use simmodule, only: store_error
291  use simvariablesmodule, only: iout
293  ! -- dummy
294  class(GwfGweExchangeType) :: this !< this exchange
295  type(GwfModelType), pointer :: gwfModel !< the flow model
296  type(GweModelType), pointer :: gweModel !< the energy transport model
297  ! -- local
298  class(SpatialModelConnectionType), pointer :: conn => null()
299  class(*), pointer :: objPtr => null()
300  class(GweGweConnectionType), pointer :: gweConn => null()
301  class(GwfGwfConnectionType), pointer :: gwfConn => null()
302  class(GwfExchangeType), pointer :: gwfEx => null()
303  integer(I4B) :: ic1, ic2, iex
304  integer(I4B) :: gwfConnIdx, gwfExIdx
305  logical(LGP) :: areEqual
306  !
307  ! loop over all connections
308  gweloop: do ic1 = 1, baseconnectionlist%Count()
309  !
310  conn => get_smc_from_list(baseconnectionlist, ic1)
311  if (.not. associated(conn%owner, gwemodel)) cycle gweloop
312  !
313  ! start with a GWE conn.
314  objptr => conn
315  gweconn => castasgwegweconnection(objptr)
316  gwfconnidx = -1
317  gwfexidx = -1
318  !
319  ! find matching GWF conn. in same list
320  gwfloop: do ic2 = 1, baseconnectionlist%Count()
321  conn => get_smc_from_list(baseconnectionlist, ic2)
322  !
323  if (associated(conn%owner, gwfmodel)) then
324  objptr => conn
325  gwfconn => castasgwfgwfconnection(objptr)
326  !
327  ! for now, connecting the same nodes nrs will be
328  ! sufficient evidence of equality
329  areequal = all(gwfconn%prim_exchange%nodem1 == &
330  gweconn%prim_exchange%nodem1)
331  areequal = areequal .and. all(gwfconn%prim_exchange%nodem2 == &
332  gweconn%prim_exchange%nodem2)
333  if (areequal) then
334  ! same DIS, same exchange: link and go to next GWE conn.
335  write (iout, '(/6a)') 'Linking exchange ', &
336  trim(gweconn%prim_exchange%name), &
337  ' to ', trim(gwfconn%prim_exchange%name), &
338  ' (using interface model) for GWE model ', &
339  trim(gwemodel%name)
340  gwfconnidx = ic2
341  call this%link_connections(gweconn, gwfconn)
342  exit gwfloop
343  end if
344  end if
345  end do gwfloop
346  !
347  ! fallback option: coupling to old gwfgwf exchange,
348  ! (this will go obsolete at some point)
349  if (gwfconnidx == -1) then
350  gwfloopexg: do iex = 1, baseexchangelist%Count()
351  gwfex => getgwfexchangefromlist(baseexchangelist, iex)
352  !
353  ! -- There is no guarantee that iex is a gwfExg, in which case
354  ! it will return as null. cycle if so.
355  if (.not. associated(gwfex)) cycle gwfloopexg
356  !
357  if (associated(gwfex%model1, gwfmodel) .or. &
358  associated(gwfex%model2, gwfmodel)) then
359 
360  ! check exchanges have same node counts
361  areequal = size(gwfex%nodem1) == size(gweconn%prim_exchange%nodem1)
362  ! then, connecting the same nodes nrs will be
363  ! sufficient evidence of equality
364  if (areequal) &
365  areequal = all(gwfex%nodem1 == gweconn%prim_exchange%nodem1)
366  if (areequal) &
367  areequal = all(gwfex%nodem2 == gweconn%prim_exchange%nodem2)
368  if (areequal) then
369  ! link exchange to connection
370  write (iout, '(/6a)') 'Linking exchange ', &
371  trim(gweconn%prim_exchange%name), &
372  ' to ', trim(gwfex%name), ' for GWE model ', &
373  trim(gwemodel%name)
374  gwfexidx = iex
375  if (gweconn%owns_exchange) then
376  gweconn%gweExchange%gwfsimvals => gwfex%simvals
377  call mem_checkin(gweconn%gweExchange%gwfsimvals, &
378  'GWFSIMVALS', gweconn%gweExchange%memoryPath, &
379  'SIMVALS', gwfex%memoryPath)
380  end if
381  !
382  !cdl link up mvt to mvr
383  if (gwfex%inmvr > 0) then
384  if (gweconn%owns_exchange) then
385  !cdl todo: check and make sure gweEx has mvt active
386  call gweconn%gweExchange%mvt%set_pointer_mvrbudobj( &
387  gwfex%mvr%budobj)
388  end if
389  end if
390  !
391  if (associated(gwfex%model2, gwfmodel)) gweconn%exgflowSign = -1
392  gweconn%gweInterfaceModel%fmi%flows_from_file = .false.
393  !
394  exit gwfloopexg
395  end if
396  end if
397  !
398  end do gwfloopexg
399  end if
400  !
401  if (gwfconnidx == -1 .and. gwfexidx == -1) then
402  ! none found, report
403  write (errmsg, '(/6a)') 'Missing GWF-GWF exchange when connecting GWE'// &
404  ' model ', trim(gwemodel%name), ' with exchange ', &
405  trim(gweconn%prim_exchange%name), ' to GWF model ', &
406  trim(gwfmodel%name)
407  call store_error(errmsg, terminate=.true.)
408  end if
409  !
410  end do gweloop
411  !
412  ! -- Return
413  return
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
integer(i4b) iout
file unit number for simulation output
Here is the call graph for this function:

◆ gwfgwe_cr()

subroutine, public gwfgweexchangemodule::gwfgwe_cr ( character(len=*), intent(in)  filename,
integer(i4b), intent(in)  id,
integer(i4b), intent(in)  m1_id,
integer(i4b), intent(in)  m2_id 
)

Definition at line 46 of file exg-gwfgwe.f90.

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  !
80  ! -- Return
81  return
integer(i4b), dimension(:), allocatable model_loc_idx
equals the local index into the basemodel list (-1 when not available)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ link_connections()

subroutine gwfgweexchangemodule::link_connections ( class(gwfgweexchangetype this,
class(gwegweconnectiontype), pointer  gweConn,
class(gwfgwfconnectiontype), pointer  gwfConn 
)
Parameters
thisthis exchange
gweconnGWE connection
gwfconnGWF connection

Definition at line 418 of file exg-gwfgwe.f90.

419  ! -- modules
421  ! -- dummy
422  class(GwfGweExchangeType) :: this !< this exchange
423  class(GweGweConnectionType), pointer :: gweConn !< GWE connection
424  class(GwfGwfConnectionType), pointer :: gwfConn !< GWF connection
425  !
426  !gweConn%exgflowja => gwfConn%exgflowja
427  if (gweconn%owns_exchange) then
428  gweconn%gweExchange%gwfsimvals => gwfconn%gwfExchange%simvals
429  call mem_checkin(gweconn%gweExchange%gwfsimvals, &
430  'GWFSIMVALS', gweconn%gweExchange%memoryPath, &
431  'SIMVALS', gwfconn%gwfExchange%memoryPath)
432  end if
433  !
434  !cdl link up mvt to mvr
435  if (gwfconn%gwfExchange%inmvr > 0) then
436  if (gweconn%owns_exchange) then
437  !cdl todo: check and make sure gweEx has mvt active
438  call gweconn%gweExchange%mvt%set_pointer_mvrbudobj( &
439  gwfconn%gwfExchange%mvr%budobj)
440  end if
441  end if
442  !
443  if (associated(gwfconn%gwfExchange%model2, gwfconn%owner)) then
444  gweconn%exgflowSign = -1
445  end if
446  !
447  ! fmi flows are not read from file
448  gweconn%gweInterfaceModel%fmi%flows_from_file = .false.
449  !
450  ! set concentration pointer for buoyancy
451  !call gwfConn%gwfInterfaceModel%buy%set_concentration_pointer( &
452  ! gweConn%gweModel%name, &
453  ! gweConn%conc, &
454  ! gweConn%icbound)
455  !
456  ! -- Return
457  return

◆ set_model_pointers()

subroutine gwfgweexchangemodule::set_model_pointers ( class(gwfgweexchangetype this)

Definition at line 86 of file exg-gwfgwe.f90.

87  ! -- dummy
88  class(GwfGweExchangeType) :: this
89  ! -- local
90  class(BaseModelType), pointer :: mb => null()
91  type(GwfModelType), pointer :: gwfmodel => null()
92  type(GweModelType), pointer :: gwemodel => null()
93  !
94  ! -- set gwfmodel
95  gwfmodel => null()
96  mb => getbasemodelfromlist(basemodellist, this%m1_idx)
97  select type (mb)
98  type is (gwfmodeltype)
99  gwfmodel => mb
100  end select
101  !
102  ! -- set gwemodel
103  gwemodel => null()
104  mb => getbasemodelfromlist(basemodellist, this%m2_idx)
105  select type (mb)
106  type is (gwemodeltype)
107  gwemodel => mb
108  end select
109  !
110  ! -- Verify that gwf model is of the correct type
111  if (.not. associated(gwfmodel)) then
112  write (errmsg, '(3a)') 'Problem with GWF-GWE exchange ', trim(this%name), &
113  '. Specified GWF Model does not appear to be of the correct type.'
114  call store_error(errmsg, terminate=.true.)
115  end if
116  !
117  ! -- Verify that gwe model is of the correct type
118  if (.not. associated(gwemodel)) then
119  write (errmsg, '(3a)') 'Problem with GWF-GWE exchange ', trim(this%name), &
120  '. Specified GWF Model does not appear to be of the correct type.'
121  call store_error(errmsg, terminate=.true.)
122  end if
123  !
124  ! -- Tell transport model fmi flows are not read from file
125  gwemodel%fmi%flows_from_file = .false.
126  !
127  ! -- Set a pointer to the GWF bndlist. This will allow the transport model
128  ! to look through the flow packages and establish a link to GWF flows
129  gwemodel%fmi%gwfbndlist => gwfmodel%bndlist
130  !
131  ! -- Return
132  return
Here is the call graph for this function: