From 5e47c9393baa1263ff023c5efad3574987fc1b5b Mon Sep 17 00:00:00 2001 From: julietbravo Date: Mon, 15 Jun 2020 12:06:48 +0200 Subject: [PATCH] Start at implemention realistic LSM --- src/modlsm.f90 | 157 +++++++++++++++++++++++++++++++++++++++++++++ src/modstartup.f90 | 5 ++ src/modsurface.f90 | 6 ++ src/program.f90 | 6 +- 4 files changed, 172 insertions(+), 2 deletions(-) create mode 100644 src/modlsm.f90 diff --git a/src/modlsm.f90 b/src/modlsm.f90 new file mode 100644 index 00000000..4f8dbb51 --- /dev/null +++ b/src/modlsm.f90 @@ -0,0 +1,157 @@ +! +! Copyright (c) 2020-2020 Wageningen University and Research (WUR) +! +! This file is part of DALES +! +! DALES is free software: you can redistribute it and/or modify +! it under the terms of the GNU General Public License as published by +! the Free Software Foundation, either version 3 of the License, or +! (at your option) any later version. +! +! DALES is distributed in the hope that it will be useful, +! but WITHOUT ANY WARRANTY; without even the implied warranty of +! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the +! GNU General Public License for more details. +! +! You should have received a copy of the GNU General Public License +! along with DALES. If not, see . +! + +module modlsm + use netcdf + implicit none + public :: initlsm, lsm, exitlsm + + ! Land-surface / van Genuchten parameters from NetCDF input table. + real, allocatable :: & + theta_res(:), theta_wp(:), theta_fc(:), theta_sat(:), gamma_sat(:), vg_a(:), vg_l(:), vg_n(:) +contains + +! +! Initialise the land-surface model +! +subroutine initlsm + implicit none + + ! Read the soil parameter table + call read_soil_table +end subroutine initlsm + +subroutine lsm + implicit none +end subroutine lsm + +! +! Cleanup (deallocate) the land-surface model +! +subroutine exitlsm + implicit none + + deallocate( theta_res, theta_wp, theta_fc, theta_sat, gamma_sat, vg_a, vg_l, vg_n ) +end subroutine exitlsm + +! +! Read the input table with the (van Genuchten) soil parameters +! +subroutine read_soil_table + implicit none + integer :: n, ncid, dimid, varid + + ! Open the NetCDF file and read the table size + print*,'Reading "van_genuchten_parameters.nc"' + call check( nf90_open('van_genuchten_parameters.nc', nf90_nowrite, ncid) ) + call check( nf90_inq_dimid(ncid, 'index', dimid) ) + call check( nf90_inquire_dimension(ncid, dimid, len=n) ) + + ! Allocate variables + allocate( & + theta_res(n), theta_wp(n), theta_fc(n), theta_sat(n), & + gamma_sat(n), vg_a(n), vg_l(n), vg_n(n) ) + + ! Read variables + call check( nf90_inq_varid(ncid, 'theta_res', varid) ) + call check( nf90_get_var(ncid, varid, theta_res) ) + + call check( nf90_inq_varid(ncid, 'theta_wp', varid) ) + call check( nf90_get_var(ncid, varid, theta_wp) ) + + call check( nf90_inq_varid(ncid, 'theta_fc', varid) ) + call check( nf90_get_var(ncid, varid, theta_fc) ) + + call check( nf90_inq_varid(ncid, 'theta_sat', varid) ) + call check( nf90_get_var(ncid, varid, theta_sat) ) + + call check( nf90_inq_varid(ncid, 'gamma_sat', varid) ) + call check( nf90_get_var(ncid, varid, gamma_sat) ) + + call check( nf90_inq_varid(ncid, 'alpha', varid) ) + call check( nf90_get_var(ncid, varid, vg_a) ) + + call check( nf90_inq_varid(ncid, 'l', varid) ) + call check( nf90_get_var(ncid, varid, vg_l) ) + + call check( nf90_inq_varid(ncid, 'n', varid) ) + call check( nf90_get_var(ncid, varid, vg_n) ) + + call check( nf90_close(ncid) ) + +end subroutine read_soil_table + +! +! Convert soil hydraulic head to soil water content, using van Genuchten parameterisation. +! +pure function psi_to_theta(theta_res, theta_sat, vg_a, vg_n, vg_m, psi) result(res) + implicit none + real, intent(in) :: theta_res, theta_sat, vg_a, vg_n, vg_m, psi + real :: res + + res = theta_res + (theta_sat - theta_res) * (1. / (1.+ abs(vg_a * psi)**vg_n))**vg_m +end function psi_to_theta + +! +! Convert soil water content to hydraulic head, using van Genuchten parameterisation. +! +pure function theta_to_psi(theta_res, theta_sat, vg_a, vg_n, vg_m, theta) result(res) + implicit none + real, intent(in) :: theta_res, theta_sat, vg_a, vg_n, vg_m, theta + real :: res + + res = -(vg_a**(-vg_n) * (-1. + ((theta_res - theta_sat)/(theta_res - theta))**(1./vg_m)))**(1./vg_n) +end function theta_to_psi + +! +! Calculate hydraulic diffusivity using van Genuchten parameterisation. +! +pure function calc_diffusivity_vg( & + theta_norm, vg_a, vg_l, vg_m, lambda_sat, theta_sat, theta_res) result(res) + implicit none + real, intent(in) :: theta_norm, vg_a, vg_l, vg_m, lambda_sat, theta_sat, theta_res + real :: res + + res = (1.-vg_m)*lambda_sat / (vg_a * vg_m * (theta_sat-theta_res)) * theta_norm**(vg_l-(1./vg_m)) * & + ( (1.-theta_norm**(1./vg_m))**(-vg_m) + (1.-theta_norm**(1/vg_m))**vg_m - 2. ) +end function calc_diffusivity_vg + +! +! Calculate hydraulic conductivity using van Genuchten parameterisation. +! +pure function calc_conductivity_vg(theta_norm, vg_l, vg_m, lambda_sat) result(res) + implicit none + real, intent(in) :: theta_norm, vg_l, vg_m, lambda_sat + real :: res + + res = lambda_sat * theta_norm**vg_l * ( 1.- (1.-theta_norm**(1./vg_m))**vg_m )**2. +end function calc_conductivity_vg + +! +! Check NetCDF calls +! +subroutine check(status) + integer, intent (in) :: status + if(status /= nf90_noerr) then + print *,'NetCDF error: ', trim(nf90_strerror(status)) + stop + end if +end subroutine check + +end module modlsm diff --git a/src/modstartup.f90 b/src/modstartup.f90 index b98534cc..1b59af6b 100644 --- a/src/modstartup.f90 +++ b/src/modstartup.f90 @@ -63,6 +63,7 @@ subroutine startup(path) use modforces, only : lforce_user use modsurfdata, only : z0,ustin,wtsurf,wqsurf,wsvsurf,ps,thls,isurf use modsurface, only : initsurface + use modlsm, only : initlsm use modfields, only : initfields use modpois, only : initpois use modradiation, only : initradiation @@ -268,6 +269,7 @@ subroutine startup(path) call initradiation call initchem call initsurface + call initlsm call initsubgrid call initmicrophysics @@ -348,6 +350,7 @@ subroutine checkinitvalues !isurf if (myid == 0) then select case (isurf) + case(-1) case(1) case(2,10) case(3:4) @@ -1143,11 +1146,13 @@ subroutine exitmodules use modradiation, only : exitradiation use modsubgrid, only : exitsubgrid use modsurface, only : exitsurface + use modlsm, only : exitlsm use modthermodynamics, only : exitthermodynamics call exittimedep call exitthermodynamics call exitsurface + call exitlsm call exitsubgrid call exitradiation call exitpois diff --git a/src/modsurface.f90 b/src/modsurface.f90 index b9934c03..9d372a1c 100644 --- a/src/modsurface.f90 +++ b/src/modsurface.f90 @@ -320,6 +320,8 @@ subroutine initsurface close(ifinput) else select case (isurf) + case (-1) + return case (1) ! Interactive land surface open (ifinput,file='surface.interactive.inp.'//cexpnr) ierr = 0 @@ -724,6 +726,10 @@ subroutine surface patchx = 0 patchy = 0 + if (isurf==-1) then + return + end if + if (isurf==10) then call surf_user return diff --git a/src/program.f90 b/src/program.f90 index 886aaded..6db67d45 100644 --- a/src/program.f90 +++ b/src/program.f90 @@ -107,7 +107,8 @@ program DALES use modboundary, only : boundary, grwdamp! JvdD ,tqaver use modthermodynamics, only : thermodynamics use modmicrophysics, only : microsources - use modsurface, only : surface + use modsurface, only : surface + use modlsm, only : lsm use modsubgrid, only : subgrid use modforces, only : forces, coriolis, lstend use modradiation, only : radiation @@ -211,8 +212,9 @@ program DALES call samptend(tend_rad) !----------------------------------------------------- -! 3.2 THE SURFACE LAYER +! 3.2 THE SURFACE LAYER / LAND-SURFACE !----------------------------------------------------- + call lsm call surface !-----------------------------------------------------