27 public :: swfdfwtype, dfw_cr
32 integer(I4B),
pointer :: is2d => null()
33 integer(I4B),
pointer :: icentral => null()
34 integer(I4B),
pointer :: iswrcond => null()
35 real(DP),
pointer :: unitconv
36 real(DP),
pointer :: timeconv
37 real(DP),
pointer :: lengthconv
38 real(DP),
dimension(:),
pointer,
contiguous :: hnew => null()
39 real(DP),
dimension(:),
pointer,
contiguous :: manningsn => null()
40 integer(I4B),
dimension(:),
pointer,
contiguous :: idcxs => null()
41 integer(I4B),
dimension(:),
pointer,
contiguous :: ibound => null()
42 integer(I4B),
dimension(:),
pointer,
contiguous :: icelltype => null()
45 integer(I4B),
pointer :: icalcvelocity => null()
46 integer(I4B),
pointer :: isavvelocity => null()
47 real(DP),
dimension(:, :),
pointer,
contiguous :: vcomp => null()
48 real(DP),
dimension(:),
pointer,
contiguous :: vmag => null()
49 integer(I4B),
pointer :: nedges => null()
50 integer(I4B),
pointer :: lastedge => null()
51 integer(I4B),
dimension(:),
pointer,
contiguous :: nodedge => null()
52 integer(I4B),
dimension(:),
pointer,
contiguous :: ihcedge => null()
53 real(DP),
dimension(:, :),
pointer,
contiguous :: propsedge => null()
54 real(DP),
dimension(:),
pointer,
contiguous :: grad_dhds_mag => null()
55 real(DP),
dimension(:),
pointer,
contiguous :: dhdsja => null()
58 integer(I4B),
pointer :: inobspkg => null()
59 type(ObsType),
pointer :: obs => null()
62 type(SwfCxsType),
pointer :: cxs
67 procedure :: allocate_scalars
68 procedure :: allocate_arrays
70 procedure :: source_options
71 procedure :: log_options
72 procedure :: source_griddata
73 procedure :: log_griddata
78 procedure :: dfw_qnm_fc_nr
84 procedure :: dfw_save_model_flows
85 procedure :: dfw_print_model_flows
87 procedure :: dfw_df_obs
88 procedure :: dfw_rp_obs
89 procedure :: dfw_bd_obs
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
108 subroutine dfw_cr(dfwobj, name_model, input_mempath, inunit, iout, &
113 type(SwfDfwType),
pointer :: dfwobj
114 character(len=*),
intent(in) :: name_model
115 character(len=*),
intent(in) :: input_mempath
116 integer(I4B),
intent(in) :: inunit
117 integer(I4B),
intent(in) :: iout
118 type(SwfCxsType),
pointer,
intent(in) :: cxs
120 logical(LGP) :: found_fname
122 character(len=*),
parameter :: fmtheader = &
123 "(1x, /1x, 'DFW -- DIFFUSIVE WAVE (DFW) PACKAGE, VERSION 1, 9/25/2023', &
124 &' INPUT READ FROM MEMPATH: ', A, /)"
130 call dfwobj%set_names(1, name_model,
'DFW',
'DFW')
133 call dfwobj%allocate_scalars()
136 dfwobj%input_mempath = input_mempath
137 dfwobj%inunit = inunit
141 call mem_set_value(dfwobj%input_fname,
'INPUT_FNAME', dfwobj%input_mempath, &
148 call obs_cr(dfwobj%obs, dfwobj%inobspkg)
154 write (iout, fmtheader) input_mempath
158 end subroutine dfw_cr
162 subroutine dfw_df(this, dis)
164 class(SwfDfwType) :: this
165 class(DisBaseType),
pointer,
intent(inout) :: dis
171 if (this%dis%is_2d())
then
180 call this%allocate_arrays()
187 end subroutine dfw_df
195 subroutine allocate_scalars(this)
198 class(SwfDfwtype) :: this
201 call this%NumericalPackageType%allocate_scalars()
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)
220 this%lengthconv =
done
223 this%icalcvelocity = 0
224 this%isavvelocity = 0
228 end subroutine allocate_scalars
232 subroutine allocate_arrays(this)
234 class(SwfDfwType) :: this
240 'MANNINGSN', this%memoryPath)
242 'IDCXS', this%memoryPath)
244 'ICELLTYPE', this%memoryPath)
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)
252 call mem_allocate(this%vcomp, 3, 0,
'VCOMP', this%memoryPath)
253 call mem_allocate(this%vmag, 0,
'VMAG', this%memoryPath)
255 do n = 1, this%dis%nodes
256 this%manningsn(n) =
dzero
258 this%icelltype(n) = 1
262 if (this%is2d == 1)
then
264 'GRAD_DHDS_MAG', this%memoryPath)
266 'DHDSJA', this%memoryPath)
267 do n = 1, this%dis%nodes
268 this%grad_dhds_mag(n) =
dzero
270 do n = 1, this%dis%njas
271 this%dhdsja(n) =
dzero
275 end subroutine allocate_arrays
279 subroutine dfw_load(this)
281 class(SwfDfwType) :: this
284 call this%source_options()
285 call this%source_griddata()
287 end subroutine dfw_load
291 subroutine source_options(this)
299 class(SwfDfwType) :: this
301 integer(I4B) :: isize
302 type(SwfDfwParamFoundType) :: found
303 type(CharacterStringType),
dimension(:),
pointer, &
304 contiguous :: obs6_fnames
308 this%input_mempath, found%icentral)
310 this%input_mempath, found%iswrcond)
312 this%input_mempath, found%lengthconv)
314 this%input_mempath, found%timeconv)
316 this%input_mempath, found%iprflow)
318 this%input_mempath, found%ipakcb)
320 this%input_mempath, found%isavvelocity)
323 if (found%icentral) this%icentral = 1
324 if (found%ipakcb) this%ipakcb = -1
327 this%unitconv = this%lengthconv**
donethird
328 this%unitconv = this%unitconv * this%timeconv
331 if (found%isavvelocity) this%icalcvelocity = this%isavvelocity
334 call get_isize(
'OBS6_FILENAME', this%input_mempath, isize)
338 errmsg =
'Multiple OBS6 keywords detected in OPTIONS block.'// &
339 ' Only one OBS6 entry allowed.'
344 call mem_setptr(obs6_fnames,
'OBS6_FILENAME', this%input_mempath)
346 found%obs6_filename = .true.
347 this%obs%inputFilename = obs6_fnames(1)
348 this%obs%active = .true.
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()
359 if (this%iout > 0)
then
360 call this%log_options(found)
363 end subroutine source_options
367 subroutine log_options(this, found)
369 class(SwfDfwType) :: this
370 type(SwfDfwParamFoundType),
intent(in) :: found
372 write (this%iout,
'(1x,a)')
'Setting DFW Options'
374 if (found%lengthconv)
then
375 write (this%iout,
'(4x,a, G0)')
'Mannings length conversion value &
376 &specified as ', this%lengthconv
379 if (found%timeconv)
then
380 write (this%iout,
'(4x,a, G0)')
'Mannings time conversion value &
381 &specified as ', this%timeconv
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
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.'
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.'
400 if (found%obs6_filename)
then
401 write (this%iout,
'(4x,a)')
'Observation package is active.'
404 if (found%isavvelocity) &
405 write (this%iout,
'(4x,a)')
'Velocity will be calculated at cell &
406 ¢ers and written to DATA-VCOMP in budget &
407 &file when requested.'
409 if (found%iswrcond)
then
410 write (this%iout,
'(4x,a, G0)')
'Conductance will be calculated using &
411 &the SWR development option.'
414 write (this%iout,
'(1x,a,/)')
'End Setting DFW Options'
416 end subroutine log_options
420 subroutine source_griddata(this)
428 class(SwfDfwType) :: this
430 character(len=LENMEMPATH) :: idmMemoryPath
431 type(SwfDfwParamFoundType) :: found
432 integer(I4B),
dimension(:),
pointer,
contiguous :: map
439 if (this%dis%nodes < this%dis%nodesuser) map => this%dis%nodeuser
443 idmmemorypath, map, found%manningsn)
444 call mem_set_value(this%idcxs,
'IDCXS', idmmemorypath, map, found%idcxs)
447 if (.not. found%manningsn)
then
448 write (errmsg,
'(a)')
'Error in GRIDDATA block: MANNINGSN not found.'
453 call store_error_filename(this%input_fname)
457 if (this%iout > 0)
then
458 call this%log_griddata(found)
461 end subroutine source_griddata
465 subroutine log_griddata(this, found)
467 class(SwfDfwType) :: this
468 type(SwfDfwParamFoundType),
intent(in) :: found
470 write (this%iout,
'(1x,a)')
'Setting DFW Griddata'
472 if (found%manningsn)
then
473 write (this%iout,
'(4x,a)')
'MANNINGSN set from input file'
476 if (found%idcxs)
then
477 write (this%iout,
'(4x,a)')
'IDCXS set from input file'
480 call this%write_cxs_tables()
482 write (this%iout,
'(1x,a,/)')
'End Setting DFW Griddata'
484 end subroutine log_griddata
486 subroutine write_cxs_tables(this)
489 class(SwfDfwType) :: this
502 end subroutine write_cxs_tables
506 subroutine dfw_ar(this, ibound, hnew)
509 class(SwfDfwType) :: this
510 integer(I4B),
dimension(:),
pointer,
contiguous :: ibound
511 real(DP),
dimension(:),
pointer,
contiguous,
intent(inout) :: hnew
516 this%ibound => ibound
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', &
526 do n = 1, this%dis%nodes
527 this%vcomp(:, n) =
dzero
533 call this%obs%obs_ar()
535 end subroutine dfw_ar
539 subroutine dfw_rp(this)
542 class(SwfDfwType) :: this
545 call this%dfw_rp_obs()
547 end subroutine dfw_rp
551 subroutine dfw_ad(this, irestore)
552 class(SwfDfwType) :: this
553 integer(I4B),
intent(in) :: irestore
556 call this%obs%obs_ad()
558 end subroutine dfw_ad
566 subroutine dfw_fc(this, kiter, matrix_sln, idxglo, rhs, stage, stage_old)
569 class(SwfDfwType) :: this
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
579 if (this%is2d == 1)
then
580 call this%calc_dhds()
584 call this%dfw_qnm_fc_nr(kiter, matrix_sln, idxglo, rhs, stage, stage_old)
586 end subroutine dfw_fc
591 subroutine dfw_qnm_fc_nr(this, kiter, matrix_sln, idxglo, rhs, stage, stage_old)
595 class(SwfDfwType) :: this
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
603 integer(I4B) :: n, m, ii, idiag
610 do n = 1, this%dis%nodes
613 idiag = this%dis%con%ia(n)
616 do ii = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
619 if (this%dis%con%mask(ii) == 0) cycle
622 m = this%dis%con%ja(ii)
625 qnm = this%qcalc(n, m, stage(n), stage(m), ii)
626 rhs(n) = rhs(n) - qnm
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)
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)
645 end subroutine dfw_qnm_fc_nr
649 subroutine dfw_fn(this, kiter, matrix_sln, idxglo, rhs, stage)
651 class(SwfDfwType) :: this
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
663 end subroutine dfw_fn
667 function qcalc(this, n, m, stage_n, stage_m, ipos)
result(qnm)
669 class(SwfDfwType) :: this
670 integer(I4B),
intent(in) :: n
671 integer(I4B),
intent(in) :: m
672 real(DP),
intent(in) :: stage_n
673 real(DP),
intent(in) :: stage_m
674 integer(I4B),
intent(in) :: ipos
676 integer(I4B) :: isympos
683 isympos = this%dis%con%jas(ipos)
685 cl1 = this%dis%con%cl1(isympos)
686 cl2 = this%dis%con%cl2(isympos)
688 cl1 = this%dis%con%cl2(isympos)
689 cl2 = this%dis%con%cl1(isympos)
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)
700 qnm = cond * (stage_m - stage_n)
709 function get_cond(this, n, m, ipos, stage_n, stage_m, cln, clm)
result(cond)
713 class(SwfDfwType) :: this
714 integer(I4B),
intent(in) :: n
715 integer(I4B),
intent(in) :: m
716 integer(I4B),
intent(in) :: ipos
717 real(DP),
intent(in) :: stage_n
718 real(DP),
intent(in) :: stage_m
719 real(DP),
intent(in) :: cln
720 real(DP),
intent(in) :: clm
728 real(DP) :: range = 1.d-6
730 real(DP) :: smooth_factor
731 real(DP) :: length_nm
739 length_nm = cln + clm
741 if (length_nm >
dprec)
then
744 depth_n = stage_n - this%dis%bot(n)
745 depth_m = stage_m - this%dis%bot(m)
748 if (this%is2d == 0)
then
749 dhds_n = abs(stage_m - stage_n) / (cln + clm)
752 dhds_n = this%grad_dhds_mag(n)
753 dhds_m = this%grad_dhds_mag(m)
757 if (this%icentral == 0)
then
759 if (stage_n > stage_m)
then
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
774 call this%dis%get_flow_width(n, m, ipos, width_n, width_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)
784 if ((cn + cm) >
dprec)
then
785 cond = cn * cm / (cn + cm)
792 end function get_cond
798 function get_cond_n(this, n, depth, dx, width, dhds)
result(c)
801 class(SwfDfwType) :: this
802 integer(I4B),
intent(in) :: n
803 real(DP),
intent(in) :: depth
804 real(DP),
intent(in) :: dx
805 real(DP),
intent(in) :: width
806 real(DP),
intent(in) :: dhds
812 real(DP) :: conveyance
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
823 c = this%unitconv * conveyance / dx / dhds_sqr
825 end function get_cond_n
833 function get_cond_swr(this, n, m, ipos, stage_n, stage_m, cln, clm)
result(cond)
837 class(SwfDfwType) :: this
838 integer(I4B),
intent(in) :: n
839 integer(I4B),
intent(in) :: m
840 integer(I4B),
intent(in) :: ipos
841 real(DP),
intent(in) :: stage_n
842 real(DP),
intent(in) :: stage_m
843 real(DP),
intent(in) :: cln
844 real(DP),
intent(in) :: clm
854 real(DP) :: range = 1.d-6
856 real(DP) :: smooth_factor
857 real(DP) :: length_nm
861 real(DP) :: area_n, area_m, area_avg
862 real(DP) :: rhn, rhm, rhavg
870 length_nm = cln + clm
872 if (length_nm >
dprec)
then
875 depth_n = stage_n - this%dis%bot(n)
876 depth_m = stage_m - this%dis%bot(m)
879 if (this%icentral == 0)
then
881 if (stage_n > stage_m)
then
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
896 call this%dis%get_flow_width(n, m, ipos, width_n, width_m)
899 weight_n = clm / length_nm
900 weight_m =
done - weight_n
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
908 if (this%is2d == 0)
then
909 rhn = this%cxs%get_hydraulic_radius(this%idcxs(n), width_n, &
911 rhm = this%cxs%get_hydraulic_radius(this%idcxs(m), width_m, &
913 rhavg = weight_n * rhn + weight_m * rhm
915 rhavg = area_avg / width_n
920 if (this%is2d == 0)
then
921 dhds_nm = abs(stage_m - stage_n) / (length_nm)
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
927 dhds_sqr = dhds_nm**
dhalf
928 if (dhds_sqr <
dem10)
then
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
944 cond = this%unitconv * rinv_avg * area_avg * rhavg / dhds_sqr / length_nm
948 end function get_cond_swr
956 function get_flow_area_nm(this, n, m, stage_n, stage_m, cln, clm, &
957 ipos)
result(area_avg)
961 class(SwfDfwType) :: this
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
978 real(DP) :: length_nm
979 real(DP) :: range = 1.d-6
981 real(DP) :: smooth_factor
986 depth_n = stage_n - this%dis%bot(n)
987 depth_m = stage_m - this%dis%bot(m)
990 if (this%icentral == 0)
then
992 if (stage_n > stage_m)
then
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
1007 call this%dis%get_flow_width(n, m, ipos, width_n, width_m)
1010 length_nm = cln + clm
1011 weight_n = clm / length_nm
1012 weight_m =
done - weight_n
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
1019 end function get_flow_area_nm
1027 subroutine calc_dhds(this)
1031 class(SwfDfwType) :: this
1035 integer(I4B) :: ipos
1036 integer(I4B) :: isympos
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)
1048 cl1 = this%dis%con%cl1(isympos)
1049 cl2 = this%dis%con%cl2(isympos)
1051 cl1 = this%dis%con%cl2(isympos)
1052 cl2 = this%dis%con%cl1(isympos)
1057 if (cl1 + cl2 >
dprec)
then
1058 this%dhdsja(isympos) = (this%hnew(m) - this%hnew(n)) / (cl1 + cl2)
1060 this%dhdsja(isympos) =
dzero
1070 end subroutine calc_dhds
1075 subroutine dfw_nur(this, neqmod, x, xtemp, dx, inewtonur, dxmax, locmax)
1077 class(SwfDfwType) :: this
1078 integer(I4B),
intent(in) :: neqmod
1079 real(DP),
dimension(neqmod),
intent(inout) :: x
1080 real(DP),
dimension(neqmod),
intent(in) :: xtemp
1081 real(DP),
dimension(neqmod),
intent(inout) :: dx
1082 integer(I4B),
intent(inout) :: inewtonur
1083 real(DP),
intent(inout) :: dxmax
1084 integer(I4B),
intent(inout) :: locmax
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)
1098 if (x(n) < botm)
then
1102 if (abs(dxx) > abs(dxmax))
then
1112 end subroutine dfw_nur
1116 subroutine dfw_cq(this, stage, stage_old, flowja)
1118 class(SwfDfwType) :: this
1119 real(DP),
intent(inout),
dimension(:) :: stage
1120 real(DP),
intent(inout),
dimension(:) :: stage_old
1121 real(DP),
intent(inout),
dimension(:) :: flowja
1123 integer(I4B) :: n, ipos, m
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)
1130 qnm = this%qcalc(n, m, stage(n), stage(m), ipos)
1132 flowja(this%dis%con%isym(ipos)) = -qnm
1136 end subroutine dfw_cq
1140 subroutine dfw_bd(this, isuppress_output, model_budget)
1144 class(SwfDfwType) :: this
1145 integer(I4B),
intent(in) :: isuppress_output
1146 type(BudgetType),
intent(inout) :: model_budget
1151 end subroutine dfw_bd
1155 subroutine dfw_save_model_flows(this, flowja, icbcfl, icbcun)
1157 class(SwfDfwType) :: this
1158 real(DP),
dimension(:),
intent(in) :: flowja
1159 integer(I4B),
intent(in) :: icbcfl
1160 integer(I4B),
intent(in) :: icbcun
1162 integer(I4B) :: ibinun
1165 if (this%ipakcb < 0)
then
1167 elseif (this%ipakcb == 0)
then
1170 ibinun = this%ipakcb
1172 if (icbcfl == 0) ibinun = 0
1175 if (ibinun /= 0)
then
1177 call this%dis%record_connection_array(flowja, ibinun, this%iout)
1181 if (this%isavvelocity /= 0)
then
1182 if (ibinun /= 0)
call this%sav_velocity(ibinun)
1185 end subroutine dfw_save_model_flows
1189 subroutine dfw_print_model_flows(this, ibudfl, flowja)
1194 class(SwfDfwType) :: this
1195 integer(I4B),
intent(in) :: ibudfl
1196 real(DP),
intent(inout),
dimension(:) :: flowja
1198 character(len=LENBIGLINE) :: line
1199 character(len=30) :: tempstr
1200 integer(I4B) :: n, ipos, m
1203 character(len=*),
parameter :: fmtiprflow = &
1204 &
"(/,4x,'CALCULATED INTERCELL FLOW FOR PERIOD ', i0, ' STEP ', i0)"
1207 if (ibudfl /= 0 .and. this%iprflow > 0)
then
1208 write (this%iout, fmtiprflow)
kper,
kstp
1209 do n = 1, this%dis%nodes
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)
1218 write (tempstr,
'(1pg15.6)') qnm
1219 line = trim(line)//
' '//trim(adjustl(tempstr))
1221 write (this%iout,
'(a)') trim(line)
1225 end subroutine dfw_print_model_flows
1229 subroutine dfw_da(this)
1235 class(SwfDfwType) :: this
1249 if (this%is2d == 1)
then
1268 call this%obs%obs_da()
1269 deallocate (this%obs)
1274 call this%NumericalPackageType%da()
1279 end subroutine dfw_da
1285 subroutine calc_velocity(this, flowja)
1288 class(SwfDfwType) :: this
1289 real(DP),
intent(in),
dimension(:) :: flowja
1293 integer(I4B) :: ipos
1294 integer(I4B) :: isympos
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.
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.)
1343 do n = 1, this%dis%nodes
1346 ic = this%dis%con%ia(n + 1) - this%dis%con%ia(n) - 1
1349 do m = 1, this%nedges
1350 if (this%nodedge(m) == n)
then
1356 if (ic > nc) nc = ic
1373 do n = 1, this%dis%nodes
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)
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)
1396 cl1 = this%dis%con%cl2(isympos)
1397 cl2 = this%dis%con%cl1(isympos)
1399 ooclsum =
done / (cl1 + cl2)
1402 di(ic) = dltot * cl1 * ooclsum
1403 area = this%get_flow_area_nm(n, m, this%hnew(n), this%hnew(m), &
1405 if (area >
dzero)
then
1406 vi(ic) = flowja(ipos) / area
1415 do m = 1, this%nedges
1416 if (this%nodedge(m) == n)
then
1419 ihc = this%ihcedge(m)
1420 area = this%propsedge(2, m)
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
1445 dsumz = dsumz + diz(iz)
1447 denom = (ncz -
done)
1449 dsumz = dsumz +
dem10 * dsumz
1451 if (dsumz >
dzero) wiz(iz) =
done - diz(iz) / dsumz
1453 wiz(iz) = wiz(iz) / denom
1461 vz = vz + wiz(iz) * viz(iz)
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)
1480 dsumx = dsumx +
dem10 * dsumx
1481 dsumy = dsumy +
dem10 * dsumy
1483 wix(ic) = (dsumx - wix(ic)) * abs(nix(ic))
1484 wiy(ic) = (dsumy - wiy(ic)) * abs(niy(ic))
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))
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)
1517 vx = vx + (bix(ic) - axy * biy(ic)) * vi(ic)
1518 vy = vy + (biy(ic) - ayx * bix(ic)) * vi(ic)
1520 denom =
done - axy * ayx
1521 if (denom /=
dzero)
then
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)
1544 end subroutine calc_velocity
1551 subroutine increase_edge_count(this, nedges)
1553 class(SwfDfwType) :: this
1554 integer(I4B),
intent(in) :: nedges
1556 this%nedges = this%nedges + nedges
1558 end subroutine increase_edge_count
1564 subroutine set_edge_properties(this, nodedge, ihcedge, q, area, nx, ny, &
1567 class(SwfDfwType) :: this
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
1576 integer(I4B) :: lastedge
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
1590 if (this%lastedge == this%nedges) this%lastedge = 0
1592 end subroutine set_edge_properties
1598 subroutine sav_velocity(this, ibinun)
1600 class(SwfDfwType) :: this
1601 integer(I4B),
intent(in) :: ibinun
1603 character(len=16) :: text
1604 character(len=16),
dimension(3) :: auxtxt
1606 integer(I4B) :: naux
1609 text =
' DATA-VCOMP'
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)
1618 do n = 1, this%dis%nodes
1619 call this%dis%record_mf6_list_entry(ibinun, n, n,
dzero, naux, &
1623 end subroutine sav_velocity
1628 subroutine dfw_df_obs(this)
1630 class(SwfDfwType) :: this
1632 integer(I4B) :: indx
1636 call this%obs%StoreObsType(
'ext-outflow', .true., indx)
1637 this%obs%obsData(indx)%ProcessIdPtr => dfwobsidprocessor
1639 end subroutine dfw_df_obs
1641 subroutine dfwobsidprocessor(obsrv, dis, inunitobs, iout)
1643 type(ObserveType),
intent(inout) :: obsrv
1644 class(DisBaseType),
intent(in) :: dis
1645 integer(I4B),
intent(in) :: inunitobs
1646 integer(I4B),
intent(in) :: iout
1649 character(len=LINELENGTH) :: string
1652 string = obsrv%IDstring
1656 obsrv%NodeNumber = n
1658 errmsg =
'Error reading data from ID string'
1663 end subroutine dfwobsidprocessor
1668 subroutine dfw_bd_obs(this)
1670 class(SwfDfwType) :: this
1676 character(len=100) :: msg
1677 type(ObserveType),
pointer :: obsrv => null()
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)
1687 select case (obsrv%ObsTypeId)
1689 msg =
'Unrecognized observation type: '//trim(obsrv%ObsTypeId)
1692 call this%obs%SaveOneSimval(obsrv, v)
1702 end subroutine dfw_bd_obs
1707 subroutine dfw_rp_obs(this)
1711 class(SwfDfwType),
intent(inout) :: this
1716 class(ObserveType),
pointer :: obsrv => null()
1722 do i = 1, this%obs%npakobs
1723 obsrv => this%obs%pakobs(i)%obsrv
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,
')'
1734 if (obsrv%indxbnds_count == 0)
then
1735 call obsrv%AddObsIndex(nn1)
1737 errmsg =
'Programming error in dfw_rp_obs'
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,
')'
1762 end subroutine dfw_rp_obs
1764 end module swfdfwmodule
This module contains the BudgetModule.
This module contains simulation constants.
integer(i4b), parameter linelength
maximum length of a standard line
real(dp), parameter dtwothirds
real constant 2/3
real(dp), parameter dp9
real constant 9/10
real(dp), parameter dem10
real constant 1e-10
real(dp), parameter donethird
real constant 1/3
integer(i4b), parameter lenbigline
maximum length of a big line
real(dp), parameter dhalf
real constant 1/2
real(dp), parameter dzero
real constant zero
real(dp), parameter dprec
real constant machine precision
real(dp), parameter dtwo
real constant 2
integer(i4b), parameter lenmempath
maximum length of the memory path
real(dp), parameter done
real constant 1
This module defines variable data types.
real(dp) function, public get_perturbation(x)
Calculate a numerical perturbation given the value of x.
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.
This module contains the derived type ObsType.
subroutine, public obs_cr(obs, inobs)
@ brief Create a new ObsType object
This module contains simulation methods.
subroutine, public store_error(msg, terminate)
Store an error message.
integer(i4b) function, public count_errors()
Return number of errors.
subroutine, public store_error_filename(filename, terminate)
Store the erroring file name.
subroutine, public store_error_unit(iunit, terminate)
Store the file unit number.
This module contains simulation variables.
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
integer(i4b), pointer, public kper
current stress period number
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.
This class is used to store a single deferred-length character string. It was designed to work in an ...