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
167 character(len=10) :: distype =
''
173 call this%dis%get_dis_type(distype)
174 if (distype ==
"DIS2D")
then
177 if (distype ==
"DISV2D")
then
186 call this%allocate_arrays()
193 end subroutine dfw_df
201 subroutine allocate_scalars(this)
204 class(SwfDfwtype) :: this
207 call this%NumericalPackageType%allocate_scalars()
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)
226 this%lengthconv =
done
229 this%icalcvelocity = 0
230 this%isavvelocity = 0
234 end subroutine allocate_scalars
238 subroutine allocate_arrays(this)
240 class(SwfDfwType) :: this
246 'MANNINGSN', this%memoryPath)
248 'IDCXS', this%memoryPath)
250 'ICELLTYPE', this%memoryPath)
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)
258 call mem_allocate(this%vcomp, 3, 0,
'VCOMP', this%memoryPath)
259 call mem_allocate(this%vmag, 0,
'VMAG', this%memoryPath)
261 do n = 1, this%dis%nodes
262 this%manningsn(n) =
dzero
264 this%icelltype(n) = 1
268 if (this%is2d == 1)
then
270 'GRAD_DHDS_MAG', this%memoryPath)
272 'DHDSJA', this%memoryPath)
273 do n = 1, this%dis%nodes
274 this%grad_dhds_mag(n) =
dzero
276 do n = 1, this%dis%njas
277 this%dhdsja(n) =
dzero
281 end subroutine allocate_arrays
285 subroutine dfw_load(this)
287 class(SwfDfwType) :: this
290 call this%source_options()
291 call this%source_griddata()
293 end subroutine dfw_load
297 subroutine source_options(this)
305 class(SwfDfwType) :: this
307 integer(I4B) :: isize
308 type(SwfDfwParamFoundType) :: found
309 type(CharacterStringType),
dimension(:),
pointer, &
310 contiguous :: obs6_fnames
314 this%input_mempath, found%icentral)
316 this%input_mempath, found%iswrcond)
318 this%input_mempath, found%lengthconv)
320 this%input_mempath, found%timeconv)
322 this%input_mempath, found%iprflow)
324 this%input_mempath, found%ipakcb)
326 this%input_mempath, found%isavvelocity)
329 if (found%icentral) this%icentral = 1
330 if (found%ipakcb) this%ipakcb = -1
333 this%unitconv = this%lengthconv**
donethird
334 this%unitconv = this%unitconv * this%timeconv
337 if (found%isavvelocity) this%icalcvelocity = this%isavvelocity
340 call get_isize(
'OBS6_FILENAME', this%input_mempath, isize)
344 errmsg =
'Multiple OBS6 keywords detected in OPTIONS block.'// &
345 ' Only one OBS6 entry allowed.'
350 call mem_setptr(obs6_fnames,
'OBS6_FILENAME', this%input_mempath)
352 found%obs6_filename = .true.
353 this%obs%inputFilename = obs6_fnames(1)
354 this%obs%active = .true.
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()
363 if (this%iout > 0)
then
364 call this%log_options(found)
367 end subroutine source_options
371 subroutine log_options(this, found)
373 class(SwfDfwType) :: this
374 type(SwfDfwParamFoundType),
intent(in) :: found
376 write (this%iout,
'(1x,a)')
'Setting DFW Options'
378 if (found%lengthconv)
then
379 write (this%iout,
'(4x,a, G0)')
'Mannings length conversion value &
380 &specified as ', this%lengthconv
383 if (found%timeconv)
then
384 write (this%iout,
'(4x,a, G0)')
'Mannings time conversion value &
385 &specified as ', this%timeconv
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
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.'
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.'
404 if (found%obs6_filename)
then
405 write (this%iout,
'(4x,a)')
'Observation package is active.'
408 if (found%isavvelocity) &
409 write (this%iout,
'(4x,a)')
'Velocity will be calculated at cell &
410 ¢ers and written to DATA-VCOMP in budget &
411 &file when requested.'
413 if (found%iswrcond)
then
414 write (this%iout,
'(4x,a, G0)')
'Conductance will be calculated using &
415 &the SWR development option.'
418 write (this%iout,
'(1x,a,/)')
'End Setting DFW Options'
420 end subroutine log_options
424 subroutine source_griddata(this)
432 class(SwfDfwType) :: this
434 character(len=LENMEMPATH) :: idmMemoryPath
435 type(SwfDfwParamFoundType) :: found
436 integer(I4B),
dimension(:),
pointer,
contiguous :: map
443 if (this%dis%nodes < this%dis%nodesuser) map => this%dis%nodeuser
447 idmmemorypath, map, found%manningsn)
448 call mem_set_value(this%idcxs,
'IDCXS', idmmemorypath, map, found%idcxs)
451 if (.not. found%manningsn)
then
452 write (errmsg,
'(a)')
'Error in GRIDDATA block: MANNINGSN not found.'
457 call store_error_filename(this%input_fname)
461 if (this%iout > 0)
then
462 call this%log_griddata(found)
465 end subroutine source_griddata
469 subroutine log_griddata(this, found)
471 class(SwfDfwType) :: this
472 type(SwfDfwParamFoundType),
intent(in) :: found
474 write (this%iout,
'(1x,a)')
'Setting DFW Griddata'
476 if (found%manningsn)
then
477 write (this%iout,
'(4x,a)')
'MANNINGSN set from input file'
480 if (found%idcxs)
then
481 write (this%iout,
'(4x,a)')
'IDCXS set from input file'
484 call this%write_cxs_tables()
486 write (this%iout,
'(1x,a,/)')
'End Setting DFW Griddata'
488 end subroutine log_griddata
490 subroutine write_cxs_tables(this)
493 class(SwfDfwType) :: this
506 end subroutine write_cxs_tables
510 subroutine dfw_ar(this, ibound, hnew)
513 class(SwfDfwType) :: this
514 integer(I4B),
dimension(:),
pointer,
contiguous :: ibound
515 real(DP),
dimension(:),
pointer,
contiguous,
intent(inout) :: hnew
520 this%ibound => ibound
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', &
530 do n = 1, this%dis%nodes
531 this%vcomp(:, n) =
dzero
537 call this%obs%obs_ar()
539 end subroutine dfw_ar
543 subroutine dfw_rp(this)
546 class(SwfDfwType) :: this
549 call this%dfw_rp_obs()
551 end subroutine dfw_rp
555 subroutine dfw_ad(this, irestore)
556 class(SwfDfwType) :: this
557 integer(I4B),
intent(in) :: irestore
560 call this%obs%obs_ad()
562 end subroutine dfw_ad
570 subroutine dfw_fc(this, kiter, matrix_sln, idxglo, rhs, stage, stage_old)
573 class(SwfDfwType) :: this
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
583 if (this%is2d == 1)
then
584 call this%calc_dhds()
588 call this%dfw_qnm_fc_nr(kiter, matrix_sln, idxglo, rhs, stage, stage_old)
590 end subroutine dfw_fc
595 subroutine dfw_qnm_fc_nr(this, kiter, matrix_sln, idxglo, rhs, stage, stage_old)
599 class(SwfDfwType) :: this
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
607 integer(I4B) :: n, m, ii, idiag
614 do n = 1, this%dis%nodes
617 idiag = this%dis%con%ia(n)
620 do ii = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
623 if (this%dis%con%mask(ii) == 0) cycle
626 m = this%dis%con%ja(ii)
629 qnm = this%qcalc(n, m, stage(n), stage(m), ii)
630 rhs(n) = rhs(n) - qnm
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)
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)
649 end subroutine dfw_qnm_fc_nr
653 subroutine dfw_fn(this, kiter, matrix_sln, idxglo, rhs, stage)
655 class(SwfDfwType) :: this
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
667 end subroutine dfw_fn
671 function qcalc(this, n, m, stage_n, stage_m, ipos)
result(qnm)
673 class(SwfDfwType) :: this
674 integer(I4B),
intent(in) :: n
675 integer(I4B),
intent(in) :: m
676 real(DP),
intent(in) :: stage_n
677 real(DP),
intent(in) :: stage_m
678 integer(I4B),
intent(in) :: ipos
680 integer(I4B) :: isympos
687 isympos = this%dis%con%jas(ipos)
689 cl1 = this%dis%con%cl1(isympos)
690 cl2 = this%dis%con%cl2(isympos)
692 cl1 = this%dis%con%cl2(isympos)
693 cl2 = this%dis%con%cl1(isympos)
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)
704 qnm = cond * (stage_m - stage_n)
713 function get_cond(this, n, m, ipos, stage_n, stage_m, cln, clm)
result(cond)
717 class(SwfDfwType) :: this
718 integer(I4B),
intent(in) :: n
719 integer(I4B),
intent(in) :: m
720 integer(I4B),
intent(in) :: ipos
721 real(DP),
intent(in) :: stage_n
722 real(DP),
intent(in) :: stage_m
723 real(DP),
intent(in) :: cln
724 real(DP),
intent(in) :: clm
732 real(DP) :: range = 1.d-6
734 real(DP) :: smooth_factor
735 real(DP) :: length_nm
743 length_nm = cln + clm
745 if (length_nm >
dprec)
then
748 depth_n = stage_n - this%dis%bot(n)
749 depth_m = stage_m - this%dis%bot(m)
752 if (this%is2d == 0)
then
753 dhds_n = abs(stage_m - stage_n) / (cln + clm)
756 dhds_n = this%grad_dhds_mag(n)
757 dhds_m = this%grad_dhds_mag(m)
761 if (this%icentral == 0)
then
763 if (stage_n > stage_m)
then
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
778 call this%dis%get_flow_width(n, m, ipos, width_n, width_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)
788 if ((cn + cm) >
dprec)
then
789 cond = cn * cm / (cn + cm)
796 end function get_cond
802 function get_cond_n(this, n, depth, dx, width, dhds)
result(c)
805 class(SwfDfwType) :: this
806 integer(I4B),
intent(in) :: n
807 real(DP),
intent(in) :: depth
808 real(DP),
intent(in) :: dx
809 real(DP),
intent(in) :: width
810 real(DP),
intent(in) :: dhds
816 real(DP) :: conveyance
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
827 c = this%unitconv * conveyance / dx / dhds_sqr
829 end function get_cond_n
837 function get_cond_swr(this, n, m, ipos, stage_n, stage_m, cln, clm)
result(cond)
841 class(SwfDfwType) :: this
842 integer(I4B),
intent(in) :: n
843 integer(I4B),
intent(in) :: m
844 integer(I4B),
intent(in) :: ipos
845 real(DP),
intent(in) :: stage_n
846 real(DP),
intent(in) :: stage_m
847 real(DP),
intent(in) :: cln
848 real(DP),
intent(in) :: clm
858 real(DP) :: range = 1.d-6
860 real(DP) :: smooth_factor
861 real(DP) :: length_nm
865 real(DP) :: area_n, area_m, area_avg
866 real(DP) :: rhn, rhm, rhavg
874 length_nm = cln + clm
876 if (length_nm >
dprec)
then
879 depth_n = stage_n - this%dis%bot(n)
880 depth_m = stage_m - this%dis%bot(m)
883 if (this%icentral == 0)
then
885 if (stage_n > stage_m)
then
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
900 call this%dis%get_flow_width(n, m, ipos, width_n, width_m)
903 weight_n = clm / length_nm
904 weight_m =
done - weight_n
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
912 if (this%is2d == 0)
then
913 rhn = this%cxs%get_hydraulic_radius(this%idcxs(n), width_n, &
915 rhm = this%cxs%get_hydraulic_radius(this%idcxs(m), width_m, &
917 rhavg = weight_n * rhn + weight_m * rhm
919 rhavg = area_avg / width_n
924 if (this%is2d == 0)
then
925 dhds_nm = abs(stage_m - stage_n) / (length_nm)
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
931 dhds_sqr = dhds_nm**
dhalf
932 if (dhds_sqr <
dem10)
then
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
948 cond = this%unitconv * rinv_avg * area_avg * rhavg / dhds_sqr / length_nm
952 end function get_cond_swr
960 function get_flow_area_nm(this, n, m, stage_n, stage_m, cln, clm, &
961 ipos)
result(area_avg)
965 class(SwfDfwType) :: this
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
982 real(DP) :: length_nm
983 real(DP) :: range = 1.d-6
985 real(DP) :: smooth_factor
990 depth_n = stage_n - this%dis%bot(n)
991 depth_m = stage_m - this%dis%bot(m)
994 if (this%icentral == 0)
then
996 if (stage_n > stage_m)
then
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
1011 call this%dis%get_flow_width(n, m, ipos, width_n, width_m)
1014 length_nm = cln + clm
1015 weight_n = clm / length_nm
1016 weight_m =
done - weight_n
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
1023 end function get_flow_area_nm
1031 subroutine calc_dhds(this)
1035 class(SwfDfwType) :: this
1039 integer(I4B) :: ipos
1040 integer(I4B) :: isympos
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)
1052 cl1 = this%dis%con%cl1(isympos)
1053 cl2 = this%dis%con%cl2(isympos)
1055 cl1 = this%dis%con%cl2(isympos)
1056 cl2 = this%dis%con%cl1(isympos)
1061 if (cl1 + cl2 >
dprec)
then
1062 this%dhdsja(isympos) = (this%hnew(m) - this%hnew(n)) / (cl1 + cl2)
1064 this%dhdsja(isympos) =
dzero
1074 end subroutine calc_dhds
1079 subroutine dfw_nur(this, neqmod, x, xtemp, dx, inewtonur, dxmax, locmax)
1081 class(SwfDfwType) :: this
1082 integer(I4B),
intent(in) :: neqmod
1083 real(DP),
dimension(neqmod),
intent(inout) :: x
1084 real(DP),
dimension(neqmod),
intent(in) :: xtemp
1085 real(DP),
dimension(neqmod),
intent(inout) :: dx
1086 integer(I4B),
intent(inout) :: inewtonur
1087 real(DP),
intent(inout) :: dxmax
1088 integer(I4B),
intent(inout) :: locmax
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)
1102 if (x(n) < botm)
then
1106 if (abs(dxx) > abs(dxmax))
then
1116 end subroutine dfw_nur
1120 subroutine dfw_cq(this, stage, stage_old, flowja)
1122 class(SwfDfwType) :: this
1123 real(DP),
intent(inout),
dimension(:) :: stage
1124 real(DP),
intent(inout),
dimension(:) :: stage_old
1125 real(DP),
intent(inout),
dimension(:) :: flowja
1127 integer(I4B) :: n, ipos, m
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)
1134 qnm = this%qcalc(n, m, stage(n), stage(m), ipos)
1136 flowja(this%dis%con%isym(ipos)) = -qnm
1140 end subroutine dfw_cq
1144 subroutine dfw_bd(this, isuppress_output, model_budget)
1148 class(SwfDfwType) :: this
1149 integer(I4B),
intent(in) :: isuppress_output
1150 type(BudgetType),
intent(inout) :: model_budget
1155 end subroutine dfw_bd
1159 subroutine dfw_save_model_flows(this, flowja, icbcfl, icbcun)
1161 class(SwfDfwType) :: this
1162 real(DP),
dimension(:),
intent(in) :: flowja
1163 integer(I4B),
intent(in) :: icbcfl
1164 integer(I4B),
intent(in) :: icbcun
1166 integer(I4B) :: ibinun
1169 if (this%ipakcb < 0)
then
1171 elseif (this%ipakcb == 0)
then
1174 ibinun = this%ipakcb
1176 if (icbcfl == 0) ibinun = 0
1179 if (ibinun /= 0)
then
1181 call this%dis%record_connection_array(flowja, ibinun, this%iout)
1185 if (this%isavvelocity /= 0)
then
1186 if (ibinun /= 0)
call this%sav_velocity(ibinun)
1189 end subroutine dfw_save_model_flows
1193 subroutine dfw_print_model_flows(this, ibudfl, flowja)
1198 class(SwfDfwType) :: this
1199 integer(I4B),
intent(in) :: ibudfl
1200 real(DP),
intent(inout),
dimension(:) :: flowja
1202 character(len=LENBIGLINE) :: line
1203 character(len=30) :: tempstr
1204 integer(I4B) :: n, ipos, m
1207 character(len=*),
parameter :: fmtiprflow = &
1208 &
"(/,4x,'CALCULATED INTERCELL FLOW FOR PERIOD ', i0, ' STEP ', i0)"
1211 if (ibudfl /= 0 .and. this%iprflow > 0)
then
1212 write (this%iout, fmtiprflow)
kper,
kstp
1213 do n = 1, this%dis%nodes
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)
1222 write (tempstr,
'(1pg15.6)') qnm
1223 line = trim(line)//
' '//trim(adjustl(tempstr))
1225 write (this%iout,
'(a)') trim(line)
1229 end subroutine dfw_print_model_flows
1233 subroutine dfw_da(this)
1239 class(SwfDfwType) :: this
1253 if (this%is2d == 1)
then
1272 call this%obs%obs_da()
1273 deallocate (this%obs)
1278 call this%NumericalPackageType%da()
1283 end subroutine dfw_da
1289 subroutine calc_velocity(this, flowja)
1292 class(SwfDfwType) :: this
1293 real(DP),
intent(in),
dimension(:) :: flowja
1297 integer(I4B) :: ipos
1298 integer(I4B) :: isympos
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.
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.)
1347 do n = 1, this%dis%nodes
1350 ic = this%dis%con%ia(n + 1) - this%dis%con%ia(n) - 1
1353 do m = 1, this%nedges
1354 if (this%nodedge(m) == n)
then
1360 if (ic > nc) nc = ic
1377 do n = 1, this%dis%nodes
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)
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)
1400 cl1 = this%dis%con%cl2(isympos)
1401 cl2 = this%dis%con%cl1(isympos)
1403 ooclsum =
done / (cl1 + cl2)
1406 di(ic) = dltot * cl1 * ooclsum
1407 area = this%get_flow_area_nm(n, m, this%hnew(n), this%hnew(m), &
1409 if (area >
dzero)
then
1410 vi(ic) = flowja(ipos) / area
1419 do m = 1, this%nedges
1420 if (this%nodedge(m) == n)
then
1423 ihc = this%ihcedge(m)
1424 area = this%propsedge(2, m)
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
1449 dsumz = dsumz + diz(iz)
1451 denom = (ncz -
done)
1453 dsumz = dsumz +
dem10 * dsumz
1455 if (dsumz >
dzero) wiz(iz) =
done - diz(iz) / dsumz
1457 wiz(iz) = wiz(iz) / denom
1465 vz = vz + wiz(iz) * viz(iz)
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)
1484 dsumx = dsumx +
dem10 * dsumx
1485 dsumy = dsumy +
dem10 * dsumy
1487 wix(ic) = (dsumx - wix(ic)) * abs(nix(ic))
1488 wiy(ic) = (dsumy - wiy(ic)) * abs(niy(ic))
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))
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)
1521 vx = vx + (bix(ic) - axy * biy(ic)) * vi(ic)
1522 vy = vy + (biy(ic) - ayx * bix(ic)) * vi(ic)
1524 denom =
done - axy * ayx
1525 if (denom /=
dzero)
then
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)
1548 end subroutine calc_velocity
1555 subroutine increase_edge_count(this, nedges)
1557 class(SwfDfwType) :: this
1558 integer(I4B),
intent(in) :: nedges
1560 this%nedges = this%nedges + nedges
1562 end subroutine increase_edge_count
1568 subroutine set_edge_properties(this, nodedge, ihcedge, q, area, nx, ny, &
1571 class(SwfDfwType) :: this
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
1580 integer(I4B) :: lastedge
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
1594 if (this%lastedge == this%nedges) this%lastedge = 0
1596 end subroutine set_edge_properties
1602 subroutine sav_velocity(this, ibinun)
1604 class(SwfDfwType) :: this
1605 integer(I4B),
intent(in) :: ibinun
1607 character(len=16) :: text
1608 character(len=16),
dimension(3) :: auxtxt
1610 integer(I4B) :: naux
1613 text =
' DATA-VCOMP'
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)
1622 do n = 1, this%dis%nodes
1623 call this%dis%record_mf6_list_entry(ibinun, n, n,
dzero, naux, &
1627 end subroutine sav_velocity
1632 subroutine dfw_df_obs(this)
1634 class(SwfDfwType) :: this
1636 integer(I4B) :: indx
1640 call this%obs%StoreObsType(
'ext-outflow', .true., indx)
1641 this%obs%obsData(indx)%ProcessIdPtr => dfwobsidprocessor
1643 end subroutine dfw_df_obs
1645 subroutine dfwobsidprocessor(obsrv, dis, inunitobs, iout)
1647 type(ObserveType),
intent(inout) :: obsrv
1648 class(DisBaseType),
intent(in) :: dis
1649 integer(I4B),
intent(in) :: inunitobs
1650 integer(I4B),
intent(in) :: iout
1653 character(len=LINELENGTH) :: string
1656 string = obsrv%IDstring
1660 obsrv%NodeNumber = n
1662 errmsg =
'Error reading data from ID string'
1667 end subroutine dfwobsidprocessor
1672 subroutine dfw_bd_obs(this)
1674 class(SwfDfwType) :: this
1680 character(len=100) :: msg
1681 type(ObserveType),
pointer :: obsrv => null()
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)
1691 select case (obsrv%ObsTypeId)
1693 msg =
'Unrecognized observation type: '//trim(obsrv%ObsTypeId)
1696 call this%obs%SaveOneSimval(obsrv, v)
1706 end subroutine dfw_bd_obs
1711 subroutine dfw_rp_obs(this)
1715 class(SwfDfwType),
intent(inout) :: this
1720 class(ObserveType),
pointer :: obsrv => null()
1726 do i = 1, this%obs%npakobs
1727 obsrv => this%obs%pakobs(i)%obsrv
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,
')'
1738 if (obsrv%indxbnds_count == 0)
then
1739 call obsrv%AddObsIndex(nn1)
1741 errmsg =
'Programming error in dfw_rp_obs'
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,
')'
1766 end subroutine dfw_rp_obs
1768 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 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 ...