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