source: trunk/LMDZ.MARS/libf/phymars/soil.F @ 1112

Last change on this file since 1112 was 1047, checked in by emillour, 11 years ago

Mars GCM:

  • IMPORTANT CHANGE: Removed all reference/use of ngridmx (dimphys.h) in routines (necessary prerequisite to using parallel dynamics); in most cases this just means adding 'ngrid' as routine argument, and making local saved variables allocatable (and allocated at first call). In the process, had to convert many *.h files to equivalent modules: yomaer.h => yomaer_h.F90 , surfdat.h => surfdat_h.F90 , comsaison.h => comsaison_h.F90 , yomlw.h => yomlw_h.F90 , comdiurn.h => comdiurn_h.F90 , dimradmars.h => dimradmars_mod.F90 , comgeomfi.h => comgeomfi_h.F90, comsoil.h => comsoil_h.F90 , slope.h => slope_mod.F90
  • Also updated EOF routines, everything is now in eofdump_mod.F90
  • Removed unused routine lectfux.F (in dyn3d)

EM

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