source: trunk/LMDZ.GENERIC/libf/phystd/soil.F @ 537

Last change on this file since 537 was 253, checked in by emillour, 13 years ago

Generic GCM

  • Massive update to version 0.7

EM+RW

File size: 6.1 KB
Line 
1      subroutine soil(ngrid,nsoil,firstcall,
2     &          therm_i,
3     &          timestep,tsurf,tsoil,
4     &          capcal,fluxgrd)
5      implicit none
6
7!-----------------------------------------------------------------------
8!  Author: Ehouarn Millour
9!
10!  Purpose: Compute soil temperature using an implict 1st order scheme
11
12!  Note: depths of layers and mid-layers, soil thermal inertia and
13!        heat capacity are commons in comsoil.h
14!-----------------------------------------------------------------------
15
16#include "dimensions.h"
17#include "dimphys.h"
18
19#include"comsoil.h"
20
21c-----------------------------------------------------------------------
22!  arguments
23!  ---------
24!  inputs:
25      integer ngrid     ! number of (horizontal) grid-points
26      integer nsoil     ! number of soil layers
27      logical firstcall ! identifier for initialization call
28      real therm_i(ngrid,nsoil) ! thermal inertia
29      real timestep         ! time step
30      real tsurf(ngrid)   ! surface temperature
31! outputs:
32      real tsoil(ngrid,nsoil) ! soil (mid-layer) temperature
33      real capcal(ngrid) ! surface specific heat
34      real fluxgrd(ngrid) ! surface diffusive heat flux
35
36! local saved variables:
37      !real,save :: mthermdiff(ngridmx,0:nsoilmx-1) ! mid-layer thermal diffusivity
38      !real,save :: thermdiff(ngridmx,nsoilmx-1) ! inter-layer thermal diffusivity
39      !real,save :: coefq(0:nsoilmx-1)          ! q_{k+1/2} coefficients
40      !real,save :: coefd(ngridmx,nsoilmx-1)    ! d_k coefficients
41      !real,save :: alph(ngridmx,nsoilmx-1)     ! alpha_k coefficients
42      !real,save :: beta(ngridmx,nsoilmx-1)     ! beta_k coefficients
43      !real,save :: mu
44
45      real mthermdiff(ngridmx,0:nsoilmx-1) ! mid-layer thermal diffusivity
46      real thermdiff(ngridmx,nsoilmx-1)    ! inter-layer thermal diffusivity
47      real coefq(0:nsoilmx-1)              ! q_{k+1/2} coefficients
48      real coefd(ngridmx,nsoilmx-1)        ! d_k coefficients
49      real alph(ngridmx,nsoilmx-1)         ! alpha_k coefficients
50      real beta(ngridmx,nsoilmx-1)         ! beta_k coefficients
51      real mu
52           
53      save mthermdiff,thermdiff,coefq,coefd,alph,beta,mu
54
55! local variables:
56      integer ig,ik
57
58
59!      print*,'tsoil=',tsoil
60!      print*,'tsurf=',tsurf
61
62! 0. Initialisations and preprocessing step
63      if (firstcall) then
64      ! note: firstcall is set to .true. or .false. by the caller
65      !       and not changed by soil.F
66! 0.1 Build mthermdiff(:), the mid-layer thermal diffusivities
67      do ig=1,ngrid
68        do ik=0,nsoil-1
69          mthermdiff(ig,ik)=therm_i(ig,ik+1)*therm_i(ig,ik+1)/volcapa
70!         write(*,*),'soil: ik: ',ik,' mthermdiff:',mthermdiff(ig,ik)
71        enddo
72      enddo
73
74! 0.2 Build thermdiff(:), the "interlayer" thermal diffusivities
75      do ig=1,ngrid
76        do ik=1,nsoil-1
77      thermdiff(ig,ik)=((layer(ik)-mlayer(ik-1))*mthermdiff(ig,ik)
78     &                +(mlayer(ik)-layer(ik))*mthermdiff(ig,ik-1))
79     &                    /(mlayer(ik)-mlayer(ik-1))
80!       write(*,*),'soil: ik: ',ik,' thermdiff:',thermdiff(ig,ik)
81        enddo
82      enddo
83
84! 0.3 Build coefficients mu, q_{k+1/2}, d_k, alpha_k and capcal
85      ! mu
86      mu=mlayer(0)/(mlayer(1)-mlayer(0))
87
88      ! q_{1/2}
89      coefq(0)=volcapa*layer(1)/timestep
90        ! q_{k+1/2}
91        do ik=1,nsoil-1
92          coefq(ik)=volcapa*(layer(ik+1)-layer(ik))
93     &                 /timestep
94        enddo
95
96      do ig=1,ngrid
97        ! d_k
98        do ik=1,nsoil-1
99          coefd(ig,ik)=thermdiff(ig,ik)/(mlayer(ik)-mlayer(ik-1))
100        enddo
101       
102        ! alph_{N-1}
103        alph(ig,nsoil-1)=coefd(ig,nsoil-1)/
104     &                  (coefq(nsoil-1)+coefd(ig,nsoil-1))
105        ! alph_k
106        do ik=nsoil-2,1,-1
107          alph(ig,ik)=coefd(ig,ik)/(coefq(ik)+coefd(ig,ik+1)*
108     &                              (1.-alph(ig,ik+1))+coefd(ig,ik))
109        enddo
110
111        ! capcal
112! Cstar
113        capcal(ig)=volcapa*layer(1)+
114     &              (thermdiff(ig,1)/(mlayer(1)-mlayer(0)))*
115     &              (timestep*(1.-alph(ig,1)))
116! Cs
117        capcal(ig)=capcal(ig)/(1.+mu*(1.0-alph(ig,1))*
118     &                         thermdiff(ig,1)/mthermdiff(ig,0))
119      !write(*,*)'soil: ig=',ig,' capcal(ig)=',capcal(ig)
120      enddo ! of do ig=1,ngrid
121           
122      else ! of if (firstcall)
123
124
125!  1. Compute soil temperatures
126! First layer:
127      do ig=1,ngrid
128        tsoil(ig,1)=(tsurf(ig)+mu*beta(ig,1)*
129     &                         thermdiff(ig,1)/mthermdiff(ig,0))/
130     &              (1.+mu*(1.0-alph(ig,1))*
131     &               thermdiff(ig,1)/mthermdiff(ig,0))
132      enddo
133! Other layers:
134      do ik=1,nsoil-1
135        do ig=1,ngrid
136          tsoil(ig,ik+1)=alph(ig,ik)*tsoil(ig,ik)+beta(ig,ik)
137        enddo
138      enddo
139     
140      endif! of if (firstcall)
141
142!  2. Compute beta coefficients (preprocessing for next time step)
143! Bottom layer, beta_{N-1}
144      do ig=1,ngrid
145        beta(ig,nsoil-1)=coefq(nsoil-1)*tsoil(ig,nsoil)
146     &                   /(coefq(nsoil-1)+coefd(ig,nsoil-1))
147      enddo
148! Other layers
149      do ik=nsoil-2,1,-1
150        do ig=1,ngrid
151          beta(ig,ik)=(coefq(ik)*tsoil(ig,ik+1)+
152     &                 coefd(ig,ik+1)*beta(ig,ik+1))/
153     &                 (coefq(ik)+coefd(ig,ik+1)*(1.0-alph(ig,ik+1))
154     &                  +coefd(ig,ik))
155        enddo
156      enddo
157
158
159!  3. Compute surface diffusive flux & calorific capacity
160      do ig=1,ngrid
161! Cstar
162!        capcal(ig)=volcapa(ig,1)*layer(ig,1)+
163!     &              (thermdiff(ig,1)/(mlayer(ig,1)-mlayer(ig,0)))*
164!     &              (timestep*(1.-alph(ig,1)))
165! Fstar
166
167!         print*,'this far in soil 1'
168!         print*,'thermdiff=',thermdiff(ig,1)
169!         print*,'mlayer=',mlayer
170!         print*,'beta=',beta(ig,1)
171!         print*,'alph=',alph(ig,1)
172!         print*,'tsoil=',tsoil(ig,1)
173
174        fluxgrd(ig)=(thermdiff(ig,1)/(mlayer(1)-mlayer(0)))*
175     &              (beta(ig,1)+(alph(ig,1)-1.0)*tsoil(ig,1))
176
177!        mu=mlayer(ig,0)/(mlayer(ig,1)-mlayer(ig,0))
178!        capcal(ig)=capcal(ig)/(1.+mu*(1.0-alph(ig,1))*
179!     &                         thermdiff(ig,1)/mthermdiff(ig,0))
180! Fs
181        fluxgrd(ig)=fluxgrd(ig)+(capcal(ig)/timestep)*
182     &              (tsoil(ig,1)*(1.+mu*(1.0-alph(ig,1))*
183     &                         thermdiff(ig,1)/mthermdiff(ig,0))
184     &               -tsurf(ig)-mu*beta(ig,1)*
185     &                          thermdiff(ig,1)/mthermdiff(ig,0))
186      enddo
187
188
189      end
190
Note: See TracBrowser for help on using the repository browser.