-
Notifications
You must be signed in to change notification settings - Fork 318
New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
Ponded water/alternative solutions #191
Closed
Closed
Changes from all commits
Commits
File filter
Filter by extension
Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
There are no files selected for viewing
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
|
@@ -28,6 +28,10 @@ module SoilHydrologyMod | |
use LandunitType , only : lun | ||
use ColumnType , only : col | ||
use PatchType , only : patch | ||
|
||
use computeFluxMod ! template for the flux module | ||
use implicitEulerMod , only : implicitEuler ! implicit Euler solver for a single state | ||
|
||
! | ||
! !PUBLIC TYPES: | ||
implicit none | ||
|
@@ -57,6 +61,12 @@ module SoilHydrologyMod | |
!----------------------------------------------------------------------- | ||
real(r8), private :: baseflow_scalar = 1.e-2_r8 | ||
|
||
! !ALTERNATIVE NUMERICAL SOLUTIONS | ||
integer, parameter :: ixExplicitEuler=1001 ! constrained Explicit Euler solution with operator splitting (original) | ||
integer, parameter :: ixImplicitEuler=1002 ! implicit Euler solution with operator splitting | ||
integer, parameter :: ixAnalytical=1003 ! analytical solution with operator splitting | ||
integer :: ixSolution=ixExplicitEuler | ||
|
||
character(len=*), parameter, private :: sourcefile = & | ||
__FILE__ | ||
|
||
|
@@ -478,6 +488,9 @@ subroutine QflxH2osfcSurf(bounds, num_hydrologyc, filter_hydrologyc, & | |
real(r8) :: dtime ! land model time step (sec) | ||
real(r8) :: frac_infclust ! fraction of submerged area that is connected | ||
real(r8) :: k_wet ! linear reservoir coefficient for h2osfc | ||
real(r8) :: dRate ! effective linear reservoir coefficient for h2osfc (account for fraction of submerged area that is connected) | ||
real(r8) :: activePond ! ponded water above threshold = h2osfc - h2osfc_thresh | ||
real(r8) :: h2osfc1 ! h2osfc at the end of the time step, given qflx_h2osfc_surf only | ||
|
||
character(len=*), parameter :: subname = 'QflxH2osfcSurf' | ||
!----------------------------------------------------------------------- | ||
|
@@ -503,11 +516,36 @@ subroutine QflxH2osfcSurf(bounds, num_hydrologyc, filter_hydrologyc, & | |
|
||
! limit runoff to value of storage above S(pc) | ||
if(h2osfc(c) > h2osfc_thresh(c) .and. h2osfcflag/=0) then | ||
|
||
! spatially variable k_wet | ||
k_wet=1.0e-4_r8 * sin((rpi/180.) * topo_slope(c)) | ||
qflx_h2osfc_surf(c) = k_wet * frac_infclust * (h2osfc(c) - h2osfc_thresh(c)) | ||
dRate=k_wet*frac_infclust | ||
|
||
! ponding above threshold | ||
activePond = h2osfc(c) - h2osfc_thresh(c) | ||
|
||
! switch between different numerical solutions | ||
select case(ixSolution) | ||
|
||
! constrained Explicit Euler solution with operator splitting (original) | ||
case(ixExplicitEuler) | ||
qflx_h2osfc_surf(c)=min(dRate*activePond,(h2osfc(c) - h2osfc_thresh(c))/dtime) | ||
|
||
! implicit Euler solution with operator splitting | ||
case(ixImplicitEuler) | ||
qflx_h2osfc_surf(c) = dRate*activePond/(1._r8 + dRate*dtime) | ||
|
||
! analytical solution with operator splitting | ||
case(ixAnalytical) | ||
h2osfc1 = activePond*exp(-dRate*dtime) + h2osfc_thresh(c) | ||
qflx_h2osfc_surf(c) = (h2osfc1 - h2osfc(c))/dtime | ||
|
||
! check case identfied | ||
case default | ||
call endrun(subname // ':: the numerical solution must be specified!') | ||
|
||
end select ! (switch between different numerical solutions) | ||
|
||
qflx_h2osfc_surf(c)=min(qflx_h2osfc_surf(c),(h2osfc(c) - h2osfc_thresh(c))/dtime) | ||
else | ||
qflx_h2osfc_surf(c)= 0._r8 | ||
endif | ||
|
@@ -547,8 +585,16 @@ subroutine QflxH2osfcDrain(bounds, num_hydrologyc, filter_hydrologyc, & | |
logical , intent(inout) :: truncate_h2osfc_to_zero( bounds%begc: ) ! whether h2osfc should be truncated to 0 to correct for roundoff errors, in order to maintain bit-for-bit the same answers as the old code | ||
! | ||
! !LOCAL VARIABLES: | ||
integer :: fc, c | ||
real(r8) :: dtime ! land model time step (sec) | ||
procedure(fluxTemplate), pointer :: funcName ! function name | ||
integer :: fc, c | ||
real(r8) :: dtime ! land model time step (sec) | ||
real(r8) :: h2osfc1 ! h2osfc at the end of the time step, given qflx_h2osfc_drain only | ||
real(r8) :: drainMax ! maximum drainage rate of ponded water | ||
real(r8) :: const ! constant in analytical integral | ||
real(r8) :: uFunc ! analytical integral | ||
real(r8), parameter :: smoothScale=0.05_r8 ! smoothing scale | ||
integer(i4b) :: err ! error code | ||
character(len=128) :: message ! error message | ||
|
||
character(len=*), parameter :: subname = 'QflxH2osfcDrain' | ||
!----------------------------------------------------------------------- | ||
|
@@ -568,7 +614,32 @@ subroutine QflxH2osfcDrain(bounds, num_hydrologyc, filter_hydrologyc, & | |
qflx_h2osfc_drain(c) = h2osfc(c)/dtime | ||
truncate_h2osfc_to_zero(c) = .true. | ||
else | ||
qflx_h2osfc_drain(c)=min(frac_h2osfc(c)*qinmax(c),h2osfc(c)/dtime) | ||
|
||
! define maximum drainage | ||
! NOTE: check fraction | ||
drainMax = frac_h2osfc(c)*qinmax(c) | ||
|
||
! switch between different numerical solutions | ||
select case(ixSolution) | ||
|
||
! constrained Explicit Euler solution with operator splitting (original) | ||
case(ixExplicitEuler) | ||
qflx_h2osfc_drain(c)=min(drainMax,h2osfc(c)/dtime) | ||
|
||
! implicit Euler solution with operator splitting | ||
case(ixImplicitEuler) | ||
funcName=>drainPond | ||
call implicitEuler(funcName, dtime, h2osfc(c), h2osfc1, err, message) | ||
if(err/=0) call endrun(subname // ':: '//trim(message)) | ||
call drainPond(h2osfc1, qflx_h2osfc_drain(c)) | ||
|
||
! analytical solution with operator splitting | ||
case(ixAnalytical) | ||
const = 1._r8 - 1._r8/(1._r8 - exp(-h2osfc(c)/smoothScale)) | ||
uFunc = -1._r8/(const*exp(drainMax*dtime/smoothScale) - 1._r8) | ||
h2osfc1 = -alog(1._r8 - uFunc)*smoothScale | ||
qflx_h2osfc_drain(c)= (h2osfc1 - h2osfc(c))/dtime | ||
|
||
if(h2osfcflag==0) then | ||
! ensure no h2osfc | ||
qflx_h2osfc_drain(c)= max(0._r8,h2osfc(c)/dtime) | ||
|
@@ -577,6 +648,22 @@ subroutine QflxH2osfcDrain(bounds, num_hydrologyc, filter_hydrologyc, & | |
end if | ||
end do | ||
|
||
contains | ||
|
||
! model-specific flux routine | ||
subroutine drainPond(storage, flux, dfdx) | ||
! dummy variables | ||
real(r8), intent(in) :: storage ! storage | ||
real(r8), intent(out) :: flux ! drainage flux | ||
real(r8), intent(out),optional :: dfdx ! derivative | ||
! local variables | ||
real(r8) :: arg ! temporary argument | ||
! procedure starts here | ||
arg = exp(-storage/smoothScale) | ||
flux = -inputRate*(1._dp - arg) | ||
if(present(dfdx)) dfdx = -inputRate*arg/smoothScale | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Comment from @martynpclark Sep 7, 2017
|
||
end subroutine drainPond | ||
|
||
end subroutine QflxH2osfcDrain | ||
|
||
|
||
|
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,28 @@ | ||
module computeFluxMod | ||
|
||
! provides interface for the flux modules | ||
|
||
USE nrtype | ||
implicit none | ||
|
||
abstract interface | ||
subroutine fluxTemplate(x,f,dfdx) | ||
USE nrtype | ||
implicit none | ||
real(dp), intent(in) :: x ! state | ||
real(dp), intent(out) :: f ! f(x) | ||
real(dp), intent(out), optional :: dfdx ! derivative | ||
end subroutine fluxTemplate | ||
end interface | ||
|
||
contains | ||
|
||
subroutine getFlux(func, x, f, dfdx) | ||
procedure(fluxTemplate), pointer :: func | ||
real(dp), intent(in) :: x ! state | ||
real(dp), intent(out) :: f ! f(x) | ||
real(dp), intent(out),optional :: dfdx ! derivative | ||
call func(x,f,dfdx) | ||
end subroutine getFlux | ||
|
||
end module computeFluxMod |
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,49 @@ | ||
module implicitEulerMod | ||
|
||
USE nrtype | ||
USE computeFluxMod | ||
implicit none | ||
|
||
contains | ||
|
||
subroutine implicitEuler(funcName, dt, xInit, xNew, err, message) | ||
! dummy variables | ||
procedure(fluxTemplate), pointer :: funcName | ||
real(dp), intent(in) :: dt ! time step | ||
real(dp), intent(in) :: xInit ! initial state | ||
real(dp), intent(out) :: xNew ! updated state | ||
integer(i4b),intent(out) :: err ! error code | ||
character(*),intent(out) :: message ! error message | ||
! local variables | ||
integer(i4b) :: iter ! iteration index | ||
integer(i4b),parameter :: maxiter=20 ! maximum number of iterations | ||
real(dp), parameter :: xTol=1.e-8_dp ! convergence tolerance | ||
real(dp) :: flux ! flux | ||
real(dp) :: dfdx ! derivative | ||
real(dp) :: xRes ! residual error | ||
real(dp) :: delX ! state update | ||
! initialiuze routine | ||
err=0; message='implicitEuler/' | ||
|
||
! initialize | ||
xNew = xInit | ||
|
||
! iterate | ||
do iter=1,maxiter | ||
! get the flux and derivative | ||
call getFlux(funcName, xNew, flux, dfdx) | ||
! get the residual and the update | ||
xRes = xNew - (xInit + dt*flux) | ||
delX = -xRes/(1._dp - dt*dfdx) | ||
! update and check for convergence | ||
xNew = xNew + delX | ||
if(abs(delX) < xTol) exit | ||
! check for non-convergence | ||
if(iter==maxiter)then | ||
message=trim(message)//'the implicit Euler solution did not converge!' | ||
err=20; return | ||
endif | ||
end do | ||
end subroutine implicitEuler | ||
|
||
end module implicitEulerMod |
Add this suggestion to a batch that can be applied as a single commit.
This suggestion is invalid because no changes were made to the code.
Suggestions cannot be applied while the pull request is closed.
Suggestions cannot be applied while viewing a subset of changes.
Only one suggestion per line can be applied in a batch.
Add this suggestion to a batch that can be applied as a single commit.
Applying suggestions on deleted lines is not supported.
You must change the existing code in this line in order to create a valid suggestion.
Outdated suggestions cannot be applied.
This suggestion has been applied or marked resolved.
Suggestions cannot be applied from pending reviews.
Suggestions cannot be applied on multi-line comments.
Suggestions cannot be applied while the pull request is queued to merge.
Suggestion cannot be applied right now. Please check back later.
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
Comment from @billsacks Sep 7, 2017
Is this a general pattern - that, to get the final flux with implicit Euler, we do a final function call like this? If so, can this final call be incorporated into the implicitEuler routine itself?
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
Comment from @martynpclark Sep 7, 2017
Yes. The final flux should be a function of the state at the end of the time step.
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
Comment from @billsacks Sep 11, 2017
Addressed in NCAR/clm-ctsm#14 in NCAR/clm-ctsm@315c14b
Update: This is now #192 , commit b872d21 . The comment referenced in that commit is the one right above: https://github.com/ESCOMP/ctsm/pull/191/files#r159001381