MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
SpatialModelConnection.f90
Go to the documentation of this file.
2  use kindmodule, only: i4b, dp, lgp
4  use sparsemodule, only: sparsematrix
6  use csrutilsmodule, only: getcsrindex
20  use listmodule, only: listtype
21  use stlvecintmodule, only: stlvecint
24 
25  implicit none
26  private
27  public :: cast_as_smc
28  public :: add_smc_to_list
29  public :: get_smc_from_list
30 
31  !> Class to manage spatial connection of a model to one
32  !! or more models of the same type. Spatial connection here
33  !! means that the model domains (spatial discretization) are
34  !! adjacent and connected via DisConnExchangeType object(s).
35  !! The connection itself is a Numerical Exchange as well,
36  !! and part of a Numerical Solution providing the amat and rhs
37  !< values for the exchange.
39 
40  class(numericalmodeltype), pointer :: owner => null() !< the model whose connection this is
41  class(numericalmodeltype), pointer :: interface_model => null() !< the interface model
42  integer(I4B), pointer :: nr_connections => null() !< total nr. of connected cells (primary)
43 
44  class(disconnexchangetype), pointer :: prim_exchange => null() !< the exchange for which the interface model is created
45  logical(LGP) :: owns_exchange !< there are two connections (in serial) for an exchange,
46  !! one of them needs to manage/own the exchange (e.g. clean up)
47  type(stlvecint), pointer :: halo_models !< models that are potentially in the halo of this interface
48  type(stlvecint), pointer :: halo_exchanges !< exchanges that are potentially part of the halo of this interface (includes primary)
49  integer(I4B), pointer :: int_stencil_depth => null() !< size of the computational stencil for the interior
50  !! default = 1, xt3d = 2, ...
51  integer(I4B), pointer :: exg_stencil_depth => null() !< size of the computational stencil at the interface
52  !! default = 1, xt3d = 2, ...
53 
54  ! The following variables are equivalent to those in Numerical Solution:
55  integer(I4B), pointer :: neq => null() !< nr. of equations in matrix system
56  class(sparsematrixtype), pointer :: matrix => null() !< system matrix for the interface
57  real(dp), dimension(:), pointer, contiguous :: rhs => null() !< rhs of interface system
58  real(dp), dimension(:), pointer, contiguous :: x => null() !< dependent variable of interface system
59  integer(I4B), dimension(:), pointer, contiguous :: active => null() !< cell status (c.f. ibound) of interface system
60 
61  ! these are not in the memory manager
62  class(gridconnectiontype), pointer :: ig_builder => null() !< facility to build the interface grid connection structure
63  integer(I4B), dimension(:), pointer :: ipos_to_sln => null() !< mapping between position in the interface matrix and the solution matrix
64  type(listtype) :: iface_dist_vars !< list with distributed variables for this interface
65  type(interfacemaptype), pointer :: interface_map => null() !< a map of the interface into models and exchanges
66 
67  contains
68 
69  ! public
70  procedure, pass(this) :: spatialconnection_ctor
71  generic :: construct => spatialconnection_ctor
72 
73  ! partly overriding NumericalExchangeType:
74  procedure :: exg_df => spatialcon_df
75  procedure :: exg_ar => spatialcon_ar
76  procedure :: exg_ac => spatialcon_ac
77  procedure :: exg_mc => spatialcon_mc
78  procedure :: exg_cf => spatialcon_cf
79  procedure :: exg_fc => spatialcon_fc
80  procedure :: exg_da => spatialcon_da
81 
82  ! protected
83  procedure, pass(this) :: spatialcon_df
84  procedure, pass(this) :: spatialcon_ar
85  procedure, pass(this) :: spatialcon_ac
86  procedure, pass(this) :: spatialcon_cf
87  procedure, pass(this) :: spatialcon_fc
88  procedure, pass(this) :: spatialcon_da
89  procedure, pass(this) :: spatialcon_setmodelptrs
90  procedure, pass(this) :: spatialcon_connect
91  procedure, pass(this) :: validateconnection
92  procedure, pass(this) :: cfg_dv
93  procedure, pass(this) :: createmodelhalo
94 
95  ! private
96  procedure, private, pass(this) :: setupgridconnection
97  procedure, private, pass(this) :: getnrofconnections
98  procedure, private, pass(this) :: allocatescalars
99  procedure, private, pass(this) :: allocatearrays
100  procedure, private, pass(this) :: createcoefficientmatrix
101  procedure, private, pass(this) :: maskownerconnections
102  procedure, private, pass(this) :: addmodelneighbors
103 
105 
106 contains ! module procedures
107 
108  !> @brief Construct the spatial connection base
109  !!
110  !! This constructor is typically called from a derived class.
111  !<
112  subroutine spatialconnection_ctor(this, model, exchange, name)
113  class(spatialmodelconnectiontype) :: this !< the connection
114  class(numericalmodeltype), intent(in), pointer :: model !< the model that owns the connection
115  class(disconnexchangetype), intent(in), pointer :: exchange !< the primary exchange from which
116  !! the connection is created
117  character(len=*), intent(in) :: name !< the connection name (for memory management mostly)
118 
119  this%name = name
120  this%memoryPath = create_mem_path(this%name)
121 
122  this%owner => model
123  this%prim_exchange => exchange
124 
125  allocate (this%ig_builder)
126  allocate (this%halo_models)
127  allocate (this%halo_exchanges)
128  allocate (this%matrix)
129  call this%allocateScalars()
130 
131  this%int_stencil_depth = 1
132  this%exg_stencil_depth = 1
133  this%nr_connections = 0
134 
135  ! this should be set in derived ctor
136  this%interface_model => null()
137 
138  end subroutine spatialconnection_ctor
139 
140  !> @brief Find all models that might participate in this interface
141  !<
142  subroutine createmodelhalo(this)
143  class(spatialmodelconnectiontype) :: this !< this connection
144 
145  call this%halo_models%init()
146  call this%halo_exchanges%init()
147 
148  call this%addModelNeighbors(this%owner%id, virtual_exchange_list, &
149  this%exg_stencil_depth, .true.)
150 
151  end subroutine createmodelhalo
152 
153  !> @brief Add neighbors and nbrs-of-nbrs to the model tree
154  !<
155  recursive subroutine addmodelneighbors(this, model_id, &
156  virtual_exchanges, &
157  depth, is_root, mask)
159  class(spatialmodelconnectiontype) :: this !< this connection
160  integer(I4B) :: model_id !< the model (id) to add neighbors for
161  type(listtype) :: virtual_exchanges !< list with all virtual exchanges
162  integer(I4B), value :: depth !< the maximal number of exchanges between
163  logical(LGP) :: is_root !< true when called for neighbor from primary exchange
164  integer(I4B), optional :: mask !< don't add this one as a neighbor
165  ! local
166  integer(I4B) :: i, n
167  class(virtualexchangetype), pointer :: v_exg
168  integer(I4B) :: neighbor_id
169  integer(I4B) :: model_mask
170  type(stlvecint) :: models_at_depth !< model ids at a certain depth, to
171  !! recurse on for nbrs-of-nbrs search
172 
173  if (.not. present(mask)) then
174  model_mask = 0
175  else
176  model_mask = mask
177  end if
178 
179  call models_at_depth%init()
180 
181  if (is_root) then
182  ! first layer in the recursive search
183  call models_at_depth%push_back_unique(model_id)
184 
185  ! fetch primary neighbor
186  if (this%prim_exchange%v_model1%id == this%owner%id) then
187  neighbor_id = this%prim_exchange%v_model2%id
188  else
189  neighbor_id = this%prim_exchange%v_model1%id
190  end if
191  ! add
192  call models_at_depth%push_back_unique(neighbor_id)
193  call this%halo_models%push_back_unique(neighbor_id)
194  call this%halo_exchanges%push_back_unique(this%prim_exchange%id)
195  else
196  ! find all direct neighbors of the model and add them,
197  ! avoiding duplicates
198  do i = 1, virtual_exchanges%Count()
199  neighbor_id = -1
200  v_exg => get_virtual_exchange_from_list(virtual_exchanges, i)
201  if (v_exg%v_model1%id == model_id) then
202  neighbor_id = v_exg%v_model2%id
203  else if (v_exg%v_model2%id == model_id) then
204  neighbor_id = v_exg%v_model1%id
205  end if
206 
207  ! check if there is a neighbor, and it is not masked
208  ! (to prevent back-and-forth connections)
209  if (neighbor_id > 0) then
210  ! check if masked
211  if (neighbor_id == model_mask) cycle
212  call models_at_depth%push_back_unique(neighbor_id)
213  call this%halo_models%push_back_unique(neighbor_id)
214  call this%halo_exchanges%push_back_unique(v_exg%id)
215  end if
216  end do
217  end if
218 
219  depth = depth - 1
220  if (depth == 0) then
221  ! and we're done with this branch
222  call models_at_depth%destroy()
223  return
224  end if
225 
226  ! now recurse on the neighbors up to the specified depth
227  do n = 1, models_at_depth%size
228  call this%addModelNeighbors(models_at_depth%at(n), virtual_exchanges, &
229  depth, .false., model_id)
230  end do
231 
232  ! we're done with the tree
233  call models_at_depth%destroy()
234 
235  end subroutine addmodelneighbors
236 
237  !> @brief Define this connection, this is where the
238  !! discretization (DISU) for the interface model is
239  !< created!
240  subroutine spatialcon_df(this)
241  class(spatialmodelconnectiontype) :: this !< this connection
242  ! local
243  integer(I4B) :: i
244  class(virtualmodeltype), pointer :: v_model
245 
246  ! create the grid connection data structure
247  this%nr_connections = this%getNrOfConnections()
248  call this%ig_builder%construct(this%owner, &
249  this%nr_connections, &
250  this%name)
251  this%ig_builder%internalStencilDepth = this%int_stencil_depth
252  this%ig_builder%exchangeStencilDepth = this%exg_stencil_depth
253  this%ig_builder%haloExchanges => this%halo_exchanges
254  do i = 1, this%halo_models%size
255  v_model => get_virtual_model(this%halo_models%at(i))
256  call this%ig_builder%addToRegionalModels(v_model)
257  end do
258  call this%setupGridConnection()
259 
260  this%neq = this%ig_builder%nrOfCells
261  call this%allocateArrays()
262 
263  end subroutine spatialcon_df
264 
265  !> @brief Allocate the connection,
266  !<
267  subroutine spatialcon_ar(this)
268  class(spatialmodelconnectiontype) :: this !< this connection
269  ! local
270  integer(I4B) :: iface_idx, glob_idx
271  class(gridconnectiontype), pointer :: gc
272 
273  ! fill mapping to global index (which can be
274  ! done now because moffset is set in sln_df)
275  gc => this%ig_builder
276  do iface_idx = 1, gc%nrOfCells
277  glob_idx = gc%idxToGlobal(iface_idx)%index + &
278  gc%idxToGlobal(iface_idx)%v_model%moffset%get()
279  gc%idxToGlobalIdx(iface_idx) = glob_idx
280  end do
281 
282  end subroutine spatialcon_ar
283 
284  !> @brief set model pointers to connection
285  !<
286  subroutine spatialcon_setmodelptrs(this)
287  class(spatialmodelconnectiontype) :: this !< this connection
288 
289  ! point x, ibound, and rhs to connection
290  this%interface_model%x => this%x
291  call mem_checkin(this%interface_model%x, 'X', &
292  this%interface_model%memoryPath, 'X', &
293  this%memoryPath)
294  this%interface_model%rhs => this%rhs
295  call mem_checkin(this%interface_model%rhs, 'RHS', &
296  this%interface_model%memoryPath, 'RHS', &
297  this%memoryPath)
298  this%interface_model%ibound => this%active
299  call mem_checkin(this%interface_model%ibound, 'IBOUND', &
300  this%interface_model%memoryPath, 'IBOUND', &
301  this%memoryPath)
302 
303  end subroutine spatialcon_setmodelptrs
304 
305  !> @brief map interface model connections to our sparse matrix,
306  !< analogously to what happens in sln_connect.
307  subroutine spatialcon_connect(this)
308  class(spatialmodelconnectiontype) :: this !< this connection
309  ! local
310  type(sparsematrix) :: sparse
311  class(matrixbasetype), pointer :: matrix_base
312 
313  call sparse%init(this%neq, this%neq, 7)
314  call this%interface_model%model_ac(sparse)
315 
316  ! create amat from sparse
317  call this%createCoefficientMatrix(sparse)
318  call sparse%destroy()
319 
320  ! map connections
321  matrix_base => this%matrix
322  call this%interface_model%model_mc(matrix_base)
323  call this%maskOwnerConnections()
324 
325  end subroutine spatialcon_connect
326 
327  !> @brief Mask the owner's connections
328  !!
329  !! Determine which connections are handled by the interface model
330  !! (using the connections object in its discretization) and
331  !< set their mask to zero for the owning model.
332  subroutine maskownerconnections(this)
333  use csrutilsmodule, only: getcsrindex
334  class(spatialmodelconnectiontype) :: this !< the connection
335  ! local
336  integer(I4B) :: ipos, n, m, nloc, mloc, csr_idx
337  type(connectionstype), pointer :: conn
338 
339  ! set the mask on connections that are calculated by the interface model
340  conn => this%interface_model%dis%con
341  do n = 1, conn%nodes
342  ! only for connections internal to the owning model
343  if (.not. this%ig_builder%idxToGlobal(n)%v_model == this%owner) then
344  cycle
345  end if
346  nloc = this%ig_builder%idxToGlobal(n)%index
347 
348  do ipos = conn%ia(n) + 1, conn%ia(n + 1) - 1
349  m = conn%ja(ipos)
350  if (.not. this%ig_builder%idxToGlobal(m)%v_model == this%owner) then
351  cycle
352  end if
353  mloc = this%ig_builder%idxToGlobal(m)%index
354 
355  if (conn%mask(ipos) > 0) then
356  ! calculated by interface model, set local model's mask to zero
357  csr_idx = getcsrindex(nloc, mloc, this%owner%ia, this%owner%ja)
358  if (csr_idx == -1) then
359  ! this can only happen with periodic boundary conditions,
360  ! then there is no need to set the mask
361  if (this%ig_builder%isPeriodic(nloc, mloc)) cycle
362 
363  write (*, *) 'Error: cannot find cell connection in global system'
364  call ustop()
365  end if
366 
367  if (this%owner%dis%con%mask(csr_idx) > 0) then
368  call this%owner%dis%con%set_mask(csr_idx, 0)
369  else
370  ! edge case, this connection is already being calculated
371  ! so we ignore it here. This can happen in the overlap
372  ! between two different exchanges when a larger stencil
373  ! (XT3D) is applied.
374  call conn%set_mask(ipos, 0)
375  end if
376  end if
377  end do
378  end do
379 
380  end subroutine maskownerconnections
381 
382  !> @brief Add connections, handled by the interface model,
383  !< to the global system's sparse
384  subroutine spatialcon_ac(this, sparse)
385  class(spatialmodelconnectiontype) :: this !< this connection
386  type(sparsematrix), intent(inout) :: sparse !< sparse matrix to store the connections
387  ! local
388  integer(I4B) :: n, m, ipos
389  integer(I4B) :: icol_start, icol_end
390  integer(I4B) :: nglo, mglo
391 
392  do n = 1, this%neq
393  if (.not. this%ig_builder%idxToGlobal(n)%v_model == this%owner) then
394  ! only add connections for own model to global matrix
395  cycle
396  end if
397 
398  nglo = this%ig_builder%idxToGlobal(n)%index + &
399  this%ig_builder%idxToGlobal(n)%v_model%moffset%get()
400 
401  icol_start = this%matrix%get_first_col_pos(n)
402  icol_end = this%matrix%get_last_col_pos(n)
403  do ipos = icol_start, icol_end
404  m = this%matrix%get_column(ipos)
405  if (m == n) cycle
406  mglo = this%ig_builder%idxToGlobal(m)%index + &
407  this%ig_builder%idxToGlobal(m)%v_model%moffset%get()
408  call sparse%addconnection(nglo, mglo, 1)
409  end do
410 
411  end do
412 
413  end subroutine spatialcon_ac
414 
415  !> @brief Creates the mapping from the local system
416  !< matrix to the global one
417  subroutine spatialcon_mc(this, matrix_sln)
418  use simmodule, only: ustop
419  class(spatialmodelconnectiontype) :: this !< this connection
420  class(matrixbasetype), pointer :: matrix_sln !< global matrix
421  ! local
422  integer(I4B) :: i, m, n, mglo, nglo, ipos, ipos_sln
423  logical(LGP) :: is_owned
424 
425  allocate (this%ipos_to_sln(this%matrix%nja))
426  do i = 1, this%matrix%nja
427  this%ipos_to_sln(i) = -1
428  end do
429 
430  do n = 1, this%neq
431  is_owned = (this%ig_builder%idxToGlobal(n)%v_model == this%owner)
432  if (.not. is_owned) cycle
433 
434  do ipos = this%matrix%ia(n), this%matrix%ia(n + 1) - 1
435  m = this%matrix%ja(ipos)
436  nglo = this%ig_builder%idxToGlobal(n)%index + &
437  this%ig_builder%idxToGlobal(n)%v_model%moffset%get()
438  mglo = this%ig_builder%idxToGlobal(m)%index + &
439  this%ig_builder%idxToGlobal(m)%v_model%moffset%get()
440 
441  ipos_sln = matrix_sln%get_position(nglo, mglo)
442  if (ipos_sln == -1) then
443  ! this should not be possible
444  write (*, *) 'Error: cannot find cell connection in global system'
445  call ustop()
446  end if
447  this%ipos_to_sln(ipos) = ipos_sln
448 
449  end do
450  end do
451 
452  end subroutine spatialcon_mc
453 
454  !> @brief Calculate (or adjust) matrix coefficients,
455  !! in this case those which are determined or affected
456  !< by the connection of a GWF model with its neighbors
457  subroutine spatialcon_cf(this, kiter)
458  class(spatialmodelconnectiontype) :: this !< this connection
459  integer(I4B), intent(in) :: kiter !< the iteration counter
460  ! local
461  integer(I4B) :: i
462 
463  ! reset interface system
464  call this%matrix%zero_entries()
465  do i = 1, this%neq
466  this%rhs(i) = 0.0_dp
467  end do
468 
469  ! calculate the interface model
470  call this%interface_model%model_cf(kiter)
471 
472  end subroutine spatialcon_cf
473 
474  !> @brief Formulate coefficients from interface model
475  !<
476  subroutine spatialcon_fc(this, kiter, matrix_sln, rhs_sln, inwtflag)
477  class(spatialmodelconnectiontype) :: this !< this connection
478  integer(I4B), intent(in) :: kiter !< the iteration counter
479  class(matrixbasetype), pointer :: matrix_sln !< the system matrix
480  real(DP), dimension(:), intent(inout) :: rhs_sln !< global right-hand-side
481  integer(I4B), optional, intent(in) :: inwtflag !< newton-raphson flag
482  ! local
483  integer(I4B) :: n, nglo
484  integer(I4B) :: icol_start, icol_end, ipos
485  class(matrixbasetype), pointer :: matrix_base
486 
487  matrix_base => this%matrix
488  call this%interface_model%model_fc(kiter, matrix_base, inwtflag)
489 
490  ! map back to solution matrix
491  do n = 1, this%neq
492  ! We only need the coefficients for our own model
493  ! (i.e. rows in the matrix that belong to this%owner):
494  if (.not. this%ig_builder%idxToGlobal(n)%v_model == this%owner) then
495  cycle
496  end if
497 
498  nglo = this%ig_builder%idxToGlobal(n)%index + &
499  this%ig_builder%idxToGlobal(n)%v_model%moffset%get() - &
500  matrix_sln%get_row_offset()
501  rhs_sln(nglo) = rhs_sln(nglo) + this%rhs(n)
502 
503  icol_start = this%matrix%get_first_col_pos(n)
504  icol_end = this%matrix%get_last_col_pos(n)
505  do ipos = icol_start, icol_end
506  call matrix_sln%add_value_pos(this%ipos_to_sln(ipos), &
507  this%matrix%get_value_pos(ipos))
508  end do
509  end do
510 
511  end subroutine spatialcon_fc
512 
513  !> @brief Deallocation
514  !<
515  subroutine spatialcon_da(this)
516  class(spatialmodelconnectiontype) :: this !< this connection
517 
518  call mem_deallocate(this%neq)
519  call mem_deallocate(this%int_stencil_depth)
520  call mem_deallocate(this%exg_stencil_depth)
521  call mem_deallocate(this%nr_connections)
522 
523  call mem_deallocate(this%x)
524  call mem_deallocate(this%rhs)
525  call mem_deallocate(this%active)
526 
527  call this%halo_models%destroy()
528  call this%halo_exchanges%destroy()
529  deallocate (this%halo_models)
530  deallocate (this%halo_exchanges)
531  call this%matrix%destroy()
532  deallocate (this%matrix)
533 
534  call this%ig_builder%destroy()
535  call this%iface_dist_vars%Clear(destroy=.true.)
536  deallocate (this%ig_builder)
537  deallocate (this%interface_map)
538  deallocate (this%ipos_to_sln)
539 
540  end subroutine spatialcon_da
541 
542  !> @brief Creates the connection structure for the
543  !! interface grid, starting from primary exchanges,
544  !! then extending inward and outward, possibly across
545  !! model boundaries.
546  !<
547  subroutine setupgridconnection(this)
548  class(spatialmodelconnectiontype) :: this !< this connection
549  ! local
550 
551  ! connect cells from primary exchange
552  call this%ig_builder%connectPrimaryExchange(this%prim_exchange)
553 
554  ! now scan for nbr-of-nbrs and create final data structures
555  call this%ig_builder%extendConnection()
556 
557  ! construct the interface map
558  call this%ig_builder%buildInterfaceMap()
559  this%interface_map => this%ig_builder%interfaceMap
560 
561  end subroutine setupgridconnection
562 
563  !> @brief Allocation of scalars
564  !<
565  subroutine allocatescalars(this)
567  class(spatialmodelconnectiontype) :: this !< this connection
568 
569  call mem_allocate(this%neq, 'NEQ', this%memoryPath)
570  call mem_allocate(this%int_stencil_depth, 'INTSTDEPTH', this%memoryPath)
571  call mem_allocate(this%exg_stencil_depth, 'EXGSTDEPTH', this%memoryPath)
572  call mem_allocate(this%nr_connections, 'NROFCONNS', this%memoryPath)
573 
574  end subroutine allocatescalars
575 
576  !> @brief Allocation of arrays
577  !<
578  subroutine allocatearrays(this)
580  class(spatialmodelconnectiontype) :: this !< this connection
581  ! local
582  integer(I4B) :: i
583 
584  call mem_allocate(this%x, this%neq, 'X', this%memoryPath)
585  call mem_allocate(this%rhs, this%neq, 'RHS', this%memoryPath)
586  call mem_allocate(this%active, this%neq, 'IACTIVE', this%memoryPath)
587 
588  ! c.f. NumericalSolution
589  do i = 1, this%neq
590  this%x(i) = dzero
591  this%active(i) = 1 ! default is active
592  this%rhs(i) = dzero
593  end do
594 
595  end subroutine allocatearrays
596 
597  !> @brief Returns total nr. of primary connections
598  !<
599  function getnrofconnections(this) result(nrConns)
600  class(spatialmodelconnectiontype) :: this !< this connection
601  integer(I4B) :: nrConns
602  !local
603 
604  nrconns = this%prim_exchange%nexg
605 
606  end function getnrofconnections
607 
608  !> @brief Create connection's matrix (ia,ja,amat) from sparse
609  !<
610  subroutine createcoefficientmatrix(this, sparse)
611  use simmodule, only: ustop
612  class(spatialmodelconnectiontype) :: this !< this connection
613  type(sparsematrix), intent(inout) :: sparse !< the sparse matrix with the cell connections
614 
615  call sparse%sort()
616  call this%matrix%init(sparse, this%memoryPath)
617 
618  end subroutine createcoefficientmatrix
619 
620  !> @brief Validate this connection
621  !<
622  subroutine validateconnection(this)
623  class(spatialmodelconnectiontype) :: this !< this connection
624  ! local
625  class(disconnexchangetype), pointer :: conEx => null()
626  character(len=LINELENGTH) :: errmsg
627 
628  conex => this%prim_exchange
629  if (conex%ixt3d > 0) then
630  ! if XT3D, we need these angles:
631  if (conex%v_model1%con_ianglex%get() == 0) then
632  write (errmsg, '(a,a,a,a,a)') 'XT3D configured on the exchange ', &
633  trim(conex%name), ' but the discretization in model ', &
634  trim(conex%v_model1%name), ' has no ANGLDEGX specified'
635  call store_error(errmsg)
636  end if
637  if (conex%v_model2%con_ianglex%get() == 0) then
638  write (errmsg, '(a,a,a,a,a)') 'XT3D configured on the exchange ', &
639  trim(conex%name), ' but the discretization in model ', &
640  trim(conex%v_model2%name), ' has no ANGLDEGX specified'
641  call store_error(errmsg)
642  end if
643  end if
644 
645  if (count_errors() > 0) then
646  call ustop()
647  end if
648 
649  end subroutine validateconnection
650 
651  !> @brief Add a variable from the interface model to be
652  !! synchronized at the configured stages by copying from
653  !! the source memory in the models/exchanges that are part
654  !< of this interface.
655  subroutine cfg_dv(this, var_name, subcomp_name, map_type, &
656  sync_stages, exg_var_name)
657  class(spatialmodelconnectiontype) :: this !< this connection
658  character(len=*) :: var_name !< name of variable, e.g. "K11"
659  character(len=*) :: subcomp_name !< subcomponent, e.g. "NPF"
660  integer(I4B) :: map_type !< type of variable map
661  integer(I4B), dimension(:) :: sync_stages !< stages to sync
662  character(len=*), optional :: exg_var_name !< needed for exchange variables, e.g. SIMVALS
663  ! local
664  type(distvartype), pointer :: dist_var => null()
665  class(*), pointer :: obj
666 
667  if (.not. present(exg_var_name)) exg_var_name = ''
668 
669  allocate (dist_var)
670  dist_var%var_name = var_name
671  dist_var%subcomp_name = subcomp_name
672  dist_var%comp_name = this%interface_model%name
673  dist_var%map_type = map_type
674  dist_var%sync_stages = sync_stages
675  dist_var%exg_var_name = exg_var_name
676 
677  obj => dist_var
678  call this%iface_dist_vars%Add(obj)
679 
680  end subroutine cfg_dv
681 
682  !> @brief Cast to SpatialModelConnectionType
683  !<
684  function cast_as_smc(obj) result(res)
685  implicit none
686  class(*), pointer, intent(inout) :: obj !< object to be cast
687  class(spatialmodelconnectiontype), pointer :: res !< the instance of SpatialModelConnectionType
688  !
689  res => null()
690  if (.not. associated(obj)) return
691  !
692  select type (obj)
693  class is (spatialmodelconnectiontype)
694  res => obj
695  end select
696  end function cast_as_smc
697 
698  !> @brief Add connection to a list
699  !<
700  subroutine add_smc_to_list(list, conn)
701  implicit none
702  ! -- dummy
703  type(listtype), intent(inout) :: list !< the list
704  class(spatialmodelconnectiontype), pointer, intent(in) :: conn !< the connection
705  ! -- local
706  class(*), pointer :: obj
707  !
708  obj => conn
709  call list%Add(obj)
710  end subroutine add_smc_to_list
711 
712  !> @brief Get the connection from a list
713  !<
714  function get_smc_from_list(list, idx) result(res)
715  type(listtype), intent(inout) :: list !< the list
716  integer(I4B), intent(in) :: idx !< the index of the connection
717  class(spatialmodelconnectiontype), pointer :: res !< the returned connection
718 
719  ! local
720  class(*), pointer :: obj
721  obj => list%GetItem(idx)
722  res => cast_as_smc(obj)
723  end function get_smc_from_list
724 
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
integer(i4b) function, public getcsrindex(i, j, ia, ja)
Return index for element i,j in CSR storage,.
Definition: CsrUtils.f90:13
Refactoring issues towards parallel:
subroutine destroy(this)
This module defines variable data types.
Definition: kind.f90:8
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public ustop(stopmess, ioutlocal)
Stop the simulation.
Definition: Sim.f90:312
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
integer(i4b) function, public count_errors()
Return number of errors.
Definition: Sim.f90:59
subroutine validateconnection(this)
Validate this connection.
subroutine spatialcon_connect(this)
map interface model connections to our sparse matrix,
class(spatialmodelconnectiontype) function, pointer, public get_smc_from_list(list, idx)
Get the connection from a list.
subroutine spatialcon_df(this)
Define this connection, this is where the discretization (DISU) for the interface model is.
subroutine spatialcon_setmodelptrs(this)
set model pointers to connection
subroutine cfg_dv(this, var_name, subcomp_name, map_type, sync_stages, exg_var_name)
Add a variable from the interface model to be synchronized at the configured stages by copying from t...
subroutine spatialconnection_ctor(this, model, exchange, name)
Construct the spatial connection base.
class(spatialmodelconnectiontype) function, pointer, public cast_as_smc(obj)
Cast to SpatialModelConnectionType.
subroutine, public add_smc_to_list(list, conn)
Add connection to a list.
recursive subroutine addmodelneighbors(this, model_id, virtual_exchanges, depth, is_root, mask)
Add neighbors and nbrs-of-nbrs to the model tree.
subroutine createcoefficientmatrix(this, sparse)
Add connections, handled by the interface model,.
subroutine createmodelhalo(this)
Find all models that might participate in this interface.
subroutine maskownerconnections(this)
Mask the owner's connections.
subroutine spatialcon_ar(this)
Allocate the connection,.
type(listtype), public virtual_exchange_list
class(virtualexchangetype) function, pointer, public get_virtual_exchange_from_list(list, idx)
class(virtualexchangetype) function, pointer, public get_virtual_exchange(exg_id)
Returns a virtual exchange with the specified id.
Exchange based on connection between discretizations of DisBaseType. The data specifies the connectio...
This class is used to construct the connections object for the interface model's spatial discretizati...
A generic heterogeneous doubly-linked list.
Definition: List.f90:14
Class to manage spatial connection of a model to one or more models of the same type....
The Virtual Exchange is based on two Virtual Models and is therefore not always strictly local or rem...