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