From b6249221bb89e298237f8e46e8fc6db1712b1512 Mon Sep 17 00:00:00 2001 From: julietbravo Date: Tue, 7 Jul 2020 11:53:21 +0200 Subject: [PATCH] Added restart functionality to new LSM --- src/modlsm.f90 | 29 ++++++++++++++-------- src/modstartup.f90 | 62 ++++++++++++++++++++++++++++++++++++++++++++-- 2 files changed, 78 insertions(+), 13 deletions(-) diff --git a/src/modlsm.f90 b/src/modlsm.f90 index 9f6dd228..9b1cbd06 100644 --- a/src/modlsm.f90 +++ b/src/modlsm.f90 @@ -105,7 +105,6 @@ module modlsm contains subroutine lsm - use modsurfdata, only : thlflux, qtflux, G0 implicit none if (.not. llsm) return @@ -1130,7 +1129,7 @@ end subroutine init_lsm_tiles ! Initialise the LSM homogeneous from namelist input ! subroutine init_homogeneous - use modglobal, only : ifnamopt, fname_options, checknamelisterror + use modglobal, only : ifnamopt, fname_options, checknamelisterror, lwarmstart use modmpi, only : myid, comm3d, mpierr, mpi_logical, my_real, mpi_integer use modsurfdata, only : tsoil, tsoilm, phiw, phiwm, wl, wlm, wmax implicit none @@ -1248,23 +1247,31 @@ subroutine init_homogeneous gD(:,:) = gD_high + if (.not. lwarmstart) then + ! Init prognostic variables in case of cold start. + ! For a warm start, these are read from restartfiles + ! in modstartup. + do k=1, kmax_soil + tsoil(:,:,k) = t_soil_p(k) + phiw (:,:,k) = theta_soil_p(k) + end do + + tsoilm(:,:,:) = tsoil(:,:,:) + phiwm (:,:,:) = phiw (:,:,:) + + wl(:,:) = 0. + wlm(:,:) = 0. + end if + do k=1, kmax_soil - tsoil(:,:,k) = t_soil_p(k) - phiw (:,:,k) = theta_soil_p(k) soil_index(:,:,k) = soil_index_p(k) end do - tsoilm(:,:,:) = tsoil(:,:,:) - phiwm (:,:,:) = phiw (:,:,:) - ! Set properties wet skin tile - wl(:,:) = 0 - wlm(:,:) = 0 - tile_ws % z0m(:,:) = c_low*z0m_low + c_high*z0m_high + c_bare*z0m_bare tile_ws % z0h(:,:) = c_low*z0h_low + c_high*z0h_high + c_bare*z0h_bare - tile_ws % lambda_stable(:,:) = c_low*lambda_s_low + c_high*lambda_s_high + c_bare*lambda_s_bare + tile_ws % lambda_stable(:,:) = c_low*lambda_s_low + c_high*lambda_s_high + c_bare*lambda_s_bare tile_ws % lambda_unstable(:,:) = c_low*lambda_us_low + c_high*lambda_us_high + c_bare*lambda_us_bare ! Max liquid water per grid point, accounting for LAI diff --git a/src/modstartup.f90 b/src/modstartup.f90 index e5b91eb5..beb2549d 100644 --- a/src/modstartup.f90 +++ b/src/modstartup.f90 @@ -848,6 +848,7 @@ subroutine readrestartfiles tres,ifinput,nsv,dt use modmpi, only : cmyid use modsubgriddata, only : ekm,ekh + use modlsm, only : kmax_soil, tile_lv, tile_hv, tile_bs, tile_ws character(50) :: name @@ -956,7 +957,35 @@ subroutine readrestartfiles end if read(ifinput) timee close(ifinput) + + else if (isurf == 11) then + name(5:5) = 'l' + write(6,*) 'loading ',name + open(unit=ifinput,file=name,form='unformatted') + read(ifinput) (((tsoil(i,j,k), i=1,i2), j=1,j2), k=1,kmax_soil) + read(ifinput) (((phiw (i,j,k), i=1,i2), j=1,j2), k=1,kmax_soil) + read(ifinput) ((tskin (i,j), i=1,i2), j=1,j2) + read(ifinput) ((Wl (i,j), i=1,i2), j=1,j2) + + read(ifinput) ((tile_lv%thlskin(i,j), i=1,i2), j=1,j2) + read(ifinput) ((tile_hv%thlskin(i,j), i=1,i2), j=1,j2) + read(ifinput) ((tile_bs%thlskin(i,j), i=1,i2), j=1,j2) + read(ifinput) ((tile_ws%thlskin(i,j), i=1,i2), j=1,j2) + + read(ifinput) ((tile_lv%qtskin(i,j), i=1,i2), j=1,j2) + read(ifinput) ((tile_hv%qtskin(i,j), i=1,i2), j=1,j2) + read(ifinput) ((tile_bs%qtskin(i,j), i=1,i2), j=1,j2) + read(ifinput) ((tile_ws%qtskin(i,j), i=1,i2), j=1,j2) + + read(ifinput) ((tile_lv%obuk(i,j), i=1,i2), j=1,j2) + read(ifinput) ((tile_hv%obuk(i,j), i=1,i2), j=1,j2) + read(ifinput) ((tile_bs%obuk(i,j), i=1,i2), j=1,j2) + read(ifinput) ((tile_ws%obuk(i,j), i=1,i2), j=1,j2) + + read(ifinput) timee + close(ifinput) end if + end subroutine readrestartfiles !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! @@ -996,13 +1025,13 @@ subroutine do_writerestartfiles use modglobal, only : i1,i2,ih,j1,j2,jh,k1,dsv,cexpnr,ifoutput,timee,rtimee,tres,nsv,dtheta,dqt,dt use modmpi, only : cmyid,myid use modsubgriddata, only : ekm,ekh + use modlsm, only : kmax_soil, tile_lv, tile_hv, tile_bs, tile_ws implicit none integer imin,ihour integer i,j,k,n character(50) name,linkname - ihour = floor(rtimee/3600) imin = floor((rtimee-ihour * 3600) /3600. * 60.) name = 'initdXXXhXXmXXXXXXXX.XXX' @@ -1114,8 +1143,37 @@ subroutine do_writerestartfiles linkname = name linkname(6:11) = "latest" call system("cp "//name //" "//linkname) - end if + else if (isurf == 11) then + name(5:5)='l' + open (ifoutput,file=name,form='unformatted') + write(ifoutput) (((tsoil(i,j,k), i=1,i2), j=1,j2), k=1,kmax_soil) + write(ifoutput) (((phiw (i,j,k), i=1,i2), j=1,j2), k=1,kmax_soil) + write(ifoutput) ((tskin (i,j), i=1,i2), j=1,j2) + write(ifoutput) ((Wl (i,j), i=1,i2), j=1,j2) + + ! Sub-grid tiles + write(ifoutput) ((tile_lv%thlskin(i,j), i=1,i2), j=1,j2) + write(ifoutput) ((tile_hv%thlskin(i,j), i=1,i2), j=1,j2) + write(ifoutput) ((tile_bs%thlskin(i,j), i=1,i2), j=1,j2) + write(ifoutput) ((tile_ws%thlskin(i,j), i=1,i2), j=1,j2) + + write(ifoutput) ((tile_lv%qtskin(i,j), i=1,i2), j=1,j2) + write(ifoutput) ((tile_hv%qtskin(i,j), i=1,i2), j=1,j2) + write(ifoutput) ((tile_bs%qtskin(i,j), i=1,i2), j=1,j2) + write(ifoutput) ((tile_ws%qtskin(i,j), i=1,i2), j=1,j2) + + write(ifoutput) ((tile_lv%obuk(i,j), i=1,i2), j=1,j2) + write(ifoutput) ((tile_hv%obuk(i,j), i=1,i2), j=1,j2) + write(ifoutput) ((tile_bs%obuk(i,j), i=1,i2), j=1,j2) + write(ifoutput) ((tile_ws%obuk(i,j), i=1,i2), j=1,j2) + + write(ifoutput) timee + close (ifoutput) + linkname = name + linkname(6:11) = "latest" + call system("cp "//name //" "//linkname) + end if if (myid==0) then write(*,'(A,F15.7,A,I4)') 'dump at time = ',rtimee,' unit = ',ifoutput