Skip to content

Commit

Permalink
Added CO2 respiration/assimilation to cross-section and time series.
Browse files Browse the repository at this point in the history
  • Loading branch information
bartvstratum committed Sep 7, 2021
1 parent 4dc4a0b commit 8c7881f
Show file tree
Hide file tree
Showing 3 changed files with 33 additions and 3 deletions.
2 changes: 1 addition & 1 deletion src/modlsm.f90
Original file line number Diff line number Diff line change
Expand Up @@ -21,7 +21,7 @@ module modlsm
use netcdf
implicit none

public :: initlsm, lsm, exitlsm, init_lsm_tiles
public :: initlsm, lsm, exitlsm, init_lsm_tiles, lags

logical :: llsm ! On/off switch LSM
logical :: lfreedrainage ! Free drainage bottom BC for soil moisture
Expand Down
13 changes: 12 additions & 1 deletion src/modlsmcrosssection.f90
Original file line number Diff line number Diff line change
Expand Up @@ -70,6 +70,7 @@ subroutine initlsmcrosssection
use modglobal, only : imax,jmax,ifnamopt,fname_options,dtmax,dtav_glob,ladaptive,j1,dt_lim,cexpnr,tres,btime,checknamelisterror
use modstat_nc, only : lnetcdf,open_nc, define_nc,ncinfo,nctiminfo,writestat_dims_nc
use modsurfdata, only : isurf
use modlsm, only : lags
implicit none

integer :: ierr, ii
Expand Down Expand Up @@ -147,6 +148,7 @@ subroutine initlsmcrosssection
nvar3 = 12
else if (isurf == 11) then
nvar3 = 16
if (lags) nvar3 = nvar3 + 2
end if

allocate(ncname3(nvar3,4))
Expand Down Expand Up @@ -189,6 +191,11 @@ subroutine initlsmcrosssection
call ncinfo(ncname3(14,:),'f2_hv', 'f2(theta) function high vegetation resistance', 's/m', 'tt0t')
call ncinfo(ncname3(15,:),'f2_b', 'f2(theta) function soil resistance', 's/m', 'tt0t')
call ncinfo(ncname3(16,:),'f3', 'f3(VPD) function vegetation resistance', 's/m', 'tt0t')
if (lags) then
call ncinfo(ncname3(17,:),'an_co2', 'Net CO2 assimilation', 'ppb m s-1', 'tt0t')
call ncinfo(ncname3(18,:),'resp_co2', 'CO2 respiration soil', 'ppb m s-1', 'tt0t')
end if

call open_nc(fname3, ncid3, nrec3, n1=imax, n2=jmax)
if (nrec3==0) then
call define_nc(ncid3, 1, tncname3)
Expand Down Expand Up @@ -290,7 +297,7 @@ subroutine wrtsurf
use modglobal, only : imax,jmax,i1,j1,cexpnr,ifoutput,rtimee
use modsurfdata, only : Qnet, H, LE, G0, rs, ra, tskin, tendskin, &
cliq, rsveg, rssoil, Wl, isurf, obl, ustar
use modlsm, only : f1, f2_lv, f2_hv, f2b, f3
use modlsm, only : f1, f2_lv, f2_hv, f2b, f3, lags, an_co2, resp_co2
use modstat_nc, only : lnetcdf, writestat_nc
implicit none

Expand Down Expand Up @@ -379,6 +386,10 @@ subroutine wrtsurf
vars(:,:,14) = f2_hv(2:i1,2:j1)
vars(:,:,15) = f2b(2:i1,2:j1)
vars(:,:,16) = f3(2:i1,2:j1)
if (lags) then
vars(:,:,17) = an_co2(2:i1,2:j1)
vars(:,:,18) = resp_co2(2:i1,2:j1)
endif
end if

call writestat_nc(ncid3, 1, tncname3, (/rtimee/), nrec3, .true.)
Expand Down
21 changes: 20 additions & 1 deletion src/modtimestat.f90
Original file line number Diff line number Diff line change
Expand Up @@ -89,6 +89,7 @@ subroutine inittimestat
use modfields, only : thlprof,qtprof,svprof
use modsurfdata, only : isurf, lhetero, xpatches, ypatches
use modstat_nc, only : lnetcdf, open_nc, define_nc, ncinfo, nctiminfo
use modlsm, only : lags
implicit none
integer :: ierr,k,location = 1
real :: gradient = 0.0
Expand Down Expand Up @@ -234,6 +235,7 @@ subroutine inittimestat
nvar = 32
else if (isurf == 11) then
nvar = 75
if (lags) nvar = nvar + 2
else
nvar = 21
end if
Expand Down Expand Up @@ -399,6 +401,12 @@ subroutine inittimestat
call ncinfo(ncname(vi,:),'qtskin_ws','Skin specific humidity water','kg kg-1','time')
vi = vi+1

if (lags) then
call ncinfo(ncname(vi,:),'an_co2','Net CO2 assimilation','ppb m s-1','time')
vi = vi+1
call ncinfo(ncname(vi,:),'resp_co2','CO2 respiration soil','ppb m s-1','time')
vi = vi+1
end if
end if

call open_nc(fname, ncid,nrec)
Expand Down Expand Up @@ -476,7 +484,7 @@ subroutine timestat
lhetero, xpatches, ypatches, qts_patch, wt_patch, wq_patch, &
thls_patch,obl,z0mav_patch, wco2av, Anav, Respav,gcco2av
use modsurface, only : patchxnr,patchynr
use modlsm, only : tile_lv, tile_hv, tile_bs, tile_ws, tile_aq, f1, f2_lv, f2_hv, f3, f2b
use modlsm, only : tile_lv, tile_hv, tile_bs, tile_ws, tile_aq, f1, f2_lv, f2_hv, f3, f2b, lags, an_co2, resp_co2
use mpi
use modmpi, only : my_real,mpi_sum,mpi_max,mpi_min,comm3d,mpierr,myid
use modstat_nc, only : lnetcdf, writestat_nc,nc_fillvalue
Expand Down Expand Up @@ -508,6 +516,7 @@ subroutine timestat
real :: G_lv_av, G_hv_av, G_bs_av, G_ws_av
real :: thlskin_lv_av, thlskin_hv_av, thlskin_bs_av, thlskin_ws_av, thlskin_aq_av
real :: qtskin_lv_av, qtskin_hv_av, qtskin_bs_av, qtskin_ws_av, qtskin_aq_av
real :: an_co2_av, resp_co2_av

integer:: i, j, k, vi

Expand Down Expand Up @@ -977,6 +986,11 @@ subroutine timestat
qtskin_ws_av = mean_2d(tile_ws%qtskin)
qtskin_aq_av = mean_2d(tile_aq%qtskin)

if (lags) then
an_co2_av = mean_2d(an_co2)
resp_co2_av = mean_2d(resp_co2)
endif

end if

! 9.8 write the results to output file
Expand Down Expand Up @@ -1146,6 +1160,11 @@ subroutine timestat
vars(vi) = qtskin_ws_av; vi = vi+1
vars(vi) = qtskin_aq_av; vi = vi+1

if (lags) then
vars(vi) = an_co2_av; vi = vi+1
vars(vi) = resp_co2_av; vi = vi+1
end if

end if

call writestat_nc(ncid,nvar,ncname,vars,nrec,.true.)
Expand Down

0 comments on commit 8c7881f

Please sign in to comment.