MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
swf-dfw.f90
Go to the documentation of this file.
1 !> @brief Stream Network Flow (SWF) Diffusive Wave (DFW) Module
2 !!
3 !! This module solves one-dimensional flow routing using a diffusive
4 !< wave approach.
5 module swfdfwmodule
6 
7  use kindmodule, only: dp, i4b, lgp
9  dzero, dhalf, done, dtwo, dthree, &
11  dprec, dem10
19  use basedismodule, only: disbasetype
20  use swfcxsmodule, only: swfcxstype
21  use obsmodule, only: obstype, obs_cr
23  use observemodule, only: observetype
25 
26  implicit none
27  private
28  public :: swfdfwtype, dfw_cr
29 
30  type, extends(numericalpackagetype) :: swfdfwtype
31 
32  ! user-provided input
33  integer(I4B), pointer :: is2d => null() !< flag to indicate this model is 2D overland flow and not 1d channel flow
34  integer(I4B), pointer :: icentral => null() !< flag to use central in space weighting (default is upstream weighting)
35  integer(I4B), pointer :: iswrcond => null() !< flag to activate the dev SWR conductance formulation
36  real(DP), pointer :: unitconv !< conversion factor used in mannings equation; calculated from timeconv and lengthconv
37  real(DP), pointer :: timeconv !< conversion factor from model length units to meters (1.0 if model uses meters for length)
38  real(DP), pointer :: lengthconv !< conversion factor from model time units to seconds (1.0 if model uses seconds for time)
39  real(DP), dimension(:), pointer, contiguous :: hnew => null() !< pointer to model xnew
40  real(DP), dimension(:), pointer, contiguous :: manningsn => null() !< mannings roughness for each reach
41  integer(I4B), dimension(:), pointer, contiguous :: idcxs => null() !< cross section id for each reach
42  integer(I4B), dimension(:), pointer, contiguous :: ibound => null() !< pointer to model ibound
43  integer(I4B), dimension(:), pointer, contiguous :: icelltype => null() !< set to 1 and is accessed by chd for checking
44 
45  ! velocity
46  integer(I4B), pointer :: icalcvelocity => null() !< flag to indicate velocity will be calculated (always on)
47  integer(I4B), pointer :: isavvelocity => null() !< flag to indicate velocity will be saved
48  real(DP), dimension(:, :), pointer, contiguous :: vcomp => null() !< velocity components: vx, vy, vz (nodes, 3)
49  real(DP), dimension(:), pointer, contiguous :: vmag => null() !< velocity magnitude (of size nodes)
50  integer(I4B), pointer :: nedges => null() !< number of cell edges
51  integer(I4B), pointer :: lastedge => null() !< last edge number
52  integer(I4B), dimension(:), pointer, contiguous :: nodedge => null() !< array of node numbers that have edges
53  integer(I4B), dimension(:), pointer, contiguous :: ihcedge => null() !< edge type (horizontal or vertical)
54  real(DP), dimension(:, :), pointer, contiguous :: propsedge => null() !< edge properties (Q, area, nx, ny, distance)
55  real(DP), dimension(:), pointer, contiguous :: grad_dhds_mag => null() !< magnitude of the gradient (of size nodes)
56  real(DP), dimension(:), pointer, contiguous :: dhdsja => null() !< gradient for each connection (of size njas)
57 
58  ! observation data
59  integer(I4B), pointer :: inobspkg => null() !< unit number for obs package
60  type(ObsType), pointer :: obs => null() !< observation package
61 
62  ! pointer to cross section data
63  type(SwfCxsType), pointer :: cxs
64 
65  contains
66 
67  procedure :: dfw_df
68  procedure :: allocate_scalars
69  procedure :: allocate_arrays
70  procedure :: dfw_load
71  procedure :: source_options
72  procedure :: log_options
73  procedure :: source_griddata
74  procedure :: log_griddata
75  procedure :: dfw_ar
76  procedure :: dfw_rp
77  procedure :: dfw_ad
78  procedure :: dfw_fc
79  procedure :: dfw_qnm_fc_nr
80  !procedure :: dfw_qnm_fc
81  procedure :: dfw_fn
82  procedure :: dfw_nur
83  procedure :: dfw_cq
84  procedure :: dfw_bd
85  procedure :: dfw_save_model_flows
86  procedure :: dfw_print_model_flows
87  procedure :: dfw_da
88  procedure :: dfw_df_obs
89  procedure :: dfw_rp_obs
90  procedure :: dfw_bd_obs
91  procedure :: qcalc
92  procedure :: get_cond
93  procedure :: get_cond_swr
94  procedure :: get_cond_n
95  procedure :: get_flow_area_nm
96  procedure :: calc_velocity
97  procedure :: sav_velocity
98  procedure, public :: increase_edge_count
99  procedure, public :: set_edge_properties
100  procedure :: calc_dhds
101  procedure :: write_cxs_tables
102 
103  end type swfdfwtype
104 
105 contains
106 
107  !> @brief create package
108  !<
109  subroutine dfw_cr(dfwobj, name_model, input_mempath, inunit, iout, &
110  cxs)
111  ! modules
113  ! dummy
114  type(SwfDfwType), pointer :: dfwobj !< object to create
115  character(len=*), intent(in) :: name_model !< name of the SWF model
116  character(len=*), intent(in) :: input_mempath !< memory path
117  integer(I4B), intent(in) :: inunit !< flag to indicate if package is active
118  integer(I4B), intent(in) :: iout !< unit number for output
119  type(SwfCxsType), pointer, intent(in) :: cxs !< the pointer to the cxs package
120  ! locals
121  logical(LGP) :: found_fname
122  ! formats
123  character(len=*), parameter :: fmtheader = &
124  "(1x, /1x, 'DFW -- DIFFUSIVE WAVE (DFW) PACKAGE, VERSION 1, 9/25/2023', &
125  &' INPUT READ FROM MEMPATH: ', A, /)"
126  !
127  ! Create the object
128  allocate (dfwobj)
129 
130  ! create name and memory path
131  call dfwobj%set_names(1, name_model, 'DFW', 'DFW')
132 
133  ! Allocate scalars
134  call dfwobj%allocate_scalars()
135 
136  ! Set variables
137  dfwobj%input_mempath = input_mempath
138  dfwobj%inunit = inunit
139  dfwobj%iout = iout
140 
141  ! set name of input file
142  call mem_set_value(dfwobj%input_fname, 'INPUT_FNAME', dfwobj%input_mempath, &
143  found_fname)
144 
145  ! Set a pointers to passed in objects
146  dfwobj%cxs => cxs
147 
148  ! create obs package
149  call obs_cr(dfwobj%obs, dfwobj%inobspkg)
150 
151  ! check if dfw is enabled
152  if (inunit > 0) then
153 
154  ! Print a message identifying the package.
155  write (iout, fmtheader) input_mempath
156 
157  end if
158 
159  end subroutine dfw_cr
160 
161  !> @brief load data from IDM to package
162  !<
163  subroutine dfw_df(this, dis)
164  ! dummy
165  class(SwfDfwType) :: this !< this instance
166  class(DisBaseType), pointer, intent(inout) :: dis !< the pointer to the discretization
167  ! locals
168  character(len=10) :: distype = ''
169 
170  ! Set a pointers to passed in objects
171  this%dis => dis
172 
173  ! Set the distype (either DISV1D or DIS2D)
174  call this%dis%get_dis_type(distype)
175  if (distype == "DIS2D") then
176  this%is2d = 1
177  end if
178  if (distype == "DISV2D") then
179  this%is2d = 1
180  end if
181 
182  ! check if dfw is enabled
183  ! this will need to become if (.not. present(dfw_options)) then
184  !if (inunit > 0) then
185 
186  ! allocate arrays
187  call this%allocate_arrays()
188 
189  ! load dfw
190  call this%dfw_load()
191 
192  !end if
193 
194  end subroutine dfw_df
195 
196  !> @ brief Allocate scalars
197  !!
198  !! Allocate and initialize scalars for the package. The base model
199  !! allocate scalars method is also called.
200  !!
201  !<
202  subroutine allocate_scalars(this)
203  ! modules
204  ! dummy
205  class(SwfDfwtype) :: this !< this instance
206  !
207  ! allocate scalars in NumericalPackageType
208  call this%NumericalPackageType%allocate_scalars()
209  !
210  ! Allocate scalars
211  call mem_allocate(this%is2d, 'IS2D', this%memoryPath)
212  call mem_allocate(this%icentral, 'ICENTRAL', this%memoryPath)
213  call mem_allocate(this%iswrcond, 'ISWRCOND', this%memoryPath)
214  call mem_allocate(this%unitconv, 'UNITCONV', this%memoryPath)
215  call mem_allocate(this%lengthconv, 'LENGTHCONV', this%memoryPath)
216  call mem_allocate(this%timeconv, 'TIMECONV', this%memoryPath)
217  call mem_allocate(this%inobspkg, 'INOBSPKG', this%memoryPath)
218  call mem_allocate(this%icalcvelocity, 'ICALCVELOCITY', this%memoryPath)
219  call mem_allocate(this%isavvelocity, 'ISAVVELOCITY', this%memoryPath)
220  call mem_allocate(this%nedges, 'NEDGES', this%memoryPath)
221  call mem_allocate(this%lastedge, 'LASTEDGE', this%memoryPath)
222 
223  this%is2d = 0
224  this%icentral = 0
225  this%iswrcond = 0
226  this%unitconv = done
227  this%lengthconv = done
228  this%timeconv = done
229  this%inobspkg = 0
230  this%icalcvelocity = 0
231  this%isavvelocity = 0
232  this%nedges = 0
233  this%lastedge = 0
234 
235  end subroutine allocate_scalars
236 
237  !> @brief allocate memory for arrays
238  !<
239  subroutine allocate_arrays(this)
240  ! dummy
241  class(SwfDfwType) :: this !< this instance
242  ! locals
243  integer(I4B) :: n
244  !
245  ! user-provided input
246  call mem_allocate(this%manningsn, this%dis%nodes, &
247  'MANNINGSN', this%memoryPath)
248  call mem_allocate(this%idcxs, this%dis%nodes, &
249  'IDCXS', this%memoryPath)
250  call mem_allocate(this%icelltype, this%dis%nodes, &
251  'ICELLTYPE', this%memoryPath)
252 
253  ! optional arrays
254  call mem_allocate(this%nodedge, 0, 'NODEDGE', this%memoryPath)
255  call mem_allocate(this%ihcedge, 0, 'IHCEDGE', this%memoryPath)
256  call mem_allocate(this%propsedge, 0, 0, 'PROPSEDGE', this%memoryPath)
257 
258  ! Specific discharge is (re-)allocated when nedges is known
259  call mem_allocate(this%vcomp, 3, 0, 'VCOMP', this%memoryPath)
260  call mem_allocate(this%vmag, 0, 'VMAG', this%memoryPath)
261 
262  do n = 1, this%dis%nodes
263  this%manningsn(n) = dzero
264  this%idcxs(n) = 0
265  this%icelltype(n) = 1
266  end do
267 
268  ! for 2d models, need to calculate and store dhds magnitude
269  if (this%is2d == 1) then
270  call mem_allocate(this%grad_dhds_mag, this%dis%nodes, &
271  'GRAD_DHDS_MAG', this%memoryPath)
272  call mem_allocate(this%dhdsja, this%dis%njas, &
273  'DHDSJA', this%memoryPath)
274  do n = 1, this%dis%nodes
275  this%grad_dhds_mag(n) = dzero
276  end do
277  do n = 1, this%dis%njas
278  this%dhdsja(n) = dzero
279  end do
280  end if
281 
282  end subroutine allocate_arrays
283 
284  !> @brief load data from IDM to package
285  !<
286  subroutine dfw_load(this)
287  ! dummy
288  class(SwfDfwType) :: this !< this instance
289 
290  ! source input data
291  call this%source_options()
292  call this%source_griddata()
293 
294  end subroutine dfw_load
295 
296  !> @brief Copy options from IDM into package
297  !<
298  subroutine source_options(this)
299  ! modules
300  use kindmodule, only: lgp
305  ! dummy
306  class(SwfDfwType) :: this !< this instance
307  ! locals
308  integer(I4B) :: isize
309  type(SwfDfwParamFoundType) :: found
310  type(CharacterStringType), dimension(:), pointer, &
311  contiguous :: obs6_fnames
312 
313  ! update defaults with idm sourced values
314  call mem_set_value(this%icentral, 'ICENTRAL', &
315  this%input_mempath, found%icentral)
316  call mem_set_value(this%iswrcond, 'ISWRCOND', &
317  this%input_mempath, found%iswrcond)
318  call mem_set_value(this%lengthconv, 'LENGTHCONV', &
319  this%input_mempath, found%lengthconv)
320  call mem_set_value(this%timeconv, 'TIMECONV', &
321  this%input_mempath, found%timeconv)
322  call mem_set_value(this%iprflow, 'IPRFLOW', &
323  this%input_mempath, found%iprflow)
324  call mem_set_value(this%ipakcb, 'IPAKCB', &
325  this%input_mempath, found%ipakcb)
326  call mem_set_value(this%isavvelocity, 'ISAVVELOCITY', &
327  this%input_mempath, found%isavvelocity)
328 
329  ! save flows option active
330  if (found%icentral) this%icentral = 1
331  if (found%ipakcb) this%ipakcb = -1
332 
333  ! calculate unit conversion
334  this%unitconv = this%lengthconv**donethird
335  this%unitconv = this%unitconv * this%timeconv
336 
337  ! save velocity active
338  if (found%isavvelocity) this%icalcvelocity = this%isavvelocity
339 
340  ! check for obs6_filename
341  call get_isize('OBS6_FILENAME', this%input_mempath, isize)
342  if (isize > 0) then
343  !
344  if (isize /= 1) then
345  errmsg = 'Multiple OBS6 keywords detected in OPTIONS block.'// &
346  ' Only one OBS6 entry allowed.'
347  call store_error(errmsg)
348  call store_error_filename(this%input_fname)
349  end if
350 
351  call mem_setptr(obs6_fnames, 'OBS6_FILENAME', this%input_mempath)
352 
353  found%obs6_filename = .true.
354  this%obs%inputFilename = obs6_fnames(1)
355  this%obs%active = .true.
356  this%inobspkg = getunit()
357  this%obs%inUnitObs = this%inobspkg
358  call openfile(this%inobspkg, this%iout, this%obs%inputFilename, 'OBS')
359  call this%obs%obs_df(this%iout, this%packName, this%filtyp, this%dis)
360  call this%dfw_df_obs()
361  end if
362 
363  ! log values to list file
364  if (this%iout > 0) then
365  call this%log_options(found)
366  end if
367 
368  end subroutine source_options
369 
370  !> @brief Write user options to list file
371  !<
372  subroutine log_options(this, found)
374  class(SwfDfwType) :: this !< this instance
375  type(SwfDfwParamFoundType), intent(in) :: found
376 
377  write (this%iout, '(1x,a)') 'Setting DFW Options'
378 
379  if (found%lengthconv) then
380  write (this%iout, '(4x,a, G0)') 'Mannings length conversion value &
381  &specified as ', this%lengthconv
382  end if
383 
384  if (found%timeconv) then
385  write (this%iout, '(4x,a, G0)') 'Mannings time conversion value &
386  &specified as ', this%timeconv
387  end if
388 
389  if (found%lengthconv .or. found%timeconv) then
390  write (this%iout, '(4x,a, G0)') 'Mannings conversion value calculated &
391  &from user-provided length_conversion and &
392  &time_conversion is ', this%unitconv
393  end if
394 
395  if (found%iprflow) then
396  write (this%iout, '(4x,a)') 'Cell-by-cell flow information will be printed &
397  &to listing file whenever ICBCFL is not zero.'
398  end if
399 
400  if (found%ipakcb) then
401  write (this%iout, '(4x,a)') 'Cell-by-cell flow information will be printed &
402  &to listing file whenever ICBCFL is not zero.'
403  end if
404 
405  if (found%obs6_filename) then
406  write (this%iout, '(4x,a)') 'Observation package is active.'
407  end if
408 
409  if (found%isavvelocity) &
410  write (this%iout, '(4x,a)') 'Velocity will be calculated at cell &
411  &centers and written to DATA-VCOMP in budget &
412  &file when requested.'
413 
414  if (found%iswrcond) then
415  write (this%iout, '(4x,a, G0)') 'Conductance will be calculated using &
416  &the SWR development option.'
417  end if
418 
419  write (this%iout, '(1x,a,/)') 'End Setting DFW Options'
420 
421  end subroutine log_options
422 
423  !> @brief copy griddata from IDM to package
424  !<
425  subroutine source_griddata(this)
426  ! modules
432  ! dummy
433  class(SwfDfwType) :: this !< this instance
434  ! locals
435  character(len=LENMEMPATH) :: idmMemoryPath
436  type(SwfDfwParamFoundType) :: found
437  integer(I4B), dimension(:), pointer, contiguous :: map
438 
439  ! set memory path
440  idmmemorypath = create_mem_path(this%name_model, 'DFW', idm_context)
441 
442  ! set map to convert user input data into reduced data
443  map => null()
444  if (this%dis%nodes < this%dis%nodesuser) map => this%dis%nodeuser
445 
446  ! update defaults with idm sourced values
447  call mem_set_value(this%manningsn, 'MANNINGSN', &
448  idmmemorypath, map, found%manningsn)
449  call mem_set_value(this%idcxs, 'IDCXS', idmmemorypath, map, found%idcxs)
450 
451  ! ensure MANNINGSN was found
452  if (.not. found%manningsn) then
453  write (errmsg, '(a)') 'Error in GRIDDATA block: MANNINGSN not found.'
454  call store_error(errmsg)
455  end if
456 
457  if (count_errors() > 0) then
458  call store_error_filename(this%input_fname)
459  end if
460 
461  ! log griddata
462  if (this%iout > 0) then
463  call this%log_griddata(found)
464  end if
465 
466  end subroutine source_griddata
467 
468  !> @brief log griddata to list file
469  !<
470  subroutine log_griddata(this, found)
472  class(SwfDfwType) :: this !< this instance
473  type(SwfDfwParamFoundType), intent(in) :: found
474 
475  write (this%iout, '(1x,a)') 'Setting DFW Griddata'
476 
477  if (found%manningsn) then
478  write (this%iout, '(4x,a)') 'MANNINGSN set from input file'
479  end if
480 
481  if (found%idcxs) then
482  write (this%iout, '(4x,a)') 'IDCXS set from input file'
483  end if
484 
485  call this%write_cxs_tables()
486 
487  write (this%iout, '(1x,a,/)') 'End Setting DFW Griddata'
488 
489  end subroutine log_griddata
490 
491  subroutine write_cxs_tables(this)
492  ! modules
493  ! dummy
494  class(SwfDfwType) :: this !< this instance
495  ! local
496  ! integer(I4B) :: idcxs
497  ! integer(I4B) :: n
498 
499  !-- TODO: write cross section tables
500  ! do n = 1, this%dis%nodes
501  ! idcxs = this%idcxs(n)
502  ! if (idcxs > 0) then
503  ! call this%cxs%write_cxs_table(idcxs, this%width(n), this%slope(n), &
504  ! this%manningsn(n), this%unitconv)
505  ! end if
506  ! end do
507  end subroutine write_cxs_tables
508 
509  !> @brief allocate memory
510  !<
511  subroutine dfw_ar(this, ibound, hnew)
512  ! modules
513  ! dummy
514  class(SwfDfwType) :: this !< this instance
515  integer(I4B), dimension(:), pointer, contiguous :: ibound !< model ibound array
516  real(DP), dimension(:), pointer, contiguous, intent(inout) :: hnew !< pointer to model head array
517  ! local
518  integer(I4B) :: n
519 
520  ! store pointer to ibound
521  this%ibound => ibound
522  this%hnew => hnew
523 
524  if (this%icalcvelocity == 1) then
525  call mem_reallocate(this%vcomp, 3, this%dis%nodes, 'VCOMP', this%memoryPath)
526  call mem_reallocate(this%vmag, this%dis%nodes, 'VMAG', this%memoryPath)
527  call mem_reallocate(this%nodedge, this%nedges, 'NODEDGE', this%memoryPath)
528  call mem_reallocate(this%ihcedge, this%nedges, 'IHCEDGE', this%memoryPath)
529  call mem_reallocate(this%propsedge, 5, this%nedges, 'PROPSEDGE', &
530  this%memoryPath)
531  do n = 1, this%dis%nodes
532  this%vcomp(:, n) = dzero
533  this%vmag(n) = dzero
534  end do
535  end if
536 
537  ! observation data
538  call this%obs%obs_ar()
539 
540  end subroutine dfw_ar
541 
542  !> @brief allocate memory
543  !<
544  subroutine dfw_rp(this)
545  ! modules
546  ! dummy
547  class(SwfDfwType) :: this !< this instance
548 
549  ! read observations
550  call this%dfw_rp_obs()
551 
552  end subroutine dfw_rp
553 
554  !> @brief advance
555  !<
556  subroutine dfw_ad(this, irestore)
557  class(SwfDfwType) :: this !< this instance
558  integer(I4B), intent(in) :: irestore !< ATS flag for retrying time step (1) or advancing (0)
559 
560  ! Push simulated values to preceding time/subtime step
561  call this%obs%obs_ad()
562 
563  end subroutine dfw_ad
564 
565  !> @brief fill coefficients
566  !!
567  !! The DFW Package is entirely Newton based. All matrix and rhs terms
568  !! are added from thish routine.
569  !!
570  !<
571  subroutine dfw_fc(this, kiter, matrix_sln, idxglo, rhs, stage, stage_old)
572  ! modules
573  ! dummy
574  class(SwfDfwType) :: this !< this instance
575  integer(I4B) :: kiter
576  class(MatrixBaseType), pointer :: matrix_sln
577  integer(I4B), intent(in), dimension(:) :: idxglo
578  real(DP), intent(inout), dimension(:) :: rhs
579  real(DP), intent(inout), dimension(:) :: stage
580  real(DP), intent(inout), dimension(:) :: stage_old
581  ! local
582 
583  ! calculate dhds at cell center for 2d case
584  if (this%is2d == 1) then
585  call this%calc_dhds()
586  end if
587 
588  ! add qnm contributions to matrix equations
589  call this%dfw_qnm_fc_nr(kiter, matrix_sln, idxglo, rhs, stage, stage_old)
590 
591  end subroutine dfw_fc
592 
593  !> @brief fill coefficients
594  !!
595  !< Add qnm contributions to matrix equations
596  subroutine dfw_qnm_fc_nr(this, kiter, matrix_sln, idxglo, rhs, stage, stage_old)
597  ! modules
599  ! dummy
600  class(SwfDfwType) :: this !< this instance
601  integer(I4B) :: kiter
602  class(MatrixBaseType), pointer :: matrix_sln
603  integer(I4B), intent(in), dimension(:) :: idxglo
604  real(DP), intent(inout), dimension(:) :: rhs
605  real(DP), intent(inout), dimension(:) :: stage
606  real(DP), intent(inout), dimension(:) :: stage_old
607  ! local
608  integer(I4B) :: n, m, ii, idiag
609  real(DP) :: qnm
610  real(DP) :: qeps
611  real(DP) :: eps
612  real(DP) :: derv
613 
614  ! Calculate conductance and put into amat
615  do n = 1, this%dis%nodes
616 
617  ! Find diagonal position for row n
618  idiag = this%dis%con%ia(n)
619 
620  ! Loop through connections adding matrix terms
621  do ii = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
622 
623  ! skip for masked cells
624  if (this%dis%con%mask(ii) == 0) cycle
625 
626  ! connection variables
627  m = this%dis%con%ja(ii)
628 
629  ! Fill the qnm term on the right-hand side
630  qnm = this%qcalc(n, m, stage(n), stage(m), ii)
631  rhs(n) = rhs(n) - qnm
632 
633  ! Derivative calculation and fill of n terms
634  eps = get_perturbation(stage(n))
635  qeps = this%qcalc(n, m, stage(n) + eps, stage(m), ii)
636  derv = (qeps - qnm) / eps
637  call matrix_sln%add_value_pos(idxglo(idiag), derv)
638  rhs(n) = rhs(n) + derv * stage(n)
639 
640  ! Derivative calculation and fill of m terms
641  eps = get_perturbation(stage(m))
642  qeps = this%qcalc(n, m, stage(n), stage(m) + eps, ii)
643  derv = (qeps - qnm) / eps
644  call matrix_sln%add_value_pos(idxglo(ii), derv)
645  rhs(n) = rhs(n) + derv * stage(m)
646 
647  end do
648  end do
649 
650  end subroutine dfw_qnm_fc_nr
651 
652  !> @brief fill newton
653  !<
654  subroutine dfw_fn(this, kiter, matrix_sln, idxglo, rhs, stage)
655  ! dummy
656  class(SwfDfwType) :: this !< this instance
657  integer(I4B) :: kiter
658  class(MatrixBaseType), pointer :: matrix_sln
659  integer(I4B), intent(in), dimension(:) :: idxglo
660  real(DP), intent(inout), dimension(:) :: rhs
661  real(DP), intent(inout), dimension(:) :: stage
662  ! local
663 
664  ! add newton terms to solution matrix
665  ! todo: add newton terms here instead?
666  ! this routine is probably not necessary as method is fully newton
667 
668  end subroutine dfw_fn
669 
670  !> @brief calculate flow between cells n and m
671  !<
672  function qcalc(this, n, m, stage_n, stage_m, ipos) result(qnm)
673  ! dummy
674  class(SwfDfwType) :: this !< this instance
675  integer(I4B), intent(in) :: n !< number for cell n
676  integer(I4B), intent(in) :: m !< number for cell m
677  real(DP), intent(in) :: stage_n !< stage in reach n
678  real(DP), intent(in) :: stage_m !< stage in reach m
679  integer(I4B), intent(in) :: ipos !< connection number
680  ! local
681  integer(I4B) :: isympos
682  real(DP) :: qnm
683  real(DP) :: cond
684  real(DP) :: cl1
685  real(DP) :: cl2
686 
687  ! Set connection lengths
688  isympos = this%dis%con%jas(ipos)
689  if (n < m) then
690  cl1 = this%dis%con%cl1(isympos)
691  cl2 = this%dis%con%cl2(isympos)
692  else
693  cl1 = this%dis%con%cl2(isympos)
694  cl2 = this%dis%con%cl1(isympos)
695  end if
696 
697  ! Calculate conductance
698  if (this%iswrcond == 0) then
699  cond = this%get_cond(n, m, ipos, stage_n, stage_m, cl1, cl2)
700  else if (this%iswrcond == 1) then
701  cond = this%get_cond_swr(n, m, ipos, stage_n, stage_m, cl1, cl2)
702  end if
703 
704  ! calculate flow between n and m
705  qnm = cond * (stage_m - stage_n)
706 
707  end function qcalc
708 
709  !> @brief calculate effective conductance between cells n and m
710  !!
711  !! Calculate half-cell conductances for cell n and cell m and then use
712  !! harmonic averaging to calculate the effective conductance between the
713  !< two cells.
714  function get_cond(this, n, m, ipos, stage_n, stage_m, cln, clm) result(cond)
715  ! modules
716  use smoothingmodule, only: squadratic
717  ! dummy
718  class(SwfDfwType) :: this !< this instance
719  integer(I4B), intent(in) :: n !< number for cell n
720  integer(I4B), intent(in) :: m !< number for cell m
721  integer(I4B), intent(in) :: ipos !< connection number
722  real(DP), intent(in) :: stage_n !< stage in reach n
723  real(DP), intent(in) :: stage_m !< stage in reach m
724  real(DP), intent(in) :: cln !< distance from cell n to shared face with m
725  real(DP), intent(in) :: clm !< distance from cell m to shared face with n
726  ! local
727  real(DP) :: depth_n
728  real(DP) :: depth_m
729  real(DP) :: dhds_n
730  real(DP) :: dhds_m
731  real(DP) :: width_n
732  real(DP) :: width_m
733  real(DP) :: range = 1.d-6
734  real(DP) :: dydx
735  real(DP) :: smooth_factor
736  real(DP) :: length_nm
737  real(DP) :: cond
738  real(DP) :: cn
739  real(DP) :: cm
740 
741  ! we are using a harmonic conductance approach here; however
742  ! the SWR Process for MODFLOW-2005/NWT uses length-weighted
743  ! average areas and hydraulic radius instead.
744  length_nm = cln + clm
745  cond = dzero
746  if (length_nm > dprec) then
747 
748  ! Calculate depth in each reach
749  depth_n = stage_n - this%dis%bot(n)
750  depth_m = stage_m - this%dis%bot(m)
751 
752  ! assign gradients
753  if (this%is2d == 0) then
754  dhds_n = abs(stage_m - stage_n) / (cln + clm)
755  dhds_m = dhds_n
756  else
757  dhds_n = this%grad_dhds_mag(n)
758  dhds_m = this%grad_dhds_mag(m)
759  end if
760 
761  ! Assign upstream depth, if not central
762  if (this%icentral == 0) then
763  ! use upstream weighting
764  if (stage_n > stage_m) then
765  depth_m = depth_n
766  else
767  depth_n = depth_m
768  end if
769  end if
770 
771  ! Calculate a smoothed depth that goes to zero over
772  ! the specified range
773  call squadratic(depth_n, range, dydx, smooth_factor)
774  depth_n = depth_n * smooth_factor
775  call squadratic(depth_m, range, dydx, smooth_factor)
776  depth_m = depth_m * smooth_factor
777 
778  ! Get the flow widths for n and m from dis package
779  call this%dis%get_flow_width(n, m, ipos, width_n, width_m)
780 
781  ! Calculate half-cell conductance for reach
782  ! n and m
783  cn = this%get_cond_n(n, depth_n, cln, width_n, dhds_n)
784  cm = this%get_cond_n(m, depth_m, clm, width_m, dhds_m)
785 
786  ! Use harmonic mean to calculate weighted
787  ! conductance between the centers of reaches
788  ! n and m
789  if ((cn + cm) > dprec) then
790  cond = cn * cm / (cn + cm)
791  else
792  cond = dzero
793  end if
794 
795  end if
796 
797  end function get_cond
798 
799  !> @brief Calculate half cell conductance
800  !!
801  !! Calculate half-cell conductance for cell n
802  !< using conveyance and Manning's equation
803  function get_cond_n(this, n, depth, dx, width, dhds) result(c)
804  ! modules
805  ! dummy
806  class(SwfDfwType) :: this !< this instance
807  integer(I4B), intent(in) :: n !< reach number
808  real(DP), intent(in) :: depth !< simulated depth (stage - elevation) in reach n for this iteration
809  real(DP), intent(in) :: dx !< half-cell distance
810  real(DP), intent(in) :: width !< width of the reach perpendicular to flow
811  real(DP), intent(in) :: dhds !< gradient
812  ! return
813  real(DP) :: c
814  ! local
815  real(DP) :: rough
816  real(DP) :: dhds_sqr
817  real(DP) :: conveyance
818 
819  ! Calculate conveyance, which is a * r**DTWOTHIRDS / roughc
820  rough = this%manningsn(n)
821  conveyance = this%cxs%get_conveyance(this%idcxs(n), width, depth, rough)
822  dhds_sqr = dhds**dhalf
823  if (dhds_sqr < dem10) then
824  dhds_sqr = dem10
825  end if
826 
827  ! Multiply by unitconv and divide conveyance by sqrt of friction slope and dx
828  c = this%unitconv * conveyance / dx / dhds_sqr
829 
830  end function get_cond_n
831 
832  !> @brief Calculate effective conductance for cells n and m using SWR method
833  !!
834  !! The SWR Process for MODFLOW uses average cell parameters from cell n and
835  !! m to calculate an effective conductance. This is different from the
836  !! default approach used in SWF, which uses harmonic averaging on two half-
837  !< cell conductances.
838  function get_cond_swr(this, n, m, ipos, stage_n, stage_m, cln, clm) result(cond)
839  ! modules
840  use smoothingmodule, only: squadratic
841  ! dummy
842  class(SwfDfwType) :: this !< this instance
843  integer(I4B), intent(in) :: n !< number for cell n
844  integer(I4B), intent(in) :: m !< number for cell m
845  integer(I4B), intent(in) :: ipos !< connection number
846  real(DP), intent(in) :: stage_n !< stage in reach n
847  real(DP), intent(in) :: stage_m !< stage in reach m
848  real(DP), intent(in) :: cln !< distance from cell n to shared face with m
849  real(DP), intent(in) :: clm !< distance from cell m to shared face with n
850  ! local
851  real(DP) :: depth_n
852  real(DP) :: depth_m
853  real(DP) :: dhds_n
854  real(DP) :: dhds_m
855  real(DP) :: dhds_nm
856  real(DP) :: dhds_sqr
857  real(DP) :: width_n
858  real(DP) :: width_m
859  real(DP) :: range = 1.d-6
860  real(DP) :: dydx
861  real(DP) :: smooth_factor
862  real(DP) :: length_nm
863  real(DP) :: cond
864  real(DP) :: ravg
865  real(DP) :: rinv_avg
866  real(DP) :: area_n, area_m, area_avg
867  real(DP) :: rhn, rhm, rhavg
868  real(DP) :: weight_n
869  real(DP) :: weight_m
870  real(DP) :: rough_n
871  real(DP) :: rough_m
872 
873  ! Use harmonic weighting for 1/manningsn, but using length-weighted
874  ! averaging for other terms
875  length_nm = cln + clm
876  cond = dzero
877  if (length_nm > dprec) then
878 
879  ! Calculate depth in each reach
880  depth_n = stage_n - this%dis%bot(n)
881  depth_m = stage_m - this%dis%bot(m)
882 
883  ! Assign upstream depth, if not central
884  if (this%icentral == 0) then
885  ! use upstream weighting
886  if (stage_n > stage_m) then
887  depth_m = depth_n
888  else
889  depth_n = depth_m
890  end if
891  end if
892 
893  ! Calculate a smoothed depth that goes to zero over
894  ! the specified range
895  call squadratic(depth_n, range, dydx, smooth_factor)
896  depth_n = depth_n * smooth_factor
897  call squadratic(depth_m, range, dydx, smooth_factor)
898  depth_m = depth_m * smooth_factor
899 
900  ! Get the flow widths for n and m from dis package
901  call this%dis%get_flow_width(n, m, ipos, width_n, width_m)
902 
903  ! linear weight toward node closer to shared face
904  weight_n = clm / length_nm
905  weight_m = done - weight_n
906 
907  ! average cross sectional flow area
908  area_n = this%cxs%get_area(this%idcxs(n), width_n, depth_n)
909  area_m = this%cxs%get_area(this%idcxs(m), width_m, depth_m)
910  area_avg = weight_n * area_n + weight_m * area_m
911 
912  ! average hydraulic radius
913  if (this%is2d == 0) then
914  rhn = this%cxs%get_hydraulic_radius(this%idcxs(n), width_n, &
915  depth_n, area_n)
916  rhm = this%cxs%get_hydraulic_radius(this%idcxs(m), width_m, &
917  depth_m, area_m)
918  rhavg = weight_n * rhn + weight_m * rhm
919  else
920  rhavg = area_avg / width_n
921  end if
922  rhavg = rhavg**dtwothirds
923 
924  ! average gradient
925  if (this%is2d == 0) then
926  dhds_nm = abs(stage_m - stage_n) / (length_nm)
927  else
928  dhds_n = this%grad_dhds_mag(n)
929  dhds_m = this%grad_dhds_mag(m)
930  dhds_nm = weight_n * dhds_n + weight_m * dhds_m
931  end if
932  dhds_sqr = dhds_nm**dhalf
933  if (dhds_sqr < dem10) then
934  dhds_sqr = dem10
935  end if
936 
937  ! weighted harmonic mean for inverse mannings value
938  weight_n = cln / length_nm
939  weight_m = done - weight_n
940  rough_n = this%cxs%get_roughness(this%idcxs(n), width_n, depth_n, &
941  this%manningsn(n), dhds_nm)
942  rough_m = this%cxs%get_roughness(this%idcxs(m), width_m, depth_m, &
943  this%manningsn(m), dhds_nm)
944  ravg = (weight_n + weight_m) / &
945  (weight_n / rough_n + weight_m / rough_m)
946  rinv_avg = done / ravg
947 
948  ! calculate conductance using averaged values
949  cond = this%unitconv * rinv_avg * area_avg * rhavg / dhds_sqr / length_nm
950 
951  end if
952 
953  end function get_cond_swr
954 
955  !> @brief Calculate flow area between cell n and m
956  !!
957  !! Calculate an average flow area between cell n and m.
958  !! First calculate a flow area for cell n and then for
959  !! cell m and linearly weight the areas using the connection
960  !< distances.
961  function get_flow_area_nm(this, n, m, stage_n, stage_m, cln, clm, &
962  ipos) result(area_avg)
963  ! module
964  use smoothingmodule, only: squadratic
965  ! dummy
966  class(SwfDfwType) :: this !< this instance
967  integer(I4B), intent(in) :: n
968  integer(I4B), intent(in) :: m
969  real(DP), intent(in) :: stage_n
970  real(DP), intent(in) :: stage_m
971  real(DP), intent(in) :: cln
972  real(DP), intent(in) :: clm
973  integer(I4B), intent(in) :: ipos
974  ! local
975  real(DP) :: depth_n
976  real(DP) :: depth_m
977  real(DP) :: width_n
978  real(DP) :: width_m
979  real(DP) :: area_n
980  real(DP) :: area_m
981  real(DP) :: weight_n
982  real(DP) :: weight_m
983  real(DP) :: length_nm
984  real(DP) :: range = 1.d-6
985  real(DP) :: dydx
986  real(DP) :: smooth_factor
987  ! return
988  real(DP) :: area_avg
989 
990  ! depths
991  depth_n = stage_n - this%dis%bot(n)
992  depth_m = stage_m - this%dis%bot(m)
993 
994  ! Assign upstream depth, if not central
995  if (this%icentral == 0) then
996  ! use upstream weighting
997  if (stage_n > stage_m) then
998  depth_m = depth_n
999  else
1000  depth_n = depth_m
1001  end if
1002  end if
1003 
1004  ! Calculate a smoothed depth that goes to zero over
1005  ! the specified range
1006  call squadratic(depth_n, range, dydx, smooth_factor)
1007  depth_n = depth_n * smooth_factor
1008  call squadratic(depth_m, range, dydx, smooth_factor)
1009  depth_m = depth_m * smooth_factor
1010 
1011  ! Get the flow widths for n and m from dis package
1012  call this%dis%get_flow_width(n, m, ipos, width_n, width_m)
1013 
1014  ! linear weight toward node closer to shared face
1015  length_nm = cln + clm
1016  weight_n = clm / length_nm
1017  weight_m = done - weight_n
1018 
1019  ! average cross sectional flow area
1020  area_n = this%cxs%get_area(this%idcxs(n), width_n, depth_n)
1021  area_m = this%cxs%get_area(this%idcxs(m), width_m, depth_m)
1022  area_avg = weight_n * area_n + weight_m * area_m
1023 
1024  end function get_flow_area_nm
1025 
1026  !> @brief Calculate average hydraulic gradient magnitude for each cell
1027  !!
1028  !! Go through each cell and calculate the average hydraulic gradient using
1029  !! an xt3d-style gradient interpolation. This is used for 2D grids in the
1030  !! calculation of an effective conductance. For 1D grids, gradients are
1031  !< calculated between cell centers.
1032  subroutine calc_dhds(this)
1033  ! modules
1035  ! dummy
1036  class(SwfDfwType) :: this !< this instance
1037  ! local
1038  integer(I4B) :: n
1039  integer(I4B) :: m
1040  integer(I4B) :: ipos
1041  integer(I4B) :: isympos
1042  real(DP) :: cl1
1043  real(DP) :: cl2
1044 
1045  do n = 1, this%dis%nodes
1046  this%grad_dhds_mag(n) = dzero
1047  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
1048  m = this%dis%con%ja(ipos)
1049  isympos = this%dis%con%jas(ipos)
1050 
1051  ! determine cl1 and cl2
1052  if (n < m) then
1053  cl1 = this%dis%con%cl1(isympos)
1054  cl2 = this%dis%con%cl2(isympos)
1055  else
1056  cl1 = this%dis%con%cl2(isympos)
1057  cl2 = this%dis%con%cl1(isympos)
1058  end if
1059 
1060  ! store for n < m in upper right triangular part of symmetric dhdsja array
1061  if (n < m) then
1062  if (cl1 + cl2 > dprec) then
1063  this%dhdsja(isympos) = (this%hnew(m) - this%hnew(n)) / (cl1 + cl2)
1064  else
1065  this%dhdsja(isympos) = dzero
1066  end if
1067  end if
1068  end do
1069  end do
1070 
1071  ! pass dhdsja into the vector interpolation to get the components
1072  ! of the gradient at the cell center
1073  call vector_interpolation_2d(this%dis, this%dhdsja, vmag=this%grad_dhds_mag)
1074 
1075  end subroutine calc_dhds
1076 
1077  !> @ brief Newton under relaxation
1078  !!
1079  !< If stage is below the bottom, then pull it up a bit
1080  subroutine dfw_nur(this, neqmod, x, xtemp, dx, inewtonur, dxmax, locmax)
1081  ! dummy
1082  class(SwfDfwType) :: this !< this instance
1083  integer(I4B), intent(in) :: neqmod !< number of equations
1084  real(DP), dimension(neqmod), intent(inout) :: x !< dependent variable
1085  real(DP), dimension(neqmod), intent(in) :: xtemp !< temporary dependent variable
1086  real(DP), dimension(neqmod), intent(inout) :: dx !< change in dependent variable
1087  integer(I4B), intent(inout) :: inewtonur !< flag to indication relaxation was applied
1088  real(DP), intent(inout) :: dxmax !< max change in x
1089  integer(I4B), intent(inout) :: locmax !< location of max change
1090  ! local
1091  integer(I4B) :: n
1092  real(DP) :: botm
1093  real(DP) :: xx
1094  real(DP) :: dxx
1095 
1096  ! Newton-Raphson under-relaxation
1097  do n = 1, this%dis%nodes
1098  if (this%ibound(n) < 1) cycle
1099  if (this%icelltype(n) > 0) then
1100  botm = this%dis%bot(n)
1101  ! only apply Newton-Raphson under-relaxation if
1102  ! solution head is below the bottom of the model
1103  if (x(n) < botm) then
1104  inewtonur = 1
1105  xx = xtemp(n) * (done - dp9) + botm * dp9
1106  dxx = x(n) - xx
1107  if (abs(dxx) > abs(dxmax)) then
1108  locmax = n
1109  dxmax = dxx
1110  end if
1111  x(n) = xx
1112  dx(n) = dzero
1113  end if
1114  end if
1115  end do
1116 
1117  end subroutine dfw_nur
1118 
1119  !> @ brief Calculate flow for each connection and store in flowja
1120  !<
1121  subroutine dfw_cq(this, stage, stage_old, flowja)
1122  ! dummy
1123  class(SwfDfwType) :: this !< this instance
1124  real(DP), intent(inout), dimension(:) :: stage !< calculated head
1125  real(DP), intent(inout), dimension(:) :: stage_old !< calculated head from previous time step
1126  real(DP), intent(inout), dimension(:) :: flowja !< vector of flows in CSR format
1127  ! local
1128  integer(I4B) :: n, ipos, m
1129  real(DP) :: qnm
1130 
1131  do n = 1, this%dis%nodes
1132  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
1133  m = this%dis%con%ja(ipos)
1134  if (m < n) cycle
1135  qnm = this%qcalc(n, m, stage(n), stage(m), ipos)
1136  flowja(ipos) = qnm
1137  flowja(this%dis%con%isym(ipos)) = -qnm
1138  end do
1139  end do
1140 
1141  end subroutine dfw_cq
1142 
1143  !> @ brief Model budget calculation for package
1144  !<
1145  subroutine dfw_bd(this, isuppress_output, model_budget)
1146  ! modules
1147  use budgetmodule, only: budgettype
1148  ! dummy variables
1149  class(SwfDfwType) :: this !< this instance
1150  integer(I4B), intent(in) :: isuppress_output !< flag to suppress model output
1151  type(BudgetType), intent(inout) :: model_budget !< model budget object
1152 
1153  ! Add any DFW budget terms
1154  ! none
1155 
1156  end subroutine dfw_bd
1157 
1158  !> @ brief save flows for package
1159  !<
1160  subroutine dfw_save_model_flows(this, flowja, icbcfl, icbcun)
1161  ! dummy
1162  class(SwfDfwType) :: this !< this instance
1163  real(DP), dimension(:), intent(in) :: flowja !< vector of flows in CSR format
1164  integer(I4B), intent(in) :: icbcfl !< flag to indicate if flows should be saved
1165  integer(I4B), intent(in) :: icbcun !< unit number for flow output
1166  ! local
1167  integer(I4B) :: ibinun
1168 
1169  ! Set unit number for binary output
1170  if (this%ipakcb < 0) then
1171  ibinun = icbcun
1172  elseif (this%ipakcb == 0) then
1173  ibinun = 0
1174  else
1175  ibinun = this%ipakcb
1176  end if
1177  if (icbcfl == 0) ibinun = 0
1178 
1179  ! Write the face flows if requested
1180  if (ibinun /= 0) then
1181  ! flowja
1182  call this%dis%record_connection_array(flowja, ibinun, this%iout)
1183  end if
1184 
1185  ! Calculate velocities at cell centers and write, if requested
1186  if (this%isavvelocity /= 0) then
1187  if (ibinun /= 0) call this%sav_velocity(ibinun)
1188  end if
1189 
1190  end subroutine dfw_save_model_flows
1191 
1192  !> @ brief print flows for package
1193  !<
1194  subroutine dfw_print_model_flows(this, ibudfl, flowja)
1195  ! modules
1196  use tdismodule, only: kper, kstp
1197  use constantsmodule, only: lenbigline
1198  ! dummy
1199  class(SwfDfwType) :: this !< this instance
1200  integer(I4B), intent(in) :: ibudfl !< print flag
1201  real(DP), intent(inout), dimension(:) :: flowja !< vector of flows in CSR format
1202  ! local
1203  character(len=LENBIGLINE) :: line
1204  character(len=30) :: tempstr
1205  integer(I4B) :: n, ipos, m
1206  real(DP) :: qnm
1207  ! formats
1208  character(len=*), parameter :: fmtiprflow = &
1209  &"(/,4x,'CALCULATED INTERCELL FLOW FOR PERIOD ', i0, ' STEP ', i0)"
1210 
1211  ! Write flowja to list file if requested
1212  if (ibudfl /= 0 .and. this%iprflow > 0) then
1213  write (this%iout, fmtiprflow) kper, kstp
1214  do n = 1, this%dis%nodes
1215  line = ''
1216  call this%dis%noder_to_string(n, tempstr)
1217  line = trim(tempstr)//':'
1218  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
1219  m = this%dis%con%ja(ipos)
1220  call this%dis%noder_to_string(m, tempstr)
1221  line = trim(line)//' '//trim(tempstr)
1222  qnm = flowja(ipos)
1223  write (tempstr, '(1pg15.6)') qnm
1224  line = trim(line)//' '//trim(adjustl(tempstr))
1225  end do
1226  write (this%iout, '(a)') trim(line)
1227  end do
1228  end if
1229 
1230  end subroutine dfw_print_model_flows
1231 
1232  !> @brief deallocate memory
1233  !<
1234  subroutine dfw_da(this)
1235  ! modules
1238  use simvariablesmodule, only: idm_context
1239  ! dummy
1240  class(SwfDfwType) :: this !< this instance
1241 
1242  ! Deallocate input memory
1243  call memorystore_remove(this%name_model, 'DFW', idm_context)
1244 
1245  ! Deallocate arrays
1246  call mem_deallocate(this%manningsn)
1247  call mem_deallocate(this%idcxs)
1248  call mem_deallocate(this%icelltype)
1249  call mem_deallocate(this%nodedge)
1250  call mem_deallocate(this%ihcedge)
1251  call mem_deallocate(this%propsedge)
1252  call mem_deallocate(this%vcomp)
1253  call mem_deallocate(this%vmag)
1254  if (this%is2d == 1) then
1255  call mem_deallocate(this%grad_dhds_mag)
1256  call mem_deallocate(this%dhdsja)
1257  end if
1258 
1259  ! Scalars
1260  call mem_deallocate(this%is2d)
1261  call mem_deallocate(this%icentral)
1262  call mem_deallocate(this%iswrcond)
1263  call mem_deallocate(this%unitconv)
1264  call mem_deallocate(this%lengthconv)
1265  call mem_deallocate(this%timeconv)
1266  call mem_deallocate(this%isavvelocity)
1267  call mem_deallocate(this%icalcvelocity)
1268  call mem_deallocate(this%nedges)
1269  call mem_deallocate(this%lastedge)
1270 
1271  ! obs package
1272  call mem_deallocate(this%inobspkg)
1273  call this%obs%obs_da()
1274  deallocate (this%obs)
1275  nullify (this%obs)
1276  nullify (this%cxs)
1277 
1278  ! deallocate parent
1279  call this%NumericalPackageType%da()
1280 
1281  ! pointers
1282  this%hnew => null()
1283 
1284  end subroutine dfw_da
1285 
1286  !> @brief Calculate the 3 components of velocity at the cell center
1287  !!
1288  !! todo: duplicated from NPF; should consolidate
1289  !<
1290  subroutine calc_velocity(this, flowja)
1291  ! modules
1292  ! dummy
1293  class(SwfDfwType) :: this !< this instance
1294  real(DP), intent(in), dimension(:) :: flowja !< vector of flows in CSR format
1295  ! local
1296  integer(I4B) :: n
1297  integer(I4B) :: m
1298  integer(I4B) :: ipos
1299  integer(I4B) :: isympos
1300  integer(I4B) :: ihc
1301  integer(I4B) :: ic
1302  integer(I4B) :: iz
1303  integer(I4B) :: nc
1304  integer(I4B) :: ncz
1305  real(DP) :: vx
1306  real(DP) :: vy
1307  real(DP) :: vz
1308  real(DP) :: xn
1309  real(DP) :: yn
1310  real(DP) :: zn
1311  real(DP) :: xc
1312  real(DP) :: yc
1313  real(DP) :: zc
1314  real(DP) :: cl1
1315  real(DP) :: cl2
1316  real(DP) :: dltot
1317  real(DP) :: ooclsum
1318  real(DP) :: dsumx
1319  real(DP) :: dsumy
1320  real(DP) :: dsumz
1321  real(DP) :: denom
1322  real(DP) :: area
1323  real(DP) :: axy
1324  real(DP) :: ayx
1325  real(DP), allocatable, dimension(:) :: vi
1326  real(DP), allocatable, dimension(:) :: di
1327  real(DP), allocatable, dimension(:) :: viz
1328  real(DP), allocatable, dimension(:) :: diz
1329  real(DP), allocatable, dimension(:) :: nix
1330  real(DP), allocatable, dimension(:) :: niy
1331  real(DP), allocatable, dimension(:) :: wix
1332  real(DP), allocatable, dimension(:) :: wiy
1333  real(DP), allocatable, dimension(:) :: wiz
1334  real(DP), allocatable, dimension(:) :: bix
1335  real(DP), allocatable, dimension(:) :: biy
1336  logical :: nozee = .true.
1337 
1338  ! Ensure dis has necessary information
1339  ! todo: do we need this for SWF?
1340  if (this%icalcvelocity /= 0 .and. this%dis%con%ianglex == 0) then
1341  call store_error('Error. ANGLDEGX not provided in '// &
1342  'discretization file. ANGLDEGX required for '// &
1343  'calculation of velocity.', terminate=.true.)
1344  end if
1345 
1346  ! Find max number of connections and allocate weight arrays
1347  nc = 0
1348  do n = 1, this%dis%nodes
1349 
1350  ! Count internal model connections
1351  ic = this%dis%con%ia(n + 1) - this%dis%con%ia(n) - 1
1352 
1353  ! Count edge connections
1354  do m = 1, this%nedges
1355  if (this%nodedge(m) == n) then
1356  ic = ic + 1
1357  end if
1358  end do
1359 
1360  ! Set max number of connections for any cell
1361  if (ic > nc) nc = ic
1362  end do
1363 
1364  ! Allocate storage arrays needed for cell-centered calculation
1365  allocate (vi(nc))
1366  allocate (di(nc))
1367  allocate (viz(nc))
1368  allocate (diz(nc))
1369  allocate (nix(nc))
1370  allocate (niy(nc))
1371  allocate (wix(nc))
1372  allocate (wiy(nc))
1373  allocate (wiz(nc))
1374  allocate (bix(nc))
1375  allocate (biy(nc))
1376 
1377  ! Go through each cell and calculate specific discharge
1378  do n = 1, this%dis%nodes
1379 
1380  ! first calculate geometric properties for x and y directions and
1381  ! the specific discharge at a face (vi)
1382  ic = 0
1383  iz = 0
1384  vi(:) = dzero
1385  di(:) = dzero
1386  viz(:) = dzero
1387  diz(:) = dzero
1388  nix(:) = dzero
1389  niy(:) = dzero
1390  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
1391  m = this%dis%con%ja(ipos)
1392  isympos = this%dis%con%jas(ipos)
1393  ihc = this%dis%con%ihc(isympos)
1394  ic = ic + 1
1395  call this%dis%connection_normal(n, m, ihc, xn, yn, zn, ipos)
1396  call this%dis%connection_vector(n, m, nozee, done, done, &
1397  ihc, xc, yc, zc, dltot)
1398  cl1 = this%dis%con%cl1(isympos)
1399  cl2 = this%dis%con%cl2(isympos)
1400  if (m < n) then
1401  cl1 = this%dis%con%cl2(isympos)
1402  cl2 = this%dis%con%cl1(isympos)
1403  end if
1404  ooclsum = done / (cl1 + cl2)
1405  nix(ic) = -xn
1406  niy(ic) = -yn
1407  di(ic) = dltot * cl1 * ooclsum
1408  area = this%get_flow_area_nm(n, m, this%hnew(n), this%hnew(m), &
1409  cl1, cl2, ipos)
1410  if (area > dzero) then
1411  vi(ic) = flowja(ipos) / area
1412  else
1413  vi(ic) = dzero
1414  end if
1415 
1416  end do
1417 
1418  ! Look through edge flows that may have been provided by an exchange
1419  ! and incorporate them into the averaging arrays
1420  do m = 1, this%nedges
1421  if (this%nodedge(m) == n) then
1422 
1423  ! propsedge: (Q, area, nx, ny, distance)
1424  ihc = this%ihcedge(m)
1425  area = this%propsedge(2, m)
1426 
1427  ic = ic + 1
1428  nix(ic) = -this%propsedge(3, m)
1429  niy(ic) = -this%propsedge(4, m)
1430  di(ic) = this%propsedge(5, m)
1431  if (area > dzero) then
1432  vi(ic) = this%propsedge(1, m) / area
1433  else
1434  vi(ic) = dzero
1435  end if
1436 
1437  end if
1438  end do
1439 
1440  ! Assign number of vertical and horizontal connections
1441  ncz = iz
1442  nc = ic
1443 
1444  ! calculate z weight (wiz) and z velocity
1445  if (ncz == 1) then
1446  wiz(1) = done
1447  else
1448  dsumz = dzero
1449  do iz = 1, ncz
1450  dsumz = dsumz + diz(iz)
1451  end do
1452  denom = (ncz - done)
1453  if (denom < dzero) denom = dzero
1454  dsumz = dsumz + dem10 * dsumz
1455  do iz = 1, ncz
1456  if (dsumz > dzero) wiz(iz) = done - diz(iz) / dsumz
1457  if (denom > 0) then
1458  wiz(iz) = wiz(iz) / denom
1459  else
1460  wiz(iz) = dzero
1461  end if
1462  end do
1463  end if
1464  vz = dzero
1465  do iz = 1, ncz
1466  vz = vz + wiz(iz) * viz(iz)
1467  end do
1468 
1469  ! distance-based weighting
1470  nc = ic
1471  dsumx = dzero
1472  dsumy = dzero
1473  dsumz = dzero
1474  do ic = 1, nc
1475  wix(ic) = di(ic) * abs(nix(ic))
1476  wiy(ic) = di(ic) * abs(niy(ic))
1477  dsumx = dsumx + wix(ic)
1478  dsumy = dsumy + wiy(ic)
1479  end do
1480 
1481  ! Finish computing omega weights. Add a tiny bit
1482  ! to dsum so that the normalized omega weight later
1483  ! evaluates to (essentially) 1 in the case of a single
1484  ! relevant connection, avoiding 0/0.
1485  dsumx = dsumx + dem10 * dsumx
1486  dsumy = dsumy + dem10 * dsumy
1487  do ic = 1, nc
1488  wix(ic) = (dsumx - wix(ic)) * abs(nix(ic))
1489  wiy(ic) = (dsumy - wiy(ic)) * abs(niy(ic))
1490  end do
1491 
1492  ! compute B weights
1493  dsumx = dzero
1494  dsumy = dzero
1495  do ic = 1, nc
1496  bix(ic) = wix(ic) * sign(done, nix(ic))
1497  biy(ic) = wiy(ic) * sign(done, niy(ic))
1498  dsumx = dsumx + wix(ic) * abs(nix(ic))
1499  dsumy = dsumy + wiy(ic) * abs(niy(ic))
1500  end do
1501  if (dsumx > dzero) dsumx = done / dsumx
1502  if (dsumy > dzero) dsumy = done / dsumy
1503  axy = dzero
1504  ayx = dzero
1505  do ic = 1, nc
1506  bix(ic) = bix(ic) * dsumx
1507  biy(ic) = biy(ic) * dsumy
1508  axy = axy + bix(ic) * niy(ic)
1509  ayx = ayx + biy(ic) * nix(ic)
1510  end do
1511 
1512  ! Calculate specific discharge. The divide by zero checking below
1513  ! is problematic for cells with only one flow, such as can happen
1514  ! with triangular cells in corners. In this case, the resulting
1515  ! cell velocity will be calculated as zero. The method should be
1516  ! improved so that edge flows of zero are included in these
1517  ! calculations. But this needs to be done with consideration for LGR
1518  ! cases in which flows are submitted from an exchange.
1519  vx = dzero
1520  vy = dzero
1521  do ic = 1, nc
1522  vx = vx + (bix(ic) - axy * biy(ic)) * vi(ic)
1523  vy = vy + (biy(ic) - ayx * bix(ic)) * vi(ic)
1524  end do
1525  denom = done - axy * ayx
1526  if (denom /= dzero) then
1527  vx = vx / denom
1528  vy = vy / denom
1529  end if
1530 
1531  this%vcomp(1, n) = vx
1532  this%vcomp(2, n) = vy
1533  this%vcomp(3, n) = vz
1534  this%vmag(n) = sqrt(vx**2 + vy**2 + vz**2)
1535 
1536  end do
1537 
1538  ! cleanup
1539  deallocate (vi)
1540  deallocate (di)
1541  deallocate (nix)
1542  deallocate (niy)
1543  deallocate (wix)
1544  deallocate (wiy)
1545  deallocate (wiz)
1546  deallocate (bix)
1547  deallocate (biy)
1548 
1549  end subroutine calc_velocity
1550 
1551  !> @brief Reserve space for nedges cells that have an edge on them.
1552  !!
1553  !! todo: duplicated from NPF; should consolidate
1554  !! This must be called before the swf%allocate_arrays routine, which is
1555  !< called from swf%ar.
1556  subroutine increase_edge_count(this, nedges)
1557  ! dummy
1558  class(SwfDfwType) :: this !< this instance
1559  integer(I4B), intent(in) :: nedges
1560 
1561  this%nedges = this%nedges + nedges
1562 
1563  end subroutine increase_edge_count
1564 
1565  !> @brief Provide the swf package with edge properties
1566  !!
1567  !! todo: duplicated from NPF; should consolidate
1568  !<
1569  subroutine set_edge_properties(this, nodedge, ihcedge, q, area, nx, ny, &
1570  distance)
1571  ! dummy
1572  class(SwfDfwType) :: this !< this instance
1573  integer(I4B), intent(in) :: nodedge
1574  integer(I4B), intent(in) :: ihcedge
1575  real(DP), intent(in) :: q
1576  real(DP), intent(in) :: area
1577  real(DP), intent(in) :: nx
1578  real(DP), intent(in) :: ny
1579  real(DP), intent(in) :: distance
1580  ! local
1581  integer(I4B) :: lastedge
1582 
1583  this%lastedge = this%lastedge + 1
1584  lastedge = this%lastedge
1585  this%nodedge(lastedge) = nodedge
1586  this%ihcedge(lastedge) = ihcedge
1587  this%propsedge(1, lastedge) = q
1588  this%propsedge(2, lastedge) = area
1589  this%propsedge(3, lastedge) = nx
1590  this%propsedge(4, lastedge) = ny
1591  this%propsedge(5, lastedge) = distance
1592 
1593  ! If this is the last edge, then the next call must be starting a new
1594  ! edge properties assignment loop, so need to reset lastedge to 0
1595  if (this%lastedge == this%nedges) this%lastedge = 0
1596 
1597  end subroutine set_edge_properties
1598 
1599  !> @brief Save specific discharge in binary format to ibinun
1600  !!
1601  !! todo: should write 2d velocity; what about for 1D channel?
1602  !<
1603  subroutine sav_velocity(this, ibinun)
1604  ! dummy
1605  class(SwfDfwType) :: this !< this instance
1606  integer(I4B), intent(in) :: ibinun
1607  ! local
1608  character(len=16) :: text
1609  character(len=16), dimension(3) :: auxtxt
1610  integer(I4B) :: n
1611  integer(I4B) :: naux
1612 
1613  ! Write the header
1614  text = ' DATA-VCOMP'
1615  naux = 3
1616  auxtxt(:) = [' vx', ' vy', ' vz']
1617  call this%dis%record_srcdst_list_header(text, this%name_model, &
1618  this%packName, this%name_model, &
1619  this%packName, naux, auxtxt, ibinun, &
1620  this%dis%nodes, this%iout)
1621 
1622  ! Write a zero for Q, and then write qx, qy, qz as aux variables
1623  do n = 1, this%dis%nodes
1624  call this%dis%record_mf6_list_entry(ibinun, n, n, dzero, naux, &
1625  this%vcomp(:, n))
1626  end do
1627 
1628  end subroutine sav_velocity
1629 
1630  !> @brief Define the observation types available in the package
1631  !!
1632  !< Method to define the observation types available in the package.
1633  subroutine dfw_df_obs(this)
1634  ! dummy variables
1635  class(SwfDfwType) :: this !< this instance
1636  ! local variables
1637  integer(I4B) :: indx
1638 
1639  ! Store obs type and assign procedure pointer
1640  ! for ext-outflow observation type.
1641  call this%obs%StoreObsType('ext-outflow', .true., indx)
1642  this%obs%obsData(indx)%ProcessIdPtr => dfwobsidprocessor
1643 
1644  end subroutine dfw_df_obs
1645 
1646  subroutine dfwobsidprocessor(obsrv, dis, inunitobs, iout)
1647  ! dummy
1648  type(ObserveType), intent(inout) :: obsrv
1649  class(DisBaseType), intent(in) :: dis
1650  integer(I4B), intent(in) :: inunitobs
1651  integer(I4B), intent(in) :: iout
1652  ! local
1653  integer(I4B) :: n
1654  character(len=LINELENGTH) :: string
1655 
1656  ! Initialize variables
1657  string = obsrv%IDstring
1658  read (string, *) n
1659 
1660  if (n > 0) then
1661  obsrv%NodeNumber = n
1662  else
1663  errmsg = 'Error reading data from ID string'
1664  call store_error(errmsg)
1665  call store_error_unit(inunitobs)
1666  end if
1667 
1668  end subroutine dfwobsidprocessor
1669 
1670  !> @brief Save observations for the package
1671  !!
1672  !< Method to save simulated values for the package.
1673  subroutine dfw_bd_obs(this)
1674  ! dummy variables
1675  class(SwfDfwType) :: this !< this instance
1676  ! local variables
1677  integer(I4B) :: i
1678  integer(I4B) :: j
1679  integer(I4B) :: n
1680  real(DP) :: v
1681  character(len=100) :: msg
1682  type(ObserveType), pointer :: obsrv => null()
1683 
1684  ! Write simulated values for all observations
1685  if (this%obs%npakobs > 0) then
1686  call this%obs%obs_bd_clear()
1687  do i = 1, this%obs%npakobs
1688  obsrv => this%obs%pakobs(i)%obsrv
1689  do j = 1, obsrv%indxbnds_count
1690  n = obsrv%indxbnds(j)
1691  v = dzero
1692  select case (obsrv%ObsTypeId)
1693  case default
1694  msg = 'Unrecognized observation type: '//trim(obsrv%ObsTypeId)
1695  call store_error(msg)
1696  end select
1697  call this%obs%SaveOneSimval(obsrv, v)
1698  end do
1699  end do
1700 
1701  ! write summary of package error messages
1702  if (count_errors() > 0) then
1703  call store_error_filename(this%input_fname)
1704  end if
1705  end if
1706 
1707  end subroutine dfw_bd_obs
1708 
1709  !> @brief Read and prepare observations for a package
1710  !!
1711  !< Method to read and prepare observations for a package.
1712  subroutine dfw_rp_obs(this)
1713  ! modules
1714  use tdismodule, only: kper
1715  ! dummy
1716  class(SwfDfwType), intent(inout) :: this !< this instance
1717  ! local
1718  integer(I4B) :: i
1719  integer(I4B) :: j
1720  integer(I4B) :: nn1
1721  class(ObserveType), pointer :: obsrv => null()
1722 
1723  ! process each package observation
1724  ! only done the first stress period since boundaries are fixed
1725  ! for the simulation
1726  if (kper == 1) then
1727  do i = 1, this%obs%npakobs
1728  obsrv => this%obs%pakobs(i)%obsrv
1729 
1730  ! get node number 1
1731  nn1 = obsrv%NodeNumber
1732  if (nn1 < 1 .or. nn1 > this%dis%nodes) then
1733  write (errmsg, '(a,1x,a,1x,i0,1x,a,1x,i0,a)') &
1734  trim(adjustl(obsrv%ObsTypeId)), &
1735  'reach must be greater than 0 and less than or equal to', &
1736  this%dis%nodes, '(specified value is ', nn1, ')'
1737  call store_error(errmsg)
1738  else
1739  if (obsrv%indxbnds_count == 0) then
1740  call obsrv%AddObsIndex(nn1)
1741  else
1742  errmsg = 'Programming error in dfw_rp_obs'
1743  call store_error(errmsg)
1744  end if
1745  end if
1746 
1747  ! check that node number 1 is valid; call store_error if not
1748  do j = 1, obsrv%indxbnds_count
1749  nn1 = obsrv%indxbnds(j)
1750  if (nn1 < 1 .or. nn1 > this%dis%nodes) then
1751  write (errmsg, '(a,1x,a,1x,i0,1x,a,1x,i0,a)') &
1752  trim(adjustl(obsrv%ObsTypeId)), &
1753  'reach must be greater than 0 and less than or equal to', &
1754  this%dis%nodes, '(specified value is ', nn1, ')'
1755  call store_error(errmsg)
1756  end if
1757  end do
1758  end do
1759 
1760  ! evaluate if there are any observation errors
1761  if (count_errors() > 0) then
1762  call store_error_filename(this%input_fname)
1763  end if
1764 
1765  end if
1766 
1767  end subroutine dfw_rp_obs
1768 
1769 end module swfdfwmodule
This module contains the BudgetModule.
Definition: Budget.f90:20
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 dtwothirds
real constant 2/3
Definition: Constants.f90:70
real(dp), parameter dp9
real constant 9/10
Definition: Constants.f90:72
real(dp), parameter dem10
real constant 1e-10
Definition: Constants.f90:113
real(dp), parameter donethird
real constant 1/3
Definition: Constants.f90:67
real(dp), parameter dnodata
real no data constant
Definition: Constants.f90:95
integer(i4b), parameter lenbigline
maximum length of a big line
Definition: Constants.f90:15
integer(i4b), parameter lenvarname
maximum length of a variable name
Definition: Constants.f90:17
real(dp), parameter dhalf
real constant 1/2
Definition: Constants.f90:68
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
real(dp), parameter dem5
real constant 1e-5
Definition: Constants.f90:108
real(dp), parameter dprec
real constant machine precision
Definition: Constants.f90:120
real(dp), parameter dtwo
real constant 2
Definition: Constants.f90:79
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:27
real(dp), parameter dthree
real constant 3
Definition: Constants.f90:80
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
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
real(dp) function, public get_perturbation(x)
Calculate a numerical perturbation given the value of x.
Definition: MathUtil.f90:372
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
subroutine, public memorystore_remove(component, subcomponent, context)
subroutine, public get_isize(name, mem_path, isize)
@ brief Get the number of elements for this variable
This module contains the base numerical package type.
This module contains the derived types ObserveType and ObsDataType.
Definition: Observe.f90:15
This module contains the derived type ObsType.
Definition: Obs.f90:127
subroutine, public obs_cr(obs, inobs)
@ brief Create a new ObsType object
Definition: Obs.f90:225
subroutine, public defaultobsidprocessor(obsrv, dis, inunitobs, iout)
@ brief Process IDstring provided for each observation
Definition: Obs.f90:246
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_filename(filename, terminate)
Store the erroring file name.
Definition: Sim.f90:203
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
character(len=linelength) idm_context
character(len=maxcharlen) warnmsg
warning message string
subroutine squadratic(x, range, dydx, y)
@ brief sQuadratic
integer(i4b), pointer, public kstp
current time step number
Definition: tdis.f90:24
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
subroutine, public vector_interpolation_2d(dis, flowja, nedges, nodedge, propsedge, vcomp, vmag, flowareaja)
Interpolate 2D vector components at cell center.
Derived type for the Budget object.
Definition: Budget.f90:39
This class is used to store a single deferred-length character string. It was designed to work in an ...
Definition: CharString.f90:23