MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
NumericalModel.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: dp, i4b
6  use basedismodule, only: disbasetype
7  use sparsemodule, only: sparsematrix
9  use listmodule, only: listtype
13 
14  implicit none
15  private
18 
20  character(len=LINELENGTH), pointer :: filename => null() !input file name
21  integer(I4B), pointer :: neq => null() !number of equations
22  integer(I4B), pointer :: nja => null() !number of connections
23  integer(I4B), pointer :: moffset => null() !offset of this model in the solution
24  integer(I4B), pointer :: icnvg => null() !convergence flag
25  integer(I4B), dimension(:), pointer, contiguous :: ia => null() !csr row pointer
26  integer(I4B), dimension(:), pointer, contiguous :: ja => null() !csr columns
27  real(dp), dimension(:), pointer, contiguous :: x => null() !dependent variable (head, conc, etc)
28  real(dp), dimension(:), pointer, contiguous :: rhs => null() !right-hand side vector
29  real(dp), dimension(:), pointer, contiguous :: cond => null() !conductance matrix
30  integer(I4B), dimension(:), pointer, contiguous :: idxglo => null() !pointer to position in solution matrix
31  real(dp), dimension(:), pointer, contiguous :: xold => null() !dependent variable for previous timestep
32  real(dp), dimension(:), pointer, contiguous :: flowja => null() !intercell flows
33  integer(I4B), dimension(:), pointer, contiguous :: ibound => null() !ibound array
34  !
35  ! -- Derived types
36  type(listtype), pointer :: bndlist => null() !array of boundary packages for this model
37  class(disbasetype), pointer :: dis => null() !discretization object
38 
39  contains
40  !
41  ! -- Required for all models (override procedures defined in BaseModelType)
42  procedure :: model_df
43  procedure :: model_ar
44  procedure :: model_fp
45  procedure :: model_da
46  !
47  ! -- Methods specific to a numerical model
48  procedure :: model_ac
49  procedure :: model_mc
50  procedure :: model_rp
51  procedure :: model_ad
52  procedure :: model_reset
53  procedure :: model_solve
54  procedure :: model_cf
55  procedure :: model_fc
56  procedure :: model_ptcchk
57  procedure :: model_ptc
58  procedure :: model_nr
59  procedure :: model_cc
60  procedure :: model_nur
61  procedure :: model_cq
62  procedure :: model_bd
63  procedure :: model_bdcalc
64  procedure :: model_bdsave
65  procedure :: model_ot
66  procedure :: model_bdentry
67  !
68  ! -- Utility methods
69  procedure :: allocate_scalars
70  procedure :: allocate_arrays
71  procedure :: set_moffset
72  procedure :: set_idsoln
73  procedure :: set_xptr
74  procedure :: set_rhsptr
75  procedure :: set_iboundptr
76  procedure :: get_mrange
77  procedure :: get_mcellid
78  procedure :: get_mnodeu
79  procedure :: get_iasym
80  procedure :: create_lstfile
81  end type numericalmodeltype
82 
83 contains
84  !
85  ! -- Type-bound procedures for a numerical model
86  !
87  subroutine model_df(this)
88  class(numericalmodeltype) :: this
89  end subroutine model_df
90 
91  subroutine model_ac(this, sparse)
92  class(numericalmodeltype) :: this
93  type(sparsematrix), intent(inout) :: sparse
94  end subroutine model_ac
95 
96  subroutine model_mc(this, matrix_sln)
97  class(numericalmodeltype) :: this
98  class(matrixbasetype), pointer :: matrix_sln
99  end subroutine model_mc
100 
101  subroutine model_ar(this)
102  class(numericalmodeltype) :: this
103  end subroutine model_ar
104 
105  subroutine model_rp(this)
106  class(numericalmodeltype) :: this
107  end subroutine model_rp
108 
109  subroutine model_ad(this)
110  class(numericalmodeltype) :: this
111  end subroutine model_ad
112 
113  subroutine model_reset(this)
114  use bndmodule, only: bndtype, getbndfromlist
115  class(numericalmodeltype) :: this
116  ! local
117  class(bndtype), pointer :: packobj
118  integer(I4B) :: ip
119 
120  do ip = 1, this%bndlist%Count()
121  packobj => getbndfromlist(this%bndlist, ip)
122  call packobj%bnd_reset()
123  end do
124 
125  end subroutine model_reset
126 
127  subroutine model_solve(this)
128  class(numericalmodeltype) :: this
129  end subroutine model_solve
130 
131  subroutine model_cf(this, kiter)
132  class(numericalmodeltype) :: this
133  integer(I4B), intent(in) :: kiter
134  end subroutine model_cf
135 
136  subroutine model_fc(this, kiter, matrix_sln, inwtflag)
137  class(numericalmodeltype) :: this
138  integer(I4B), intent(in) :: kiter
139  class(matrixbasetype), pointer :: matrix_sln
140  integer(I4B), intent(in) :: inwtflag
141  end subroutine model_fc
142 
143  subroutine model_ptcchk(this, iptc)
144  class(numericalmodeltype) :: this
145  integer(I4B), intent(inout) :: iptc
146  iptc = 0
147  end subroutine model_ptcchk
148 
149  subroutine model_ptc(this, vec_residual, iptc, ptcf)
150  class(numericalmodeltype) :: this
151  class(vectorbasetype), pointer :: vec_residual
152  integer(I4B), intent(inout) :: iptc
153  real(DP), intent(inout) :: ptcf
154  end subroutine model_ptc
155 
156  subroutine model_nr(this, kiter, matrix, inwtflag)
157  class(numericalmodeltype) :: this
158  integer(I4B), intent(in) :: kiter
159  class(matrixbasetype), pointer :: matrix
160  integer(I4B), intent(in) :: inwtflag
161  end subroutine model_nr
162 
163  subroutine model_cc(this, innertot, kiter, iend, icnvgmod, cpak, ipak, dpak)
164  class(numericalmodeltype) :: this
165  integer(I4B), intent(in) :: innertot
166  integer(I4B), intent(in) :: kiter
167  integer(I4B), intent(in) :: iend
168  integer(I4B), intent(in) :: icnvgmod
169  character(len=LENPAKLOC), intent(inout) :: cpak
170  integer(I4B), intent(inout) :: ipak
171  real(DP), intent(inout) :: dpak
172  end subroutine model_cc
173 
174  subroutine model_nur(this, neqmod, x, xtemp, dx, inewtonur, dxmax, locmax)
175  class(numericalmodeltype) :: this
176  integer(I4B), intent(in) :: neqmod
177  real(DP), dimension(neqmod), intent(inout) :: x
178  real(DP), dimension(neqmod), intent(in) :: xtemp
179  real(DP), dimension(neqmod), intent(inout) :: dx
180  integer(I4B), intent(inout) :: inewtonur
181  real(DP), intent(inout) :: dxmax
182  integer(I4B), intent(inout) :: locmax
183  end subroutine model_nur
184 
185  subroutine model_cq(this, icnvg, isuppress_output)
186  class(numericalmodeltype) :: this
187  integer(I4B), intent(in) :: icnvg
188  integer(I4B), intent(in) :: isuppress_output
189  end subroutine model_cq
190 
191  subroutine model_bd(this, icnvg, isuppress_output)
192  class(numericalmodeltype) :: this
193  integer(I4B), intent(in) :: icnvg
194  integer(I4B), intent(in) :: isuppress_output
195  end subroutine model_bd
196 
197  subroutine model_bdcalc(this, icnvg)
198  class(numericalmodeltype) :: this
199  integer(I4B), intent(in) :: icnvg
200  end subroutine model_bdcalc
201 
202  subroutine model_bdsave(this, icnvg)
203  class(numericalmodeltype) :: this
204  integer(I4B), intent(in) :: icnvg
205  end subroutine model_bdsave
206 
207  subroutine model_ot(this)
208  class(numericalmodeltype) :: this
209  end subroutine model_ot
210 
211  subroutine model_bdentry(this, budterm, budtxt, rowlabel)
212  class(numericalmodeltype) :: this
213  real(DP), dimension(:, :), intent(in) :: budterm
214  character(len=LENBUDTXT), dimension(:), intent(in) :: budtxt
215  character(len=*), intent(in) :: rowlabel
216  end subroutine model_bdentry
217 
218  subroutine model_fp(this)
219  class(numericalmodeltype) :: this
220  end subroutine model_fp
221 
222  subroutine model_da(this)
223  ! -- modules
225  class(numericalmodeltype) :: this
226 
227  ! -- Scalars
228  call mem_deallocate(this%neq)
229  call mem_deallocate(this%nja)
230  call mem_deallocate(this%icnvg)
231  call mem_deallocate(this%moffset)
232  deallocate (this%filename)
233  !
234  ! -- Arrays
235  call mem_deallocate(this%xold)
236  call mem_deallocate(this%flowja, 'FLOWJA', this%memoryPath)
237  call mem_deallocate(this%idxglo)
238  !
239  ! -- derived types
240  call this%bndlist%Clear()
241  deallocate (this%bndlist)
242  !
243  ! -- nullify pointers
244  call mem_deallocate(this%x, 'X', this%memoryPath)
245  call mem_deallocate(this%rhs, 'RHS', this%memoryPath)
246  call mem_deallocate(this%ibound, 'IBOUND', this%memoryPath)
247  !
248  ! -- BaseModelType
249  call this%BaseModelType%model_da()
250  !
251  end subroutine model_da
252 
253  subroutine set_moffset(this, moffset)
254  class(numericalmodeltype) :: this
255  integer(I4B), intent(in) :: moffset
256  this%moffset = moffset
257  end subroutine set_moffset
258 
259  subroutine get_mrange(this, mstart, mend)
260  class(numericalmodeltype) :: this
261  integer(I4B), intent(inout) :: mstart
262  integer(I4B), intent(inout) :: mend
263  mstart = this%moffset + 1
264  mend = mstart + this%neq - 1
265  end subroutine get_mrange
266 
267  subroutine set_idsoln(this, id)
268  class(numericalmodeltype) :: this
269  integer(I4B), intent(in) :: id
270  this%idsoln = id
271  end subroutine set_idsoln
272 
273  subroutine allocate_scalars(this, modelname)
275  class(numericalmodeltype) :: this
276  character(len=*), intent(in) :: modelname
277  !
278  ! -- allocate basetype members
279  call this%BaseModelType%allocate_scalars(modelname)
280  !
281  ! -- allocate members from this type
282  call mem_allocate(this%neq, 'NEQ', this%memoryPath)
283  call mem_allocate(this%nja, 'NJA', this%memoryPath)
284  call mem_allocate(this%icnvg, 'ICNVG', this%memoryPath)
285  call mem_allocate(this%moffset, 'MOFFSET', this%memoryPath)
286  allocate (this%filename)
287  allocate (this%bndlist)
288  !
289  this%filename = ''
290  this%neq = 0
291  this%nja = 0
292  this%icnvg = 0
293  this%moffset = 0
294  end subroutine allocate_scalars
295 
296  subroutine allocate_arrays(this)
297  use constantsmodule, only: dzero
299  class(numericalmodeltype) :: this
300  integer(I4B) :: i
301  !
302  call mem_allocate(this%xold, this%neq, 'XOLD', this%memoryPath)
303  call mem_allocate(this%flowja, this%nja, 'FLOWJA', this%memoryPath)
304  call mem_allocate(this%idxglo, this%nja, 'IDXGLO', this%memoryPath)
305  !
306  ! -- initialize
307  do i = 1, size(this%flowja)
308  this%flowja(i) = dzero
309  end do
310  end subroutine allocate_arrays
311 
312  subroutine set_xptr(this, xsln, sln_offset, varNameTgt, memPathTgt)
314  ! -- dummy
315  class(numericalmodeltype) :: this
316  real(DP), dimension(:), pointer, contiguous, intent(in) :: xsln
317  integer(I4B) :: sln_offset
318  character(len=*), intent(in) :: varNameTgt
319  character(len=*), intent(in) :: memPathTgt
320  ! -- local
321  integer(I4B) :: offset
322  ! -- code
323  offset = this%moffset - sln_offset
324  this%x => xsln(offset + 1:offset + this%neq)
325  call mem_checkin(this%x, 'X', this%memoryPath, varnametgt, mempathtgt)
326  end subroutine set_xptr
327 
328  subroutine set_rhsptr(this, rhssln, sln_offset, varNameTgt, memPathTgt)
330  ! -- dummy
331  class(numericalmodeltype) :: this
332  real(DP), dimension(:), pointer, contiguous, intent(in) :: rhssln
333  integer(I4B) :: sln_offset
334  character(len=*), intent(in) :: varNameTgt
335  character(len=*), intent(in) :: memPathTgt
336  ! -- local
337  integer(I4B) :: offset
338  ! -- code
339  offset = this%moffset - sln_offset
340  this%rhs => rhssln(offset + 1:offset + this%neq)
341  call mem_checkin(this%rhs, 'RHS', this%memoryPath, varnametgt, mempathtgt)
342  end subroutine set_rhsptr
343 
344  subroutine set_iboundptr(this, iboundsln, sln_offset, varNameTgt, memPathTgt)
346  ! -- dummy
347  class(numericalmodeltype) :: this
348  integer(I4B), dimension(:), pointer, contiguous, intent(in) :: iboundsln
349  integer(I4B) :: sln_offset
350  character(len=*), intent(in) :: varNameTgt
351  character(len=*), intent(in) :: memPathTgt
352  ! -- local
353  integer(I4B) :: offset
354  ! -- code
355  offset = this%moffset - sln_offset
356  this%ibound => iboundsln(offset + 1:offset + this%neq)
357  call mem_checkin(this%ibound, 'IBOUND', this%memoryPath, varnametgt, &
358  mempathtgt)
359  end subroutine set_iboundptr
360 
361  subroutine get_mcellid(this, node, mcellid)
362  use bndmodule, only: bndtype, getbndfromlist
363  class(numericalmodeltype) :: this
364  integer(I4B), intent(in) :: node
365  character(len=*), intent(inout) :: mcellid
366  ! -- local
367  character(len=20) :: cellid
368  integer(I4B) :: ip, ipaknode, istart, istop
369  class(bndtype), pointer :: packobj
370 
371  if (node < 1) then
372  cellid = ''
373  else if (node <= this%dis%nodes) then
374  call this%dis%noder_to_string(node, cellid)
375  else
376  cellid = '***ERROR***'
377  ipaknode = node - this%dis%nodes
378  istart = 1
379  do ip = 1, this%bndlist%Count()
380  packobj => getbndfromlist(this%bndlist, ip)
381  if (packobj%npakeq == 0) cycle
382  istop = istart + packobj%npakeq - 1
383  if (istart <= ipaknode .and. ipaknode <= istop) then
384  write (cellid, '(a, a, a, i0, a, i0, a)') '(', &
385  trim(packobj%filtyp), '_', &
386  packobj%ibcnum, '-', ipaknode - packobj%ioffset, ')'
387  exit
388  end if
389  istart = istop + 1
390  end do
391  end if
392  write (mcellid, '(i0, a, a, a, a)') this%id, '_', this%macronym, '-', &
393  trim(adjustl(cellid))
394  end subroutine get_mcellid
395 
396  subroutine get_mnodeu(this, node, nodeu)
397  use bndmodule, only: bndtype, getbndfromlist
398  class(numericalmodeltype) :: this
399  integer(I4B), intent(in) :: node
400  integer(I4B), intent(inout) :: nodeu
401  ! -- local
402  if (node <= this%dis%nodes) then
403  nodeu = this%dis%get_nodeuser(node)
404  else
405  nodeu = -(node - this%dis%nodes)
406  end if
407  end subroutine get_mnodeu
408 
409  function get_iasym(this) result(iasym)
410  class(numericalmodeltype) :: this
411  integer(I4B) :: iasym
412  iasym = 0
413  end function get_iasym
414 
415  function castasnumericalmodelclass(obj) result(res)
416  implicit none
417  class(*), pointer, intent(inout) :: obj
418  class(numericalmodeltype), pointer :: res
419  !
420  res => null()
421  if (.not. associated(obj)) return
422  !
423  select type (obj)
424  class is (numericalmodeltype)
425  res => obj
426  end select
427  end function castasnumericalmodelclass
428 
429  subroutine addnumericalmodeltolist(list, model)
430  implicit none
431  ! -- dummy
432  type(listtype), intent(inout) :: list
433  class(numericalmodeltype), pointer, intent(inout) :: model
434  ! -- local
435  class(*), pointer :: obj
436  !
437  obj => model
438  call list%Add(obj)
439  end subroutine addnumericalmodeltolist
440 
441  function getnumericalmodelfromlist(list, idx) result(res)
442  implicit none
443  ! -- dummy
444  type(listtype), intent(inout) :: list
445  integer(I4B), intent(in) :: idx
446  class(numericalmodeltype), pointer :: res
447  ! -- local
448  class(*), pointer :: obj
449  !
450  obj => list%GetItem(idx)
451  res => castasnumericalmodelclass(obj)
452  end function getnumericalmodelfromlist
453 
454  subroutine create_lstfile(this, lst_fname, model_fname, defined, headertxt)
455  ! -- modules
456  use kindmodule, only: lgp
458  ! -- dummy
459  class(numericalmodeltype) :: this
460  character(len=*), intent(inout) :: lst_fname
461  character(len=*), intent(in) :: model_fname
462  logical(LGP), intent(in) :: defined
463  character(len=*), intent(in) :: headertxt
464  ! -- local
465  integer(I4B) :: i, istart, istop
466  !
467  ! -- set list file name if not provided
468  if (.not. defined) then
469  !
470  ! -- initialize
471  lst_fname = ' '
472  istart = 0
473  istop = len_trim(model_fname)
474  !
475  ! -- identify '.' character position from back of string
476  do i = istop, 1, -1
477  if (model_fname(i:i) == '.') then
478  istart = i
479  exit
480  end if
481  end do
482  !
483  ! -- if not found start from string end
484  if (istart == 0) istart = istop + 1
485  !
486  ! -- set list file name
487  lst_fname = model_fname(1:istart)
488  istop = istart + 3
489  lst_fname(istart:istop) = '.lst'
490  end if
491  !
492  ! -- create the list file
493  this%iout = getunit()
494  call openfile(this%iout, 0, lst_fname, 'LIST', filstat_opt='REPLACE')
495  !
496  ! -- write list file header
497  call write_listfile_header(this%iout, headertxt)
498  end subroutine create_lstfile
499 
500 end module numericalmodelmodule
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 linelength
maximum length of a standard line
Definition: Constants.f90:45
integer(i4b), parameter lenpakloc
maximum length of a package location
Definition: Constants.f90:50
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
integer(i4b), parameter lenbudtxt
maximum length of a budget component names
Definition: Constants.f90:37
integer(i4b) function, public getunit()
Get a free unit number.
subroutine, public openfile(iu, iout, fname, ftype, fmtarg_opt, accarg_opt, filstat_opt, mode_opt)
Open a file.
Definition: InputOutput.f90:30
This module defines variable data types.
Definition: kind.f90:8
subroutine, public addnumericalmodeltolist(list, model)
subroutine model_ac(this, sparse)
subroutine model_df(this)
subroutine model_cc(this, innertot, kiter, iend, icnvgmod, cpak, ipak, dpak)
subroutine model_bdsave(this, icnvg)
subroutine model_ad(this)
subroutine model_solve(this)
subroutine get_mrange(this, mstart, mend)
subroutine model_mc(this, matrix_sln)
class(numericalmodeltype) function, pointer, public getnumericalmodelfromlist(list, idx)
subroutine set_rhsptr(this, rhssln, sln_offset, varNameTgt, memPathTgt)
subroutine model_ar(this)
subroutine model_nr(this, kiter, matrix, inwtflag)
subroutine set_idsoln(this, id)
subroutine get_mnodeu(this, node, nodeu)
subroutine model_nur(this, neqmod, x, xtemp, dx, inewtonur, dxmax, locmax)
subroutine model_cq(this, icnvg, isuppress_output)
subroutine model_da(this)
subroutine model_bd(this, icnvg, isuppress_output)
subroutine model_bdcalc(this, icnvg)
subroutine model_fp(this)
subroutine get_mcellid(this, node, mcellid)
subroutine model_reset(this)
class(numericalmodeltype) function, pointer castasnumericalmodelclass(obj)
subroutine model_ptc(this, vec_residual, iptc, ptcf)
subroutine set_moffset(this, moffset)
subroutine model_ptcchk(this, iptc)
subroutine create_lstfile(this, lst_fname, model_fname, defined, headertxt)
subroutine allocate_arrays(this)
integer(i4b) function get_iasym(this)
subroutine allocate_scalars(this, modelname)
subroutine set_iboundptr(this, iboundsln, sln_offset, varNameTgt, memPathTgt)
subroutine model_fc(this, kiter, matrix_sln, inwtflag)
subroutine model_ot(this)
subroutine model_cf(this, kiter)
subroutine set_xptr(this, xsln, sln_offset, varNameTgt, memPathTgt)
subroutine model_bdentry(this, budterm, budtxt, rowlabel)
subroutine model_rp(this)
This module contains version information.
Definition: version.f90:7
subroutine write_listfile_header(iout, cmodel_type, write_sys_command, write_kind_info)
@ brief Write program header
Definition: version.f90:98
Highest level model type. All models extend this parent type.
Definition: BaseModel.f90:13
@ brief BndType
A generic heterogeneous doubly-linked list.
Definition: List.f90:14