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

Last change on this file since 706 was 285, checked in by aslmd, 13 years ago

08/09/11 == AS

LMDZ.MARS + MESOSCALE.
---> Setting up a more realistic water ice source at the poles (notably outliers)

surfini.F ?
Main changes and bug fixes

  • reference to comcstfi.h was wrong. big problem because e.g. pi was not known.
  • commented about a problem to be fixed, due to surfini being called before initracer.
  • MESOSCALE: put the mesoscale north cap definition into a precompiling flag #MESOSCALE

for the moment: if [alb_mean_TES > 0.26 and lat > 70] then outliers
(previously done in meso_inc_caps.F)

inifis.F ?
Just changed a comment with wrong formatting

--> below, only MESOSCALE

soil.F ?
if somewhere IT > IT_outliers, then makes it = IT_outliers

physiq.F ?
meso_inc/meso_inc_caps.F ?
meso_inc/meso_inc_ini.F ?
meso_inc_caps no longer called. keep for reference for the moment.

meso_inc/meso_inc_var.F ?
deleted lines with *_lim variables, now useless

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