source: trunk/LMDZ.MARS/libf/phymars/yomlw_h.F90 @ 1972

Last change on this file since 1972 was 1772, checked in by aslmd, 7 years ago

LMDZ.MARS yomlw 1. add the possibility to deallocate common modules with subroutine (as in previous commit) 2. moved xi initialization to yomlw

File size: 3.2 KB
Line 
1module yomlw_h
2! Coefficients for the longwave radiation subroutines
3use dimradmars_mod, only: nir, nabsmx, npademx
4implicit none
5
6  real,save :: at(2,2)
7  real,save :: bt(2,2)
8  real,save :: tref= 200.0 ! tref: temperature dependence of the absorption
9  real,save :: xp(6,nir)
10  real,save :: tstand = 200.0 ! reference temperature for the Planck function
11  real,save :: ga(npademx,nabsmx)
12  real,save :: gb(npademx,nabsmx)
13  real,save :: cst_voigt(2,nabsmx)
14  real,save :: gcp ! = g/cpp (set in callradite)
15
16! Number of layers on which LTE calculations (in lw and sw) are performed
17! (Computed in nlthermeq) :
18  integer,save :: nlaylte
19
20  real,save,allocatable :: xi(:,:,:,:)
21  real,save,allocatable :: xi_ground(:,:)
22  real,save,allocatable :: xi_emis(:,:,:)
23
24contains
25
26! Allocate array (subroutine ini_yomlw_h is called by phys_state_var_init)
27  subroutine ini_yomlw_h(ngrid)
28   
29    use dimradmars_mod, only: nuco2,nflev
30    implicit none
31   
32    integer,intent(in) :: ngrid ! number of atmospheric columns
33   
34    allocate(xi(ngrid,nuco2,0:nflev+1,0:nflev+1))
35    allocate(xi_ground(ngrid,nuco2))
36    allocate(xi_emis(ngrid,nuco2,nflev-1))
37   
38    xi (:,:,:,:)=0. ! initialisation previously done in lwmain at firstcall
39 
40    ! initialize xp(),at(),bt(),ga(),gb(),cst_voigt()
41
42    ! What follows was taken from old "sulw.F" routine (by Jean-Jacques
43    ! Morcrette, ECMWF, further simplified by F. Forget 01/2000),
44    ! adapted to F90
45   
46    ! ROOTS AND WEIGHTS FOR THE 2-POINT GAUSSIAN QUADRATURE
47!     DATA (RT1(IG1),IG1=1,2) / -0.577350269, +0.577350269 /
48!     DATA (WG1(IG1),IG1=1,2) /  1.0        ,  1.0         /
49    ! COEFFICIENTS OF THE POLYNOMIALS GIVING THE PLANCK FUNCTIONS
50    xp = reshape ( (/ &
51           0.63849788E+01, 0.30969419E+02, 0.44790835E+02, &
52           0.52651048E+01,-0.18799237E+02, 0.92836181E+01, &
53           0.26166790E+02, 0.12348011E+03, 0.17868306E+03, &
54           0.33657659E+02,-0.66869343E+02, 0.21017507E+02, &
55           0.11101254E+02, 0.86037325E+02, 0.25892695E+03, &
56           0.35582991E+03, 0.16958020E+03,-0.41311413E+02, &
57           0.47045285E+02, 0.12234377E+03, 0.61873275E+02, &
58           -0.31971883E+02, 0.59168472E+01, 0.91927407E+01 &
59                    /) , (/ 6,nir /) )
60
61    ! temperature dependency of absorber amounts:
62    at = reshape( (/ &
63           0.694E-03, 0.272E-02, 0.275E-02, 0.178E-01  &
64                   /) , (/ 2,2 /) )
65    bt = reshape( (/ &
66           0.328E-05, 0.298E-05,-0.705E-04,-0.163E-04  &
67                   /) , (/ 2,2 /) )
68
69    ga(1:4,1) = (/ 0.288231E-04, 0.170794E-01,-0.339714E-01, 0.000000E+00 /)
70    gb(1:4,1) = (/ 0.288231E-04, 0.145426E-01, 0.543812E+00, 0.100000E+01 /)
71   
72    ga(1:4,2) = (/ 0.289299E-01, 0.190634E+01, 0.384061E+01, 0.000000E+00 /)
73    gb(1:4,2) = (/ 0.289299E-01, 0.189485E+01, 0.600363E+01, 0.100000E+01 /)
74   
75    cst_voigt = reshape( (/ &
76                  0.500E-02, 0.100E-01, 0.150E-01, 0.100E+00  &
77                         /) , (/ 2,2 /) )
78 
79  end subroutine ini_yomlw_h
80   
81  subroutine end_yomlw_h
82
83    implicit none
84
85    if (allocated(xi)) deallocate(xi)
86    if (allocated(xi_ground)) deallocate(xi_ground)
87    if (allocated(xi_emis)) deallocate(xi_emis)
88
89  end subroutine end_yomlw_h
90 
91end module yomlw_h
Note: See TracBrowser for help on using the repository browser.