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

Last change on this file since 1009 was 833, checked in by jbmadeleine, 12 years ago

Mars GCM:

Implemented the thermal inertia feedback:

  • Added a flag in callphys.def called tifeedback, set to false by default;
  • Changed physiq.F to call soil.F with a new thermal inertia if tifeedback = true;
  • Added a routine called soil_tifeedback.F that computes the new thermal inertia of the subsurface when ice is deposited;

    Modified files: soil.F, physiq.F, inifis.F, callkeys.h
    Added files: soil_tifeedback.F

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
[833]53      if (firstcall.or.tifeedback) then
[38]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           
[833]129      endif ! of if (firstcall.or.tifeedback)
[38]130
131!  1. Compute soil temperatures
[833]132      IF (.not.firstcall) THEN
[38]133! First layer:
134      do ig=1,ngrid
135        tsoil(ig,1)=(tsurf(ig)+mu*beta(ig,1)*
136     &                         thermdiff(ig,1)/mthermdiff(ig,0))/
137     &              (1.+mu*(1.0-alph(ig,1))*
138     &               thermdiff(ig,1)/mthermdiff(ig,0))
139      enddo
140! Other layers:
141      do ik=1,nsoil-1
142        do ig=1,ngrid
143          tsoil(ig,ik+1)=alph(ig,ik)*tsoil(ig,ik)+beta(ig,ik)
144        enddo
145      enddo
146     
[833]147      ENDIF! of if (.not.firstcall)
[38]148
149!  2. Compute beta coefficients (preprocessing for next time step)
150! Bottom layer, beta_{N-1}
151      do ig=1,ngrid
152        beta(ig,nsoil-1)=coefq(nsoil-1)*tsoil(ig,nsoil)
153     &                   /(coefq(nsoil-1)+coefd(ig,nsoil-1))
154      enddo
155! Other layers
156      do ik=nsoil-2,1,-1
157        do ig=1,ngrid
158          beta(ig,ik)=(coefq(ik)*tsoil(ig,ik+1)+
159     &                 coefd(ig,ik+1)*beta(ig,ik+1))/
160     &                 (coefq(ik)+coefd(ig,ik+1)*(1.0-alph(ig,ik+1))
161     &                  +coefd(ig,ik))
162        enddo
163      enddo
164
165!  3. Compute surface diffusive flux & calorific capacity
166      do ig=1,ngrid
167! Cstar
168!        capcal(ig)=volcapa(ig,1)*layer(ig,1)+
169!     &              (thermdiff(ig,1)/(mlayer(ig,1)-mlayer(ig,0)))*
170!     &              (timestep*(1.-alph(ig,1)))
171! Fstar
172        fluxgrd(ig)=(thermdiff(ig,1)/(mlayer(1)-mlayer(0)))*
173     &              (beta(ig,1)+(alph(ig,1)-1.0)*tsoil(ig,1))
174
175!        mu=mlayer(ig,0)/(mlayer(ig,1)-mlayer(ig,0))
176!        capcal(ig)=capcal(ig)/(1.+mu*(1.0-alph(ig,1))*
177!     &                         thermdiff(ig,1)/mthermdiff(ig,0))
178! Fs
179        fluxgrd(ig)=fluxgrd(ig)+(capcal(ig)/timestep)*
180     &              (tsoil(ig,1)*(1.+mu*(1.0-alph(ig,1))*
181     &                         thermdiff(ig,1)/mthermdiff(ig,0))
182     &               -tsurf(ig)-mu*beta(ig,1)*
183     &                          thermdiff(ig,1)/mthermdiff(ig,0))
184      enddo
185
186      end
187
Note: See TracBrowser for help on using the repository browser.