From db40554172aa43424117c3df69d4c791d17efec3 Mon Sep 17 00:00:00 2001 From: julietbravo Date: Mon, 29 Jun 2020 12:23:47 +0200 Subject: [PATCH] Added data structures for land-surface tiles --- src/modlsm.f90 | 107 ++++++++++++++++++++++++++++++++++++++++++--- src/modsurface.f90 | 35 +++++++-------- 2 files changed, 118 insertions(+), 24 deletions(-) diff --git a/src/modlsm.f90 b/src/modlsm.f90 index 4f8dbb51..0b9c5bea 100644 --- a/src/modlsm.f90 +++ b/src/modlsm.f90 @@ -20,26 +20,77 @@ 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(:) + + logical :: sw_lsm + logical :: sw_homogeneous + + ! Data structure for sub-grid tiles + type lsm_tile + ! Static properties: + real, allocatable :: z0m(:,:), z0h(:,:) + ! Dynamic tile fraction: + real, allocatable :: frac(:,:) + ! Monin-obukhov / surface layer: + real, allocatable :: obuk(:,:), ustar(:,:), ra(:,:) + ! Conductivity skin layer: + real, allocatable :: lambda_stable(:,:), lamda_unstable(:,:) + ! Surface fluxes: + real, allocatable :: H(:,:), LE(:,:), G(:,:), wthl(:,:), wqt(:,:) + ! Surface temperature and humidity + real, allocatable :: tskin(:,:), qskin(:,:) + end type lsm_tile + + type(lsm_tile) low_veg, high_veg, bare_soil, wet_skin + contains +subroutine lsm + implicit none +end subroutine lsm + ! ! Initialise the land-surface model ! subroutine initlsm + use modglobal, only : ifnamopt, fname_options, checknamelisterror + use modmpi, only : myid, comm3d, mpierr, mpi_logical + use modsurfdata, only : isurf implicit none - ! Read the soil parameter table - call read_soil_table -end subroutine initlsm + integer :: ierr -subroutine lsm - implicit none -end subroutine lsm + ! Read namelist + namelist /NAMLSM/ & + sw_lsm, sw_homogeneous + + if (myid == 0) then + open(ifnamopt, file=fname_options, status='old', iostat=ierr) + read(ifnamopt, NAMLSM, iostat=ierr) + write(6, NAMLSM) + close(ifnamopt) + end if + + ! Broadcast namelist values to all MPI tasks + call MPI_BCAST(sw_lsm, 1, mpi_logical, 0, comm3d, mpierr) + call MPI_BCAST(sw_homogeneous, 1, mpi_logical, 0, comm3d, mpierr) + + if (sw_lsm) then + ! Allocate required fields in modsurfacedata, + ! and arrays / tiles from this module + call allocate_fields + + ! Read the soil parameter table + call read_soil_table + + end if + +end subroutine initlsm ! ! Cleanup (deallocate) the land-surface model @@ -50,6 +101,50 @@ subroutine exitlsm deallocate( theta_res, theta_wp, theta_fc, theta_sat, gamma_sat, vg_a, vg_l, vg_n ) end subroutine exitlsm +! +! Allocate all LSM fields +! +subroutine allocate_fields + use modglobal, only : i2, j2 + implicit none + + ! Allocate the tiled variables + call allocate_tile(low_veg) + call allocate_tile(high_veg) + call allocate_tile(bare_soil) + call allocate_tile(wet_skin) + +end subroutine allocate_fields + +! +! Allocate all fields of a LSM tile +! +subroutine allocate_tile(tile) + use modglobal, only : i2, j2 + implicit none + type(lsm_tile), intent(inout) :: tile + + ! Static properties: + allocate(tile%z0m (i2, j2)) + allocate(tile%z0h (i2, j2)) + ! Dynamic tile fraction: + allocate(tile%frac (i2, j2)) + ! Monin-obukhov / surface layer: + allocate(tile%obuk (i2, j2)) + allocate(tile%ustar(i2, j2)) + allocate(tile%ra (i2, j2)) + ! Surface fluxes: + allocate(tile%H (i2, j2)) + allocate(tile%LE (i2, j2)) + allocate(tile%G (i2, j2)) + allocate(tile%wthl (i2, j2)) + allocate(tile%wqt (i2, j2)) + ! Surface temperature and humidity + allocate(tile%tskin(i2, j2)) + allocate(tile%qskin(i2, j2)) + +end subroutine allocate_tile + ! ! Read the input table with the (van Genuchten) soil parameters ! diff --git a/src/modsurface.f90 b/src/modsurface.f90 index 9d372a1c..99d83aee 100644 --- a/src/modsurface.f90 +++ b/src/modsurface.f90 @@ -80,9 +80,10 @@ subroutine initsurface integer :: i,j,k, landindex, ierr, defined_landtypes, landtype_0 = -1 integer :: tempx,tempy - character(len=1500) :: readbuffer + character(len=1500) :: readbuffer + namelist/NAMSURFACE/ & !< Soil related variables - isurf,tsoilav, tsoildeepav, phiwav, rootfav, & + isurf, tsoilav, tsoildeepav, phiwav, rootfav, & ! Land surface related variables lmostlocal, lsmoothflux, lneutral, z0mav, z0hav, rsisurf2, Cskinav, lambdaskinav, albedoav, Qnetav, cvegav, Wlav, & ! Jarvis-Steward related variables @@ -163,7 +164,7 @@ subroutine initsurface call MPI_BCAST(phiwp , 1, MY_REAL , 0, comm3d, mpierr) call MPI_BCAST(R10 , 1, MY_REAL , 0, comm3d, mpierr) call MPI_BCAST(lsplitleaf , 1, MPI_LOGICAL, 0, comm3d, mpierr) - + call MPI_BCAST(land_use(1:mpatch,1:mpatch),mpatch*mpatch, MPI_INTEGER, 0, comm3d, mpierr) if(lCO2Ags .and. (.not. lrsAgs)) then @@ -320,8 +321,6 @@ 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 @@ -696,7 +695,7 @@ subroutine initsurface if (lsplitleaf) then allocate(PARdirField (2:i1,2:j1)) allocate(PARdifField (2:i1,2:j1)) - endif + endif endif return end subroutine initsurface @@ -868,7 +867,7 @@ subroutine surface phimzf = phim(zf(1)/obl(i,j)) phihzf = phih(zf(1)/obl(i,j)) - + dudz (i,j) = ustar(i,j) * phimzf / (fkar*zf(1))*(upcu/horv) dvdz (i,j) = ustar(i,j) * phimzf / (fkar*zf(1))*(vpcv/horv) dthldz(i,j) = - thlflux(i,j) / ustar(i,j) * phihzf / (fkar*zf(1)) @@ -901,7 +900,7 @@ subroutine surface phimzf = phim(zf(1)/obl(i,j)) phihzf = phih(zf(1)/obl(i,j)) - + upcu = 0.5 * (u0(i,j,1) + u0(i+1,j,1)) + cu vpcv = 0.5 * (v0(i,j,1) + v0(i,j+1,1)) + cv horv = sqrt(upcu ** 2. + vpcv ** 2.) @@ -984,10 +983,10 @@ subroutine surface svflux(i,j,n) = wsvsurf(n) enddo endif - + phimzf = phim(zf(1)/obl(i,j)) phihzf = phih(zf(1)/obl(i,j)) - + dudz (i,j) = ustar(i,j) * phimzf / (fkar*zf(1))*(upcu/horv) dvdz (i,j) = ustar(i,j) * phimzf / (fkar*zf(1))*(vpcv/horv) dthldz(i,j) = - thlflux(i,j) / ustar(i,j) * phihzf / (fkar*zf(1)) @@ -1162,7 +1161,7 @@ subroutine getobl if(Rib > 0) L = 0.01 if(Rib < 0) L = -0.01 end if - + do while (.true.) iter = iter + 1 Lold = L @@ -1301,7 +1300,7 @@ subroutine getobl if(Rib > 0) L = 0.01 if(Rib < 0) L = -0.01 end if - + do while (.true.) iter = iter + 1 Lold = L @@ -1378,7 +1377,7 @@ end function psih ! stability function Phi for momentum. ! Many functional forms of Phi have been suggested, see e.g. Optis 2015 - ! Phi and Psi above are related by an integral and should in principle match, + ! Phi and Psi above are related by an integral and should in principle match, ! currently they do not. ! FJ 2018: For very stable situations, zeta > 1 add cap to phi - the linear expression is valid only for zeta < 1 function phim(zeta) @@ -1398,7 +1397,7 @@ function phim(zeta) return end function phim - ! stability function Phi for heat. + ! stability function Phi for heat. function phih(zeta) implicit none real :: phih @@ -1416,7 +1415,7 @@ function phih(zeta) return end function phih - + function E1(x) implicit none real :: E1 @@ -1428,7 +1427,7 @@ function E1(x) do k=1,99 !E1sum = E1sum + (-1.0) ** (k + 0.0) * x ** (k + 0.0) / ( (k + 0.0) * factorial(k) ) E1sum = E1sum + (-1.0 * x) ** k / ( k * factorial(k) ) ! FJ changed this for compilation with cray fortran - + end do E1 = -0.57721566490153286060 - log(x) - E1sum @@ -1668,7 +1667,7 @@ subroutine do_lsm real :: Ag, PARdir, PARdif !Variables for 2leaf AGS real :: MW_Air = 28.97 real :: MW_CO2 = 44 - + real :: sinbeta, kdrbl, kdf, kdr, ref, ref_dir real :: iLAI, fSL real :: PARdfU, PARdfD, PARdfT, PARdrU, PARdrD, PARdrT, dirPAR, difPAR @@ -1960,7 +1959,7 @@ subroutine do_lsm gc_inf = LAI(i,j) * sum(weight_g * gnet) else !lsplitleaf - + ! Calculate upscaling from leaf to canopy: net flow CO2 into the plant (An) AGSa1 = 1.0 / (1 - f0) Dstar = D0 / (AGSa1 * (f0 - fmin))