source: trunk/LMDZ.COMMON/libf/evolution/soil_evolution_mod.F90 @ 2825

Last change on this file since 2825 was 2794, checked in by llange, 3 years ago

MARS PEM:

  • Add a PEMETAT0 that read "startfi_pem.nc"
  • Add the soil in the model: soil temperature, thermal properties, ice table
  • Add a routine that compute CO2 + H2O adsorption
  • Minor corrections in PEM.F90

LL

File size: 4.2 KB
Line 
1      MODULE soil_evolution_mod
2
3
4      IMPLICIT NONE
5     
6      CONTAINS   
7
8
9     subroutine soil_pem_CN(ngrid,nsoil, &
10               therm_i, timestep,tsurf,tsurf_prev,tsoil)
11
12
13      use comsoil_h_PEM, only:   & fluxgeo,layer_PEM
14      use comsoil_h,only: volcapa
15      implicit none
16
17!-----------------------------------------------------------------------
18!  Author: LL
19!  Purpose: Compute soil temperature using an implict 1st order scheme
20
21!  Note: depths of layers and mid-layers, soil thermal inertia and
22!        heat capacity are commons in comsoil_PEM.h
23!        A convergence loop is added until equilibrium
24!-----------------------------------------------------------------------
25
26#include "dimensions.h"
27!#include "dimphys.h"
28
29!#include"comsoil.h"
30
31
32!-----------------------------------------------------------------------
33!  arguments
34!  ---------
35!  inputs:
36      integer,intent(in) :: ngrid       ! number of (horizontal) grid-points
37      integer,intent(in) :: nsoil       ! number of soil layers 
38      real,intent(in) :: therm_i(ngrid,nsoil) ! thermal inertia
39      real,intent(in) :: timestep           ! time step
40      real,intent(in) :: tsurf(ngrid)   ! surface temperature
41      real,intent(in) :: tsurf_prev(ngrid)   ! surface temperature at previous timestep
42
43! outputs:
44      real,intent(inout) :: tsoil(ngrid,nsoil) ! soil (mid-layer) temperature
45
46
47     
48! local variables:
49      integer ig,ik,isoil
50      REAL :: ice_inertia = 2120.
51      real :: alph_PEM(nsoil)
52      real :: beta_PEM(nsoil)
53      real :: rhoc(nsoil)
54      real :: volcapa_ice = 1.43e7
55
56      real :: k_soil(nsoil)
57      real :: d1(nsoil),d2(nsoil),d3(nsoil),re(nsoil)
58      real :: Tcol(nsoil)
59
60
61do ig = 1,ngrid
62
63   Tcol(:) = tsoil(ig,:)
64!0. Build thermal properties of the ground
65
66
67  do isoil = 1,nsoil
68    if(abs(ice_inertia-therm_i(ig,isoil)).lt.1e-6) then
69      rhoc(isoil) = volcapa_ice
70    else
71      rhoc(isoil) = volcapa
72    endif
73    k_soil(isoil) = therm_i(ig,isoil)*therm_i(ig,isoil)/rhoc(isoil)
74  enddo
75
76!1. Build Crank-Nicholson scheme for heat equation, following Schorgofer derivations
77!1.a all the points which are not impacted by boundary condition
78  do isoil = 2,nsoil-1
79   
80    alph_pem(isoil) = timestep/((rhoc(isoil) + rhoc(isoil+1))/2.)*k_soil(isoil+1)/((layer_PEM(isoil+1) - layer_PEM(isoil))*(layer_PEM(isoil+1) - layer_PEM(isoil-1)))
81    beta_pem(isoil) = timestep/((rhoc(isoil) + rhoc(isoil+1))/2.)*k_soil(isoil)/((layer_PEM(isoil) - layer_PEM(isoil-1))*(layer_PEM(isoil+1) - layer_PEM(isoil-1)))
82  enddo
83
84!1.b At the surface: T(z=0) = Tsurf.
85
86alph_pem(1) = k_soil(2)*timestep/(rhoc(1)*layer_PEM(2)*(layer_PEM(2)-layer_PEM(1)))
87beta_pem(1) = k_soil(1)*timestep/(rhoc(1)*layer_PEM(2)*layer_PEM(1))
88
89
90
91!1.c At the bottom: dT/dz = -Fgeo/k.
92beta_pem(nsoil) = timestep*k_soil(nsoil)/(2.*rhoc(nsoil)*(layer_PEM(isoil) - layer_PEM(isoil-1))**2)
93
94!1 .d set the tridiagonal with diagonals d1,d2,d3, and right element re:
95
96
97  d1(:) = -beta_pem(:)   
98  d2(:) = 1. + alph_pem(:) + beta_pem(:)
99  d3(:) = -alph_pem(:)   
100  d2(nsoil) = 1. + beta_pem(nsoil)
101 
102  re(1)= alph_pem(1)*Tcol(2) + (1.-alph_pem(1)-beta_pem(1))*Tcol(1) + beta_pem(1)*(tsurf(ig)+tsurf_prev(ig))
103  do isoil=2,nsoil-1
104     re(isoil) = beta_pem(isoil)*Tcol(isoil-1) + (1.-alph_pem(isoil)-beta_pem(isoil))*Tcol(isoil) + alph_pem(isoil)*Tcol(isoil+1)
105  enddo
106  re(nsoil) = beta_pem(nsoil)*Tcol(nsoil-1) + (1.-beta_pem(nsoil))*Tcol(nsoil) + &
107  timestep/rhoc(nsoil)*fluxgeo/(layer_PEM(nsoil)-layer_PEM(nsoil-1))
108
109
110
111!2. Update by tridiagonal inversion
112  call tridag(d1,d2,d3,re,Tcol,nsoil)
113 
114
115
116      tsoil(ig,:) = Tcol(:)
117  enddo
118
119
120      END SUBROUTINE soil_pem_CN
121
122!====================================
123
124      SUBROUTINE TRIDAG(A,B,C,R,U,N)
125! From Numerical Recipes in Fortran 77
126      INTEGER, INTENT(IN)::  N
127      INTEGER ,PARAMETER:: NMAX = 1000
128
129      REAL, INTENT(IN) :: A(N),B(N),C(N),R(N)
130      REAL,INTENT(INOUT) :: U(N)
131
132      INTEGER j
133      REAL bet,gam(NMAX)
134
135      IF(B(1).EQ.0.)PAUSE
136      BET=B(1)
137      U(1)=R(1)/BET
138      DO 11 J=2,N
139        GAM(J)=C(J-1)/BET
140        BET=B(J)-A(J)*GAM(J)
141        IF(BET.EQ.0.)PAUSE
142        U(J)=(R(J)-A(J)*U(J-1))/BET
14311    CONTINUE
144      DO 12 J=N-1,1,-1
145        U(J)=U(J)-GAM(J+1)*U(J+1)
14612    CONTINUE
147      RETURN
148      END  SUBROUTINE TRIDAG
149
150
151      END MODULE soil_evolution_mod
152
153
154     
Note: See TracBrowser for help on using the repository browser.