MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
GhostNode.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: dp, i4b
4  use constantsmodule, only: linelength
9 
10  implicit none
11 
12  private
13  public :: ghostnodetype
14  public :: gnc_cr
15 
17  logical, pointer :: smgnc => null() ! single model gnc
18  logical, pointer :: implicit => null() ! lhs or rhs
19  logical, pointer :: i2kn => null() ! not used
20  integer(I4B), pointer :: nexg => null() ! number of gncs
21  integer(I4B), pointer :: numjs => null() ! number of connecting nodes
22  class(numericalmodeltype), pointer :: m1 => null() ! pointer to model 1
23  class(numericalmodeltype), pointer :: m2 => null() ! pointer to model 2
24  integer(I4B), dimension(:), pointer, contiguous :: nodem1 => null() ! array of nodes in model 1
25  integer(I4B), dimension(:), pointer, contiguous :: nodem2 => null() ! array of nodes in model 2
26  integer(I4B), dimension(:, :), pointer, contiguous :: nodesj => null() ! array of interpolation nodes
27  real(dp), dimension(:), pointer, contiguous :: cond => null() ! array of conductance
28  integer(I4B), dimension(:), pointer, contiguous :: idxglo => null() ! connection position in amat
29  integer(I4B), dimension(:), pointer, contiguous :: idxsymglo => null() ! symmetric position in amat
30  real(dp), dimension(:, :), pointer, contiguous :: alphasj => null() ! interpolation factors
31  integer(I4B), dimension(:), pointer, contiguous :: idiagn => null() ! amat diagonal position of n
32  integer(I4B), dimension(:), pointer, contiguous :: idiagm => null() ! amat diagonal position of m
33  integer(I4B), dimension(:, :), pointer, contiguous :: jposinrown => null() ! amat j position in row n
34  integer(I4B), dimension(:, :), pointer, contiguous :: jposinrowm => null() ! amat j position in row m
35 
36  contains
37 
38  procedure :: gnc_df
39  procedure :: gnc_ac
40  procedure :: gnc_mc
41  procedure, private :: gnc_fmsav
42  procedure :: gnc_fc
43  procedure :: gnc_fn
44  procedure :: gnc_cq
45  procedure :: gnc_ot
46  procedure :: gnc_da
47  procedure :: deltaqgnc
48  procedure :: allocate_scalars
49  procedure, private :: allocate_arrays
50  procedure, private :: read_options
51  procedure, private :: read_dimensions
52  procedure, private :: read_data
53  procedure, private :: nodeu_to_noder
54  end type ghostnodetype
55 
56 contains
57 
58  !> @brief Create new GNC exchange object
59  !<
60  subroutine gnc_cr(gncobj, name_parent, inunit, iout)
61  ! -- dummy
62  type(ghostnodetype), pointer, intent(inout) :: gncobj
63  character(len=*), intent(in) :: name_parent
64  integer(I4B), intent(in) :: inunit
65  integer(I4B), intent(in) :: iout
66  !
67  ! -- Allocate the gnc exchange object
68  allocate (gncobj)
69  !
70  ! -- create name and memory path. name_parent will either be model name or the
71  ! exchange name.
72  call gncobj%set_names(1, name_parent, 'GNC', 'GNC')
73  !
74  ! -- allocate scalars
75  call gncobj%allocate_scalars()
76  !
77  ! -- Set variables
78  gncobj%inunit = inunit
79  gncobj%iout = iout
80  !
81  ! -- Return
82  return
83  end subroutine gnc_cr
84 
85  !> @brief Initialize a gnc object
86  !<
87  subroutine gnc_df(this, m1, m2)
88  ! -- modules
91  use simvariablesmodule, only: errmsg
92  ! -- dummy
93  class(ghostnodetype) :: this
94  class(numericalmodeltype), target :: m1
95  class(numericalmodeltype), target, optional :: m2
96  !
97  ! -- Point or set attributes
98  this%m1 => m1
99  this%m2 => m1
100  !
101  ! -- If m2 is present, then GNC spans two models
102  if (present(m2)) then
103  this%m2 => m2
104  this%smgnc = .false.
105  end if
106  !
107  ! -- Initialize block parser
108  call this%parser%Initialize(this%inunit, this%iout)
109  !
110  ! -- read gnc options
111  call this%read_options()
112  !
113  ! -- read gnc dimensions
114  call this%read_dimensions()
115  !
116  ! -- allocate arrays
117  call this%allocate_arrays()
118  !
119  ! -- Allocate and read the gnc entries
120  call this%read_data()
121  !
122  ! -- Trap for implicit gnc but models are in different solutions
123  if (this%m1%idsoln /= this%m2%idsoln) then
124  if (this%implicit) then
125  write (errmsg, '(a)') 'GNC is implicit but models are in '// &
126  'different solutions.'
127  call store_error(errmsg)
128  call store_error_unit(this%inunit)
129  end if
130  end if
131  !
132  ! -- Return
133  return
134  end subroutine gnc_df
135 
136  !> @brief Single or Two-Model GNC Add Connections
137  !!
138  !! For implicit GNC, expand the sparse solution matrix
139  !<
140  subroutine gnc_ac(this, sparse)
141  ! -- modules
142  use sparsemodule, only: sparsematrix
143  ! -- dummy
144  class(ghostnodetype) :: this
145  type(sparsematrix), intent(inout) :: sparse
146  ! -- local
147  integer(I4B) :: ignc, jidx, noden, nodem, nodej
148  !
149  ! -- Expand the sparse matrix for ghost node connections. No need to add
150  ! connection between n and m as they must be connected some other way
151  ! that will calculate the conductance.
152  if (this%implicit) then
153  do ignc = 1, this%nexg
154  noden = this%nodem1(ignc) + this%m1%moffset
155  nodem = this%nodem2(ignc) + this%m2%moffset
156  jloop: do jidx = 1, this%numjs
157  nodej = this%nodesj(jidx, ignc)
158  if (nodej == 0) cycle
159  nodej = nodej + this%m1%moffset
160  call sparse%addconnection(nodem, nodej, 1)
161  call sparse%addconnection(nodej, nodem, 1)
162  call sparse%addconnection(noden, nodej, 1)
163  call sparse%addconnection(nodej, noden, 1)
164  end do jloop
165  end do
166  end if
167  !
168  ! -- Return
169  return
170  end subroutine gnc_ac
171 
172  !> @brief Single or Two-Model GNC Map Connections
173  !!
174  !! Fill the following mapping arrays:
175  !! this%idiagn, this%idiagm (diagonal positions in solution amat)
176  !! this%idxglo (nm connection in solution amat)
177  !! this%idxsymglo (mn connection in solution amat)
178  !! this%jposinrown (position of j in row n of solution amat)
179  !! this%jposinrowm (position of j in row m of solution amat)
180  !<
181  subroutine gnc_mc(this, matrix_sln)
182  ! -- modules
184  use simvariablesmodule, only: errmsg
185  ! -- dummy
186  class(ghostnodetype) :: this
187  class(matrixbasetype), pointer :: matrix_sln
188  ! -- local
189  integer(I4B) :: noden, nodem, ipos, ignc, jidx, nodej
190  ! -- formats
191  character(len=*), parameter :: fmterr = &
192  "('GHOST NODE ERROR. Cell ', i0, ' in model ', a, &
193  &' is not connected to cell ', i0, ' in model ', a)"
194  !
195  ! -- Find the location of Cnm in the global solution and store it in
196  ! this%idxglo
197  do ignc = 1, this%nexg
198  noden = this%nodem1(ignc) + this%m1%moffset
199  nodem = this%nodem2(ignc) + this%m2%moffset
200  !
201  ! -- store diagonal positions in idiagn and idiagm
202  this%idiagn(ignc) = matrix_sln%get_position_diag(noden)
203  this%idiagm(ignc) = matrix_sln%get_position_diag(nodem)
204  !if(this%implicit) then
205  ! this%idiagn(ignc) = iasln(noden)
206  ! this%idiagm(ignc) = iasln(nodem)
207  !endif
208  !
209  ! -- find location of m in row n of global solution, and v.v.
210  this%idxglo(ignc) = matrix_sln%get_position(noden, nodem)
211  this%idxsymglo(ignc) = matrix_sln%get_position(nodem, noden)
212  !
213  ! -- Check to make sure idxglo is set
214  if (this%idxglo(ignc) == -1) then
215  write (errmsg, fmterr) this%nodem1(ignc), trim(this%m1%name), &
216  this%nodem2(ignc), trim(this%m2%name)
217  call store_error(errmsg)
218  end if
219  !
220  end do
221  !
222  ! -- Stop if errors
223  if (count_errors() > 0) then
224  call store_error_unit(this%inunit)
225  end if
226  !
227  ! -- find locations of j in rows n and row m of global solution
228  if (this%implicit) then
229  do ignc = 1, this%nexg
230  noden = this%nodem1(ignc) + this%m1%moffset
231  nodem = this%nodem2(ignc) + this%m2%moffset
232  !
233  do jidx = 1, this%numjs
234  nodej = this%nodesj(jidx, ignc)
235  if (nodej > 0) nodej = nodej + this%m1%moffset
236  !
237  ! -- search for nodej in row n, unless it is 0
238  if (nodej == 0) then
239  ipos = 0
240  this%jposinrown(jidx, ignc) = ipos
241  else
242  this%jposinrown(jidx, ignc) = matrix_sln%get_position(noden, nodej)
243  end if
244  !
245  ! -- search for nodej in row m
246  if (nodej == 0) then
247  ipos = 0
248  this%jposinrowm(jidx, ignc) = ipos
249  else
250  this%jposinrowm(jidx, ignc) = matrix_sln%get_position(nodem, nodej)
251  end if
252  end do
253  end do
254  end if
255  !
256  ! -- Return
257  return
258  end subroutine gnc_mc
259 
260  !> @brief Store the n-m Picard conductance in cond prior to the Newton terms
261  !! terms being added
262  !<
263  subroutine gnc_fmsav(this, kiter, matrix)
264  ! -- modules
265  use constantsmodule, only: dzero
266  ! -- dummy
267  class(ghostnodetype) :: this
268  integer(I4B), intent(in) :: kiter
269  class(matrixbasetype), pointer :: matrix
270  ! -- local
271  integer(I4B) :: ignc, ipos
272  real(DP) :: cond
273  !
274  ! -- An ipos value of zero indicates that noden is not connected to
275  ! nodem, and therefore the conductance is zero.
276  gncloop: do ignc = 1, this%nexg
277  ipos = this%idxglo(ignc)
278  if (ipos > 0) then
279  cond = matrix%get_value_pos(ipos)
280  else
281  cond = dzero
282  end if
283  this%cond(ignc) = cond
284  end do gncloop
285  !
286  ! -- Return
287  return
288  end subroutine gnc_fmsav
289 
290  !> @brief Fill matrix terms
291  !!
292  !! Add the GNC terms to the solution matrix or model rhs depending on whether
293  !! whether GNC is implicit or explicit
294  !<
295  subroutine gnc_fc(this, kiter, matrix)
296  ! -- modules
297  use constantsmodule, only: dzero
298  ! -- dummy
299  class(ghostnodetype) :: this
300  integer(I4B), intent(in) :: kiter
301  class(matrixbasetype), pointer :: matrix
302  ! -- local
303  integer(I4B) :: ignc, j, noden, nodem, ipos, jidx, iposjn, iposjm
304  real(DP) :: cond, alpha, aterm, rterm
305  !
306  ! -- If this is a single model gnc (not an exchange across models), then
307  ! pull conductances out of the system matrix and store them in this%cond
308  if (this%smgnc) call this%gnc_fmsav(kiter, matrix)
309  !
310  ! -- Add gnc terms to rhs or to the matrix depending on whether gnc is implicit
311  ! or explicit
312  gncloop: do ignc = 1, this%nexg
313  noden = this%nodem1(ignc)
314  nodem = this%nodem2(ignc)
315  if (this%m1%ibound(noden) == 0 .or. &
316  this%m2%ibound(nodem) == 0) cycle gncloop
317  ipos = this%idxglo(ignc)
318  cond = this%cond(ignc)
319  jloop: do jidx = 1, this%numjs
320  j = this%nodesj(jidx, ignc)
321  if (j == 0) cycle
322  alpha = this%alphasj(jidx, ignc)
323  if (alpha == dzero) cycle
324  aterm = alpha * cond
325  if (this%implicit) then
326  iposjn = this%jposinrown(jidx, ignc)
327  iposjm = this%jposinrowm(jidx, ignc)
328  call matrix%add_value_pos(this%idiagn(ignc), aterm)
329  call matrix%add_value_pos(iposjn, -aterm)
330  call matrix%add_value_pos(this%idxsymglo(ignc), -aterm)
331  call matrix%add_value_pos(iposjm, aterm)
332  else
333  rterm = aterm * (this%m1%x(noden) - this%m1%x(j))
334  this%m1%rhs(noden) = this%m1%rhs(noden) - rterm
335  this%m2%rhs(nodem) = this%m2%rhs(nodem) + rterm
336  end if
337  end do jloop
338  end do gncloop
339  !
340  ! -- Return
341  return
342  end subroutine gnc_fc
343 
344  !> @brief Fill GNC Newton terms
345  !!
346  !! Required arguments:
347  !! kiter : outer iteration number
348  !! matrix_sln: the solution matrix
349  !! condsat is of size(njas) if single model, otherwise nexg
350  !!
351  !! Optional arguments:
352  !! ihc_opt : an optional vector of size(nexg), which contains a horizontal
353  !! connection code (0=vertical, 1=horizontal, 2=vertically staggered)
354  !! ivarcv_opt : variable vertical conductance flag (default is 0)
355  !! ictm1_opt : icelltype for model 1 integer vector (default is 1)
356  !! ictm2_opt : icelltype for model 2 integer vector (default is 1)
357  !<
358  subroutine gnc_fn(this, kiter, matrix_sln, condsat, ihc_opt, &
359  ivarcv_opt, ictm1_opt, ictm2_opt)
360  ! -- modules
361  use constantsmodule, only: dzero
363  ! -- dummy
364  class(ghostnodetype) :: this
365  integer(I4B) :: kiter
366  class(matrixbasetype), pointer :: matrix_sln
367  real(DP), dimension(:), intent(in) :: condsat
368  integer(I4B), dimension(:), optional :: ihc_opt
369  integer(I4B), optional :: ivarcv_opt
370  integer(I4B), dimension(:), optional :: ictm1_opt
371  integer(I4B), dimension(:), optional :: ictm2_opt
372  ! -- local
373  integer(I4B) :: ignc, jidx, ipos, isympos, ihc, ivarcv
374  integer(I4B) :: nodej, noden, nodem
375  integer(I4B) :: iups, ictup
376  real(DP) :: csat, alpha, consterm, term, derv
377  real(DP) :: xup, topup, botup
378  !
379  ! -- Set the ivarcv to indicate whether or not the vertical conductance
380  ! is a function of water table
381  ivarcv = 0
382  if (present(ivarcv_opt)) ivarcv = ivarcv_opt
383  !
384  gncloop: do ignc = 1, this%nexg
385  noden = this%nodem1(ignc)
386  nodem = this%nodem2(ignc)
387  if (this%m1%ibound(noden) == 0 .or. &
388  this%m2%ibound(nodem) == 0) cycle gncloop
389  !
390  ! -- Assign variables depending on whether single model gnc or exchange
391  ! gnc
392  if (this%smgnc) then
393  ipos = this%m1%dis%con%getjaindex(noden, nodem)
394  isympos = this%m1%dis%con%jas(ipos)
395  ihc = this%m1%dis%con%ihc(isympos)
396  csat = condsat(isympos)
397  else
398  ihc = ihc_opt(ignc)
399  csat = condsat(ignc)
400  end if
401  !
402  ! If vertical connection and not variable cv, then cycle
403  if (ihc == 0 .and. ivarcv == 0) cycle
404  !
405  ! determine upstream node (0 is noden, 1 is nodem)
406  iups = 0
407  if (this%m2%x(nodem) > this%m1%x(noden)) iups = 1
408  !
409  ! -- Set the upstream top and bot, and then recalculate for a
410  ! vertically staggered horizontal connection
411  if (iups == 0) then
412  topup = this%m1%dis%top(noden)
413  botup = this%m1%dis%bot(noden)
414  ictup = 1
415  if (present(ictm1_opt)) ictup = ictm1_opt(noden)
416  xup = this%m1%x(noden)
417  else
418  topup = this%m2%dis%top(nodem)
419  botup = this%m2%dis%bot(nodem)
420  ictup = 1
421  if (present(ictm2_opt)) ictup = ictm2_opt(nodem)
422  xup = this%m2%x(nodem)
423  end if
424  !
425  ! -- No newton terms if upstream cell is confined
426  if (ictup == 0) cycle
427  !
428  ! -- Handle vertically staggered horizontal connection
429  if (ihc == 2) then
430  topup = min(this%m1%dis%top(noden), this%m2%dis%top(nodem))
431  botup = max(this%m1%dis%bot(noden), this%m2%dis%bot(nodem))
432  end if
433  !
434  ! -- Process each contributing node
435  jloop: do jidx = 1, this%numjs
436  nodej = this%nodesj(jidx, ignc)
437  if (nodej == 0) cycle
438  if (this%m1%ibound(nodej) == 0) cycle
439  alpha = this%alphasj(jidx, ignc)
440  if (alpha == dzero) cycle
441  consterm = csat * alpha * (this%m1%x(noden) - this%m1%x(nodej))
442  derv = squadraticsaturationderivative(topup, botup, xup)
443  term = consterm * derv
444  if (iups == 0) then
445  call matrix_sln%add_value_pos(this%idiagn(ignc), term)
446  if (this%m2%ibound(nodem) > 0) then
447  call matrix_sln%add_value_pos(this%idxsymglo(ignc), -term)
448  end if
449  this%m1%rhs(noden) = this%m1%rhs(noden) + term * this%m1%x(noden)
450  this%m2%rhs(nodem) = this%m2%rhs(nodem) - term * this%m1%x(noden)
451  else
452  call matrix_sln%add_value_pos(this%idiagm(ignc), -term)
453  if (this%m1%ibound(noden) > 0) then
454  call matrix_sln%add_value_pos(this%idxglo(ignc), term)
455  end if
456  this%m1%rhs(noden) = this%m1%rhs(noden) + term * this%m2%x(nodem)
457  this%m2%rhs(nodem) = this%m2%rhs(nodem) - term * this%m2%x(nodem)
458  end if
459  end do jloop
460  end do gncloop
461  !
462  ! -- Return
463  return
464  end subroutine gnc_fn
465 
466  !> @brief Single Model GNC Output
467  !!
468  !! Output GNC deltaQgnc values
469  !<
470  subroutine gnc_ot(this, ibudfl)
471  ! -- dummy
472  class(ghostnodetype) :: this
473  integer(I4B), intent(in) :: ibudfl
474  ! -- local
475  integer(I4B) :: ignc
476  real(DP) :: deltaQgnc
477  character(len=LINELENGTH) :: nodenstr, nodemstr
478  ! -- format
479  character(len=*), parameter :: fmtgnc = "(i10, 2a10, 2(1pg15.6))"
480  !
481  ! -- Process each gnc and output deltaQgnc
482  if (ibudfl /= 0 .and. this%iprflow /= 0) then
483  write (this%iout, '(//, a)') 'GHOST NODE CORRECTION RESULTS'
484  write (this%iout, '(3a10, 2a15)') 'GNC NUM', 'NODEN', 'NODEM', &
485  'DELTAQGNC', 'CONDNM'
486  do ignc = 1, this%nexg
487  deltaqgnc = this%deltaQgnc(ignc)
488  call this%m1%dis%noder_to_string(this%nodem1(ignc), nodenstr)
489  call this%m2%dis%noder_to_string(this%nodem2(ignc), nodemstr)
490  write (this%iout, fmtgnc) ignc, trim(adjustl(nodenstr)), &
491  trim(adjustl(nodemstr)), &
492  deltaqgnc, this%cond(ignc)
493  end do
494  end if
495  !
496  ! -- Return
497  return
498  end subroutine gnc_ot
499 
500  !> @brief Add GNC to flowja
501  !<
502  subroutine gnc_cq(this, flowja)
503  ! -- dummy
504  class(ghostnodetype) :: this
505  real(DP), dimension(:), intent(inout) :: flowja
506  ! -- local
507  integer(I4B) :: ignc, n1, n2, ipos, isympos
508  real(DP) :: deltaQgnc
509  !
510  ! -- go through each gnc and add deltagnc to flowja
511  do ignc = 1, this%nexg
512  !
513  ! -- calculate correction term between n1 and n2 connection
514  n1 = this%nodem1(ignc)
515  n2 = this%nodem2(ignc)
516  deltaqgnc = this%deltaQgnc(ignc)
517  !
518  ! -- find the positions of this connection in the csr array
519  ipos = this%m1%dis%con%getjaindex(n1, n2)
520  isympos = this%m1%dis%con%isym(ipos)
521  !
522  ! -- add/subtract the corrections
523  flowja(ipos) = flowja(ipos) + deltaqgnc
524  flowja(isympos) = flowja(isympos) - deltaqgnc
525  !
526  end do
527  !
528  ! -- Return
529  return
530  end subroutine gnc_cq
531 
532  !> @brief Single Model deltaQgnc (ghost node correction flux)
533  !!
534  !! Calculate the deltaQgnc value for any GNC in the GNC list
535  !<
536  function deltaqgnc(this, ignc)
537  ! -- modules
538  use constantsmodule, only: dzero
539  ! -- Return
540  real(dp) :: deltaqgnc
541  ! -- dummy
542  class(ghostnodetype) :: this
543  integer(I4B), intent(in) :: ignc
544  ! -- local
545  integer(I4B) :: noden, nodem, nodej, jidx
546  real(dp) :: sigalj, alpha, hd, aterm, cond
547  !
548  ! -- initialize values
549  deltaqgnc = dzero
550  sigalj = dzero
551  hd = dzero
552  noden = this%nodem1(ignc)
553  nodem = this%nodem2(ignc)
554  !
555  ! -- calculate deltaQgnc
556  if (this%m1%ibound(noden) /= 0 .and. this%m2%ibound(nodem) /= 0) then
557  jloop: do jidx = 1, this%numjs
558  nodej = this%nodesj(jidx, ignc)
559  if (nodej == 0) cycle jloop
560  if (this%m1%ibound(nodej) == 0) cycle jloop
561  alpha = this%alphasj(jidx, ignc)
562  sigalj = sigalj + alpha
563  hd = hd + alpha * this%m1%x(nodej)
564  end do jloop
565  aterm = sigalj * this%m1%x(noden) - hd
566  cond = this%cond(ignc)
567  deltaqgnc = aterm * cond
568  end if
569  !
570  ! -- Return
571  return
572  end function deltaqgnc
573 
574  !> @brief Allocate gnc scalar variables
575  !<
576  subroutine allocate_scalars(this)
577  ! -- modules
579  ! -- dummy
580  class(ghostnodetype) :: this
581  !
582  ! -- allocate scalars in NumericalPackageType
583  call this%NumericalPackageType%allocate_scalars()
584  !
585  call mem_allocate(this%smgnc, 'SMGNC', this%memoryPath)
586  call mem_allocate(this%implicit, 'IMPLICIT', this%memoryPath)
587  call mem_allocate(this%i2kn, 'I2KN', this%memoryPath)
588  call mem_allocate(this%nexg, 'NEXG', this%memoryPath)
589  call mem_allocate(this%numjs, 'NUMJS', this%memoryPath)
590  !
591  ! -- Initialize values
592  this%smgnc = .true.
593  this%implicit = .true.
594  this%i2kn = .false.
595  this%nexg = 0
596  this%numjs = 0
597  !
598  ! -- Return
599  return
600  end subroutine allocate_scalars
601 
602  !> @brief Allocate gnc scalar variables
603  !<
604  subroutine allocate_arrays(this)
605  ! -- modules
607  ! -- dummy
608  class(ghostnodetype) :: this
609  !
610  ! -- allocate memory for arrays
611  call mem_allocate(this%nodem1, this%nexg, 'NODEM1', this%memoryPath)
612  call mem_allocate(this%nodem2, this%nexg, 'NODEM2', this%memoryPath)
613  call mem_allocate(this%nodesj, this%numjs, this%nexg, 'NODESJ', &
614  this%memoryPath)
615  call mem_allocate(this%alphasj, this%numjs, this%nexg, 'ALPHASJ', &
616  this%memoryPath)
617  call mem_allocate(this%cond, this%nexg, 'COND', this%memoryPath)
618  call mem_allocate(this%idxglo, this%nexg, 'IDXGLO', this%memoryPath)
619  call mem_allocate(this%idiagn, this%nexg, 'IDIAGN', this%memoryPath)
620  call mem_allocate(this%idiagm, this%nexg, 'IDIAGM', this%memoryPath)
621  call mem_allocate(this%idxsymglo, this%nexg, 'IDXSYMGLO', this%memoryPath)
622  if (this%implicit) then
623  call mem_allocate(this%jposinrown, this%numjs, this%nexg, 'JPOSINROWN', &
624  this%memoryPath)
625  call mem_allocate(this%jposinrowm, this%numjs, this%nexg, 'JPOSINROWM', &
626  this%memoryPath)
627  else
628  call mem_allocate(this%jposinrown, 0, 0, 'JPOSINROWN', this%memoryPath)
629  call mem_allocate(this%jposinrowm, 0, 0, 'JPOSINROWM', this%memoryPath)
630  end if
631  !
632  ! -- Return
633  return
634  end subroutine allocate_arrays
635 
636  !> @brief Deallocate memory
637  !<
638  subroutine gnc_da(this)
639  ! -- modules
641  ! -- dummy
642  class(ghostnodetype) :: this
643  !
644  call mem_deallocate(this%smgnc)
645  call mem_deallocate(this%implicit)
646  call mem_deallocate(this%i2kn)
647  call mem_deallocate(this%nexg)
648  call mem_deallocate(this%numjs)
649  !
650  ! -- Arrays
651  if (this%inunit > 0) then
652  call mem_deallocate(this%nodem1)
653  call mem_deallocate(this%nodem2)
654  call mem_deallocate(this%nodesj)
655  call mem_deallocate(this%alphasj)
656  call mem_deallocate(this%cond)
657  call mem_deallocate(this%idxglo)
658  call mem_deallocate(this%idiagn)
659  call mem_deallocate(this%idiagm)
660  call mem_deallocate(this%idxsymglo)
661  call mem_deallocate(this%jposinrown)
662  call mem_deallocate(this%jposinrowm)
663  end if
664  !
665  ! -- deallocate NumericalPackageType
666  call this%NumericalPackageType%da()
667  !
668  ! -- Return
669  return
670  end subroutine gnc_da
671 
672  !> @brief Read a gnc options block
673  !!
674  !! Read options from input file
675  !<
676  subroutine read_options(this)
677  ! -- modules
678  use simmodule, only: store_error
679  use simvariablesmodule, only: errmsg
680  ! -- dummy
681  class(ghostnodetype) :: this
682  ! -- local
683  character(len=LINELENGTH) :: keyword
684  integer(I4B) :: ierr
685  logical :: isfound, endOfBlock
686  !
687  ! -- get options block
688  call this%parser%GetBlock('OPTIONS', isfound, ierr, &
689  supportopenclose=.true., blockrequired=.false.)
690  !
691  ! -- parse options block if detected
692  if (isfound) then
693  write (this%iout, '(1x,a)') 'PROCESSING GNC OPTIONS'
694  do
695  call this%parser%GetNextLine(endofblock)
696  if (endofblock) exit
697  call this%parser%GetStringCaps(keyword)
698  select case (keyword)
699  case ('PRINT_INPUT')
700  this%iprpak = 1
701  write (this%iout, '(4x,a)') &
702  'THE LIST OF GHOST-NODE CORRECTIONS WILL BE PRINTED.'
703  case ('PRINT_FLOWS')
704  this%iprflow = 1
705  write (this%iout, '(4x,a)') &
706  'DELTAQGNC VALUES WILL BE PRINTED TO THE LIST FILE.'
707  case ('I2KN')
708  this%i2kn = .true.
709  write (this%iout, '(4x,a)') &
710  'SECOND ORDER CORRECTION WILL BE APPLIED.'
711  case ('EXPLICIT')
712  this%implicit = .false.
713  write (this%iout, '(4x,a)') 'GHOST NODE CORRECTION IS EXPLICIT.'
714  case default
715  write (errmsg, '(a,a)') 'Unknown GNC option: ', &
716  trim(keyword)
717  call store_error(errmsg)
718  call this%parser%StoreErrorUnit()
719  end select
720  end do
721  write (this%iout, '(1x,a)') 'END OF GNC OPTIONS'
722  end if
723  !
724  ! -- Set the iasym flag if the correction is implicit
725  if (this%implicit) this%iasym = 1
726  !
727  ! -- Return
728  return
729  end subroutine read_options
730 
731  !> @brief Single Model GNC Read Dimensions
732  !!
733  !! Read dimensions (size of gnc list) from input file
734  !<
735  subroutine read_dimensions(this)
736  ! -- modules
737  use simmodule, only: store_error
738  use simvariablesmodule, only: errmsg
739  ! -- dummy
740  class(ghostnodetype) :: this
741  ! -- local
742  character(len=LINELENGTH) :: keyword
743  integer(I4B) :: ierr
744  logical :: isfound, endOfBlock
745  !
746  ! -- get options block
747  call this%parser%GetBlock('DIMENSIONS', isfound, ierr, &
748  supportopenclose=.true.)
749  !
750  ! -- parse options block if detected
751  if (isfound) then
752  write (this%iout, '(1x,a)') 'PROCESSING GNC DIMENSIONS'
753  do
754  call this%parser%GetNextLine(endofblock)
755  if (endofblock) exit
756  call this%parser%GetStringCaps(keyword)
757  select case (keyword)
758  case ('NUMGNC')
759  this%nexg = this%parser%GetInteger()
760  write (this%iout, '(4x,a,i7)') 'NUMGNC = ', this%nexg
761  case ('NUMALPHAJ')
762  this%numjs = this%parser%GetInteger()
763  write (this%iout, '(4x,a,i7)') 'NUMAPHAJ = ', this%numjs
764  case default
765  write (errmsg, '(a,a)') 'Unknown GNC dimension: ', &
766  trim(keyword)
767  call store_error(errmsg)
768  call this%parser%StoreErrorUnit()
769  end select
770  end do
771  write (this%iout, '(1x,a)') 'END OF GNC DIMENSIONS'
772  else
773  call store_error('Required DIMENSIONS block not found.', terminate=.true.)
774  end if
775  !
776  ! -- Return
777  return
778  end subroutine read_dimensions
779 
780  !> @brief Read a GNCDATA block
781  !!
782  !! Read list of GNCs from input file
783  !<
784  subroutine read_data(this)
785  ! -- modules
787  use simvariablesmodule, only: errmsg
788  ! -- dummy
789  class(ghostnodetype) :: this
790  ! -- local
791  character(len=LINELENGTH) :: line, nodestr, fmtgnc, cellid, &
792  cellidm, cellidn
793  integer(I4B) :: lloc, ierr, ival
794  integer(I4B) :: ignc, jidx, nodeun, nodeum, nerr
795  integer(I4B), dimension(:), allocatable :: nodesuj
796  logical :: isfound, endOfBlock
797  !
798  ! -- Construct the fmtgnc format
799  write (fmtgnc, '("(2i10,",i0,"i10,",i0, "(1pg15.6))")') this%numjs, &
800  this%numjs
801  !
802  ! -- Allocate the temporary nodesuj, which stores the user-based nodej
803  ! node numbers
804  allocate (nodesuj(this%numjs))
805  !
806  ! -- get GNCDATA block
807  call this%parser%GetBlock('GNCDATA', isfound, ierr, supportopenclose=.true.)
808  !
809  ! -- process GNC data
810  if (isfound) then
811  write (this%iout, '(1x,a)') 'PROCESSING GNCDATA'
812  do ignc = 1, this%nexg
813  call this%parser%GetNextLine(endofblock)
814  if (endofblock) exit
815  call this%parser%GetCurrentLine(line)
816  lloc = 1
817  !
818  ! -- cellidn (read as cellid and convert to user node)
819  call this%parser%GetCellid(this%m1%dis%ndim, cellidn)
820  nodeun = this%m1%dis%nodeu_from_cellid(cellidn, this%parser%iuactive, &
821  this%iout)
822  !
823  ! -- convert user node to reduced node number
824  call this%nodeu_to_noder(nodeun, this%nodem1(ignc), this%m1)
825  !
826  ! -- cellidm (read as cellid and convert to user node)
827  call this%parser%GetCellid(this%m2%dis%ndim, cellidm)
828  nodeum = this%m2%dis%nodeu_from_cellid(cellidm, this%parser%iuactive, &
829  this%iout)
830  !
831  ! -- convert user node to reduced node number
832  call this%nodeu_to_noder(nodeum, this%nodem2(ignc), this%m2)
833  !
834  ! -- cellidsj (read as cellid)
835  do jidx = 1, this%numjs
836  ! read cellidj as cellid of model 1
837  call this%parser%GetCellid(this%m1%dis%ndim, cellid)
838  ival = this%m1%dis%nodeu_from_cellid(cellid, this%parser%iuactive, &
839  this%iout, allow_zero=.true.)
840  nodesuj(jidx) = ival
841  if (ival > 0) then
842  call this%nodeu_to_noder(ival, this%nodesj(jidx, ignc), this%m1)
843  else
844  this%nodesj(jidx, ignc) = 0
845  end if
846  end do
847  !
848  ! -- alphaj
849  do jidx = 1, this%numjs
850  this%alphasj(jidx, ignc) = this%parser%GetDouble()
851  end do
852  !
853  ! -- Echo if requested
854  if (this%iprpak /= 0) &
855  write (this%iout, fmtgnc) nodeun, nodeum, &
856  (nodesuj(jidx), jidx=1, this%numjs), &
857  (this%alphasj(jidx, ignc), jidx=1, this%numjs)
858  !
859  ! -- Check to see if noden is outside of active domain
860  if (this%nodem1(ignc) <= 0) then
861  call this%m1%dis%nodeu_to_string(nodeun, nodestr)
862  write (errmsg, *) &
863  trim(adjustl(this%m1%name))// &
864  ' Cell is outside active grid domain: '// &
865  trim(adjustl(nodestr))
866  call store_error(errmsg)
867  end if
868  !
869  ! -- Check to see if nodem is outside of active domain
870  if (this%nodem2(ignc) <= 0) then
871  call this%m2%dis%nodeu_to_string(nodeum, nodestr)
872  write (errmsg, *) &
873  trim(adjustl(this%m2%name))// &
874  ' Cell is outside active grid domain: '// &
875  trim(adjustl(nodestr))
876  call store_error(errmsg)
877  end if
878  !
879  ! -- Check to see if any nodejs are outside of active domain
880  do jidx = 1, this%numjs
881  if (this%nodesj(jidx, ignc) < 0) then
882  call this%m1%dis%nodeu_to_string(nodesuj(jidx), nodestr)
883  write (errmsg, *) &
884  trim(adjustl(this%m1%name))// &
885  ' Cell is outside active grid domain: '// &
886  trim(adjustl(nodestr))
887  call store_error(errmsg)
888  end if
889  end do
890  !
891  end do
892  !
893  ! -- Stop if errors
894  nerr = count_errors()
895  if (nerr > 0) then
896  call store_error('Errors encountered in GNC input file.')
897  call this%parser%StoreErrorUnit()
898  end if
899  !
900  write (this%iout, '(1x,a)') 'END OF GNCDATA'
901  else
902  write (errmsg, '(a)') 'Required GNCDATA block not found.'
903  call store_error(errmsg)
904  call this%parser%StoreErrorUnit()
905  end if
906  !
907  ! -- deallocate nodesuj array
908  deallocate (nodesuj)
909  !
910  ! -- Return
911  return
912  end subroutine read_data
913 
914  !> @brief Convert the user-based node number into a reduced number
915  !<
916  subroutine nodeu_to_noder(this, nodeu, noder, model)
917  ! -- modules
919  use simmodule, only: store_error
920  use simvariablesmodule, only: errmsg
921  ! -- dummy
922  class(ghostnodetype) :: this
923  integer(I4B), intent(in) :: nodeu
924  integer(I4B), intent(inout) :: noder
925  class(numericalmodeltype), intent(in) :: model
926  !
927  if (nodeu < 1 .or. nodeu > model%dis%nodesuser) then
928  write (errmsg, *) &
929  trim(adjustl(model%name))// &
930  ' node number < 0 or > model nodes: ', nodeu
931  call store_error(errmsg)
932  else
933  noder = model%dis%get_nodenumber(nodeu, 0)
934  end if
935  !
936  ! -- Return
937  return
938  end subroutine nodeu_to_noder
939 
940 end module ghostnodemodule
This module contains block parser methods.
Definition: BlockParser.f90:7
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:44
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
subroutine read_dimensions(this)
Single Model GNC Read Dimensions.
Definition: GhostNode.f90:736
subroutine nodeu_to_noder(this, nodeu, noder, model)
Convert the user-based node number into a reduced number.
Definition: GhostNode.f90:917
subroutine gnc_da(this)
Deallocate memory.
Definition: GhostNode.f90:639
subroutine gnc_fn(this, kiter, matrix_sln, condsat, ihc_opt, ivarcv_opt, ictm1_opt, ictm2_opt)
Fill GNC Newton terms.
Definition: GhostNode.f90:360
subroutine, public gnc_cr(gncobj, name_parent, inunit, iout)
Create new GNC exchange object.
Definition: GhostNode.f90:61
subroutine read_data(this)
Read a GNCDATA block.
Definition: GhostNode.f90:785
subroutine allocate_scalars(this)
Allocate gnc scalar variables.
Definition: GhostNode.f90:577
subroutine gnc_fc(this, kiter, matrix)
Fill matrix terms.
Definition: GhostNode.f90:296
subroutine gnc_df(this, m1, m2)
Initialize a gnc object.
Definition: GhostNode.f90:88
subroutine read_options(this)
Read a gnc options block.
Definition: GhostNode.f90:677
subroutine allocate_arrays(this)
Allocate gnc scalar variables.
Definition: GhostNode.f90:605
subroutine gnc_ot(this, ibudfl)
Single Model GNC Output.
Definition: GhostNode.f90:471
subroutine gnc_cq(this, flowja)
Add GNC to flowja.
Definition: GhostNode.f90:503
subroutine gnc_mc(this, matrix_sln)
Single or Two-Model GNC Map Connections.
Definition: GhostNode.f90:182
subroutine gnc_fmsav(this, kiter, matrix)
Store the n-m Picard conductance in cond prior to the Newton terms terms being added.
Definition: GhostNode.f90:264
subroutine gnc_ac(this, sparse)
Single or Two-Model GNC Add Connections.
Definition: GhostNode.f90:141
real(dp) function deltaqgnc(this, ignc)
Single Model deltaQgnc (ghost node correction flux)
Definition: GhostNode.f90:537
This module defines variable data types.
Definition: kind.f90:8
This module contains the base numerical package type.
This module contains simulation methods.
Definition: Sim.f90:10
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, public store_error_unit(iunit, terminate)
Store the file unit number.
Definition: Sim.f90:168
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=maxcharlen) errmsg
error message string
real(dp) function squadraticsaturationderivative(top, bot, x, eps)
@ brief Derivative of the quadratic saturation function