Skip to content

Commit a9606e7

Browse files
committed
LST-DA: all temperature layers updated
1 parent d738f97 commit a9606e7

2 files changed

Lines changed: 48 additions & 24 deletions

File tree

bldsva/intf_DA/pdaf/model/clm5_0/enkf_clm_mod_5.F90

Lines changed: 24 additions & 12 deletions
Original file line numberDiff line numberDiff line change
@@ -86,6 +86,7 @@ subroutine define_clm_statevec(mype)
8686
use shr_kind_mod, only: r8 => shr_kind_r8
8787
use decompMod , only : get_proc_bounds
8888
use clm_varpar , only : nlevsoi
89+
use clm_varpar , only : nlevgrnd
8990
use clm_varcon , only : ispval
9091
use ColumnType , only : col
9192
use PatchType , only : patch
@@ -97,6 +98,7 @@ subroutine define_clm_statevec(mype)
9798
integer :: i
9899
integer :: j
99100
integer :: jj
101+
integer :: lev
100102
integer :: p
101103
integer :: c
102104
integer :: g
@@ -358,7 +360,7 @@ subroutine define_clm_statevec(mype)
358360

359361
clm_varsize = endp-begp+1
360362
! clm_paramsize = endp-begp+1 !LAI
361-
clm_statevecsize = 3* (endp-begp+1) !TSKIN, then TSOIL and TV
363+
clm_statevecsize = (1 + nlevgrnd + 1)* (endp-begp+1) !TSKIN, then nlevgrnd times TSOIL and TV
362364

363365
IF (allocated(state_pdaf2clm_p_p)) deallocate(state_pdaf2clm_p_p)
364366
allocate(state_pdaf2clm_p_p(clm_statevecsize))
@@ -374,13 +376,15 @@ subroutine define_clm_statevec(mype)
374376
state_pdaf2clm_p_p(cc) = p !TSKIN
375377
state_pdaf2clm_c_p(cc) = patch%column(p) !TSKIN
376378
state_pdaf2clm_j_p(cc) = 1
377-
state_pdaf2clm_p_p(cc+clm_varsize) = p !TSOIL
378-
state_pdaf2clm_c_p(cc+clm_varsize) = patch%column(p) !TSOIL
379-
state_pdaf2clm_j_p(cc+clm_varsize) = 1
380-
state_pdaf2clm_p_p(cc+2*clm_varsize) = p !TV
381-
state_pdaf2clm_c_p(cc+2*clm_varsize) = patch%column(p) !TV
382-
state_pdaf2clm_j_p(cc+2*clm_varsize) = 1
383-
end do
379+
do lev=1,nlevgrnd
380+
! ivar = 2-26: TSOIL
381+
state_pdaf2clm_p_p(cc + lev*clm_varsize) = p
382+
state_pdaf2clm_c_p(cc + lev*clm_varsize) = patch%column(p)
383+
state_pdaf2clm_j_p(cc + lev*clm_varsize) = lev
384+
end do
385+
state_pdaf2clm_p_p(cc+(1+nlevgrnd)*clm_varsize) = p !TV
386+
state_pdaf2clm_c_p(cc+(1+nlevgrnd)*clm_varsize) = patch%column(p) !TV
387+
state_pdaf2clm_j_p(cc+(1+nlevgrnd)*clm_varsize) = 1
384388

385389
endif
386390

@@ -423,6 +427,7 @@ subroutine set_clm_statevec(tstartcycle, mype)
423427
use clm_instMod, only : temperature_inst
424428
use clm_instMod, only : canopystate_inst
425429
use clm_varpar , only : nlevsoi
430+
use clm_varpar , only : nlevgrnd
426431
! use clm_varcon, only: nameg, namec
427432
! use GetGlobalValuesMod, only: GetGlobalWrite
428433
use ColumnType , only : col
@@ -441,6 +446,7 @@ subroutine set_clm_statevec(tstartcycle, mype)
441446
real(r8), pointer :: t_skin(:)
442447
real(r8), pointer :: tlai(:)
443448
integer :: i,j,jj,g,cc=0,offset=0
449+
integer :: lev
444450
character (len = 34) :: fn !TSMP-PDAF: function name for state vector output
445451
character (len = 34) :: fn2 !TSMP-PDAF: function name for swc output
446452

@@ -523,8 +529,10 @@ subroutine set_clm_statevec(tstartcycle, mype)
523529
do cc = 1, clm_varsize
524530
! t_skin iterated over patches
525531
clm_statevec(cc) = t_skin(state_pdaf2clm_p_p(cc))
526-
clm_statevec(cc+clm_varsize) = t_soisno(state_pdaf2clm_c_p(cc+clm_varsize), state_pdaf2clm_j_p(cc+clm_varsize))
527-
clm_statevec(cc+2*clm_varsize) = t_veg(state_pdaf2clm_p_p(cc+2*clm_varsize))
532+
do lev=1,nlevgrnd
533+
clm_statevec(cc+lev*clm_varsize) = t_soisno(state_pdaf2clm_c_p(cc+lev*clm_varsize), state_pdaf2clm_j_p(cc+lev*clm_varsize))
534+
end do
535+
clm_statevec(cc+(1+nlevgrnd)*clm_varsize) = t_veg(state_pdaf2clm_p_p(cc+(1+nlevgrnd)*clm_varsize))
528536
end do
529537
endif
530538

@@ -565,6 +573,7 @@ end subroutine set_clm_statevec
565573

566574
subroutine update_clm(tstartcycle, mype) bind(C,name="update_clm")
567575
use clm_varpar , only : nlevsoi
576+
use clm_varpar , only : nlevgrnd
568577
use clm_time_manager , only : update_DA_nstep
569578
use shr_kind_mod , only : r8 => shr_kind_r8
570579
use PatchType , only : patch
@@ -601,6 +610,7 @@ subroutine update_clm(tstartcycle, mype) bind(C,name="update_clm")
601610
real(r8) :: swc_update ! updated SWC in loop
602611

603612
integer :: i,j,jj,g,c,p,cc=0,offset=0
613+
integer :: lev
604614
character (len = 31) :: fn !TSMP-PDAF: function name for state vector outpu
605615
character (len = 31) :: fn2 !TSMP-PDAF: function name for state vector outpu
606616
character (len = 32) :: fn3 !TSMP-PDAF: function name for state vector outpu
@@ -782,8 +792,10 @@ subroutine update_clm(tstartcycle, mype) bind(C,name="update_clm")
782792
do p = clm_begp, clm_endp
783793
c = patch%column(p)
784794
t_skin(p) = clm_statevec(state_clm2pdaf_p(p,1))
785-
t_soisno(c,1) = clm_statevec(state_clm2pdaf_p(p,1) + clm_varsize)
786-
t_veg(p) = clm_statevec(state_clm2pdaf_p(p,1) + 2*clm_varsize)
795+
do lev=1,nlevgrnd
796+
t_soisno(c,lev) = clm_statevec(state_clm2pdaf_p(p,1) + (lev)*clm_varsize)
797+
end do
798+
t_veg(p) = clm_statevec(state_clm2pdaf_p(p,1) + (1+nlevgrnd)*clm_varsize)
787799
end do
788800
endif
789801

bldsva/intf_DA/pdaf/model/eclm/enkf_clm_mod_5.F90

Lines changed: 24 additions & 12 deletions
Original file line numberDiff line numberDiff line change
@@ -86,6 +86,7 @@ subroutine define_clm_statevec(mype)
8686
use shr_kind_mod, only: r8 => shr_kind_r8
8787
use decompMod , only : get_proc_bounds
8888
use clm_varpar , only : nlevsoi
89+
use clm_varpar , only : nlevgrnd
8990
use clm_varcon , only : ispval
9091
use ColumnType , only : col
9192
use PatchType , only : patch
@@ -97,6 +98,7 @@ subroutine define_clm_statevec(mype)
9798
integer :: i
9899
integer :: j
99100
integer :: jj
101+
integer :: lev
100102
integer :: p
101103
integer :: c
102104
integer :: g
@@ -358,7 +360,7 @@ subroutine define_clm_statevec(mype)
358360

359361
clm_varsize = endp-begp+1
360362
! clm_paramsize = endp-begp+1 !LAI
361-
clm_statevecsize = 3* (endp-begp+1) !TSKIN, then TSOIL and TV
363+
clm_statevecsize = (1 + nlevgrnd + 1)* (endp-begp+1) !TSKIN, then nlevgrnd times TSOIL and TV
362364

363365
IF (allocated(state_pdaf2clm_p_p)) deallocate(state_pdaf2clm_p_p)
364366
allocate(state_pdaf2clm_p_p(clm_statevecsize))
@@ -374,13 +376,15 @@ subroutine define_clm_statevec(mype)
374376
state_pdaf2clm_p_p(cc) = p !TSKIN
375377
state_pdaf2clm_c_p(cc) = patch%column(p) !TSKIN
376378
state_pdaf2clm_j_p(cc) = 1
377-
state_pdaf2clm_p_p(cc+clm_varsize) = p !TSOIL
378-
state_pdaf2clm_c_p(cc+clm_varsize) = patch%column(p) !TSOIL
379-
state_pdaf2clm_j_p(cc+clm_varsize) = 1
380-
state_pdaf2clm_p_p(cc+2*clm_varsize) = p !TV
381-
state_pdaf2clm_c_p(cc+2*clm_varsize) = patch%column(p) !TV
382-
state_pdaf2clm_j_p(cc+2*clm_varsize) = 1
383-
end do
379+
do lev=1,nlevgrnd
380+
! ivar = 2-26: TSOIL
381+
state_pdaf2clm_p_p(cc + lev*clm_varsize) = p
382+
state_pdaf2clm_c_p(cc + lev*clm_varsize) = patch%column(p)
383+
state_pdaf2clm_j_p(cc + lev*clm_varsize) = lev
384+
end do
385+
state_pdaf2clm_p_p(cc+(1+nlevgrnd)*clm_varsize) = p !TV
386+
state_pdaf2clm_c_p(cc+(1+nlevgrnd)*clm_varsize) = patch%column(p) !TV
387+
state_pdaf2clm_j_p(cc+(1+nlevgrnd)*clm_varsize) = 1
384388

385389
endif
386390

@@ -423,6 +427,7 @@ subroutine set_clm_statevec(tstartcycle, mype)
423427
use clm_instMod, only : temperature_inst
424428
use clm_instMod, only : canopystate_inst
425429
use clm_varpar , only : nlevsoi
430+
use clm_varpar , only : nlevgrnd
426431
! use clm_varcon, only: nameg, namec
427432
! use GetGlobalValuesMod, only: GetGlobalWrite
428433
use ColumnType , only : col
@@ -441,6 +446,7 @@ subroutine set_clm_statevec(tstartcycle, mype)
441446
real(r8), pointer :: t_skin(:)
442447
real(r8), pointer :: tlai(:)
443448
integer :: i,j,jj,g,cc=0,offset=0
449+
integer :: lev
444450
character (len = 34) :: fn !TSMP-PDAF: function name for state vector output
445451
character (len = 34) :: fn2 !TSMP-PDAF: function name for swc output
446452

@@ -523,8 +529,10 @@ subroutine set_clm_statevec(tstartcycle, mype)
523529
do cc = 1, clm_varsize
524530
! t_skin iterated over patches
525531
clm_statevec(cc) = t_skin(state_pdaf2clm_p_p(cc))
526-
clm_statevec(cc+clm_varsize) = t_soisno(state_pdaf2clm_c_p(cc+clm_varsize), state_pdaf2clm_j_p(cc+clm_varsize))
527-
clm_statevec(cc+2*clm_varsize) = t_veg(state_pdaf2clm_p_p(cc+2*clm_varsize))
532+
do lev=1,nlevgrnd
533+
clm_statevec(cc+lev*clm_varsize) = t_soisno(state_pdaf2clm_c_p(cc+lev*clm_varsize), state_pdaf2clm_j_p(cc+lev*clm_varsize))
534+
end do
535+
clm_statevec(cc+(1+nlevgrnd)*clm_varsize) = t_veg(state_pdaf2clm_p_p(cc+(1+nlevgrnd)*clm_varsize))
528536
end do
529537
endif
530538

@@ -565,6 +573,7 @@ end subroutine set_clm_statevec
565573

566574
subroutine update_clm(tstartcycle, mype) bind(C,name="update_clm")
567575
use clm_varpar , only : nlevsoi
576+
use clm_varpar , only : nlevgrnd
568577
use clm_time_manager , only : update_DA_nstep
569578
use shr_kind_mod , only : r8 => shr_kind_r8
570579
use PatchType , only : patch
@@ -601,6 +610,7 @@ subroutine update_clm(tstartcycle, mype) bind(C,name="update_clm")
601610
real(r8) :: swc_update ! updated SWC in loop
602611

603612
integer :: i,j,jj,g,c,p,cc=0,offset=0
613+
integer :: lev
604614
character (len = 31) :: fn !TSMP-PDAF: function name for state vector outpu
605615
character (len = 31) :: fn2 !TSMP-PDAF: function name for state vector outpu
606616
character (len = 32) :: fn3 !TSMP-PDAF: function name for state vector outpu
@@ -782,8 +792,10 @@ subroutine update_clm(tstartcycle, mype) bind(C,name="update_clm")
782792
do p = clm_begp, clm_endp
783793
c = patch%column(p)
784794
t_skin(p) = clm_statevec(state_clm2pdaf_p(p,1))
785-
t_soisno(c,1) = clm_statevec(state_clm2pdaf_p(p,1) + clm_varsize)
786-
t_veg(p) = clm_statevec(state_clm2pdaf_p(p,1) + 2*clm_varsize)
795+
do lev=1,nlevgrnd
796+
t_soisno(c,lev) = clm_statevec(state_clm2pdaf_p(p,1) + (lev)*clm_varsize)
797+
end do
798+
t_veg(p) = clm_statevec(state_clm2pdaf_p(p,1) + (1+nlevgrnd)*clm_varsize)
787799
end do
788800
endif
789801

0 commit comments

Comments
 (0)