source: trunk/LMDZ.COMMON/libf/evolution/soil_pem_compute_mod.F90 @ 3116

Last change on this file since 3116 was 3076, checked in by jbclement, 21 months ago

PEM:
Big cleaning/improvements of the PEM:

  • Conversion of "abort_pem.F" and "soil_settings_PEM.F" into Fortran 90;
  • Transformation of every PEM subroutines into module;
  • Rewriting of many subroutines with modern Fortran syntax;
  • Correction of a bug in "pem.F90" when calling 'recomp_tend_co2_slope'. The arguments were given in disorder and emissivity was missing;
  • Update of "launch_pem.sh" in deftank.

JBC

File size: 4.2 KB
RevLine 
[3076]1MODULE soil_pem_compute_mod
[2794]2
[3076]3implicit none
[2794]4
[3076]5!=======================================================================
6contains
7!=======================================================================
[2794]8
[3076]9SUBROUTINE soil_pem_compute(ngrid,nsoil,firstcall,therm_i,timestep,tsurf,tsoil)
10
11use comsoil_h_PEM, only: layer_PEM, mlayer_PEM, mthermdiff_PEM, thermdiff_PEM, coefq_PEM, coefd_PEM, mu_PEM, alph_PEM, beta_PEM, fluxgeo
12use comsoil_h,     only: volcapa
13
14implicit none
15
[2794]16!-----------------------------------------------------------------------
17!  Author: LL
18!  Purpose: Compute soil temperature using an implict 1st order scheme
19
20!  Note: depths of layers and mid-layers, soil thermal inertia and
21!        heat capacity are commons in comsoil_PEM.h
22!-----------------------------------------------------------------------
23
24#include "dimensions.h"
25
26!-----------------------------------------------------------------------
27!  arguments
28!  ---------
29!  inputs:
[3076]30integer,                      intent(in) :: ngrid     ! number of (horizontal) grid-points
31integer,                      intent(in) :: nsoil     ! number of soil layers 
32logical,                      intent(in) :: firstcall ! identifier for initialization call
33real, dimension(ngrid,nsoil), intent(in) :: therm_i   ! thermal inertia [SI]
34real,                         intent(in) :: timestep  ! time step [s]
35real, dimension(ngrid),       intent(in) :: tsurf     ! surface temperature [K]
[2794]36 
37! outputs:
[3076]38real, dimension(ngrid,nsoil), intent(inout) :: tsoil ! soil (mid-layer) temperature [K]
[2794]39! local variables:
[3076]40integer :: ig, ik   
[2794]41
42! 0. Initialisations and preprocessing step
43 if (firstcall) then
44! 0.1 Build mthermdiff_PEM(:), the mid-layer thermal diffusivities
[3076]45    do ig = 1,ngrid
46        do ik = 0,nsoil - 1
47            mthermdiff_PEM(ig,ik) = therm_i(ig,ik + 1)*therm_i(ig,ik + 1)/volcapa   
48        enddo
49    enddo
[2794]50
51! 0.2 Build thermdiff(:), the "interlayer" thermal diffusivities
[3076]52    do ig = 1,ngrid
53        do ik = 1,nsoil - 1
54            thermdiff_PEM(ig,ik) = ((layer_PEM(ik) - mlayer_PEM(ik - 1))*mthermdiff_PEM(ig,ik) &
55                                   + (mlayer_PEM(ik) - layer_PEM(ik))*mthermdiff_PEM(ig,ik - 1))/(mlayer_PEM(ik) - mlayer_PEM(ik - 1))
56        enddo
57    enddo
[2794]58
59! 0.3 Build coefficients mu_PEM, q_{k+1/2}, d_k, alph_PEMa_k and capcal
[3076]60    ! mu_PEM
61    mu_PEM = mlayer_PEM(0)/(mlayer_PEM(1) - mlayer_PEM(0))
[2794]62
[3076]63    ! q_{1/2}
64    coefq_PEM(0) = volcapa*layer_PEM(1)/timestep
65    ! q_{k+1/2}
66    do ik = 1,nsoil - 1
67        coefq_PEM(ik) = volcapa*(layer_PEM(ik + 1) - layer_PEM(ik))/timestep
68    enddo
[2794]69
[3076]70    do ig = 1,ngrid
71        ! d_k
72        do ik = 1,nsoil - 1
73            coefd_PEM(ig,ik) = thermdiff_PEM(ig,ik)/(mlayer_PEM(ik)-mlayer_PEM(ik - 1))
74        enddo
75
76        ! alph_PEM_{N-1}
77        alph_PEM(ig,nsoil - 1) = coefd_PEM(ig,nsoil-1)/(coefq_PEM(nsoil - 1) + coefd_PEM(ig,nsoil - 1))
[2794]78        ! alph_PEM_k
[3076]79        do ik = nsoil - 2,1,-1
80            alph_PEM(ig,ik) = coefd_PEM(ig,ik)/(coefq_PEM(ik) + coefd_PEM(ig,ik + 1)*(1. - alph_PEM(ig,ik + 1)) + coefd_PEM(ig,ik))
81        enddo
[2794]82      enddo ! of do ig=1,ngrid
83endif ! of if (firstcall)
84
[3076]85IF (.not. firstcall) THEN
86! 2. Compute soil temperatures
[2794]87! First layer:
[3076]88    do ig = 1,ngrid
89        tsoil(ig,1) = (tsurf(ig) + mu_PEM*beta_PEM(ig,1)*thermdiff_PEM(ig,1)/mthermdiff_PEM(ig,0))/ &
90                      (1. + mu_PEM*(1. - alph_PEM(ig,1))*thermdiff_PEM(ig,1)/mthermdiff_PEM(ig,0))
[2794]91
92! Other layers:
[3076]93        do ik = 1,nsoil - 1
94                tsoil(ig,ik + 1) = alph_PEM(ig,ik)*tsoil(ig,ik) + beta_PEM(ig,ik)
[2835]95        enddo
[3076]96    enddo
97endif
[2794]98
99!  2. Compute beta_PEM coefficients (preprocessing for next time step)
100! Bottom layer, beta_PEM_{N-1}
[3076]101do ig = 1,ngrid
102    beta_PEM(ig,nsoil - 1) = coefq_PEM(nsoil - 1)*tsoil(ig,nsoil)/(coefq_PEM(nsoil - 1) + coefd_PEM(ig,nsoil - 1)) &
103                             + fluxgeo/(coefq_PEM(nsoil - 1) + coefd_PEM(ig,nsoil - 1))
104enddo
[2794]105! Other layers
[3076]106do ik = nsoil-2,1,-1
107    do ig = 1,ngrid
108        beta_PEM(ig,ik) = (coefq_PEM(ik)*tsoil(ig,ik + 1) + coefd_PEM(ig,ik + 1)*beta_PEM(ig,ik + 1))/ &
109                          (coefq_PEM(ik) + coefd_PEM(ig,ik + 1)*(1. - alph_PEM(ig,ik + 1)) + coefd_PEM(ig,ik))
110    enddo
111enddo
[2794]112
[3076]113END SUBROUTINE soil_pem_compute
[2794]114
[3076]115END MODULE soil_pem_compute_mod
Note: See TracBrowser for help on using the repository browser.