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

Last change on this file since 2153 was 1268, checked in by aslmd, 11 years ago

LMDZ.MARS. related to r1266. forgot to remove a few now-obsolete dimensions.h includes in Mars physics.

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