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

Last change on this file since 1103 was 1047, checked in by emillour, 12 years ago

Mars GCM:

  • IMPORTANT CHANGE: Removed all reference/use of ngridmx (dimphys.h) in routines (necessary prerequisite to using parallel dynamics); in most cases this just means adding 'ngrid' as routine argument, and making local saved variables allocatable (and allocated at first call). In the process, had to convert many *.h files to equivalent modules: yomaer.h => yomaer_h.F90 , surfdat.h => surfdat_h.F90 , comsaison.h => comsaison_h.F90 , yomlw.h => yomlw_h.F90 , comdiurn.h => comdiurn_h.F90 , dimradmars.h => dimradmars_mod.F90 , comgeomfi.h => comgeomfi_h.F90, comsoil.h => comsoil_h.F90 , slope.h => slope_mod.F90
  • Also updated EOF routines, everything is now in eofdump_mod.F90
  • Removed unused routine lectfux.F (in dyn3d)

EM

File size: 2.9 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 inifis)
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    ! initialize xp(),at(),bt(),ga(),gb(),cst_voigt()
39
40    ! What follows was taken from old "sulw.F" routine (by Jean-Jacques
41    ! Morcrette, ECMWF, further simplified by F. Forget 01/2000),
42    ! adapted to F90
43   
44    ! ROOTS AND WEIGHTS FOR THE 2-POINT GAUSSIAN QUADRATURE
45!     DATA (RT1(IG1),IG1=1,2) / -0.577350269, +0.577350269 /
46!     DATA (WG1(IG1),IG1=1,2) /  1.0        ,  1.0         /
47    ! COEFFICIENTS OF THE POLYNOMIALS GIVING THE PLANCK FUNCTIONS
48    xp = reshape ( (/ &
49           0.63849788E+01, 0.30969419E+02, 0.44790835E+02, &
50           0.52651048E+01,-0.18799237E+02, 0.92836181E+01, &
51           0.26166790E+02, 0.12348011E+03, 0.17868306E+03, &
52           0.33657659E+02,-0.66869343E+02, 0.21017507E+02, &
53           0.11101254E+02, 0.86037325E+02, 0.25892695E+03, &
54           0.35582991E+03, 0.16958020E+03,-0.41311413E+02, &
55           0.47045285E+02, 0.12234377E+03, 0.61873275E+02, &
56           -0.31971883E+02, 0.59168472E+01, 0.91927407E+01 &
57                    /) , (/ 6,nir /) )
58
59    ! temperature dependency of absorber amounts:
60    at = reshape( (/ &
61           0.694E-03, 0.272E-02, 0.275E-02, 0.178E-01  &
62                   /) , (/ 2,2 /) )
63    bt = reshape( (/ &
64           0.328E-05, 0.298E-05,-0.705E-04,-0.163E-04  &
65                   /) , (/ 2,2 /) )
66
67    ga(1:4,1) = (/ 0.288231E-04, 0.170794E-01,-0.339714E-01, 0.000000E+00 /)
68    gb(1:4,1) = (/ 0.288231E-04, 0.145426E-01, 0.543812E+00, 0.100000E+01 /)
69   
70    ga(1:4,2) = (/ 0.289299E-01, 0.190634E+01, 0.384061E+01, 0.000000E+00 /)
71    gb(1:4,2) = (/ 0.289299E-01, 0.189485E+01, 0.600363E+01, 0.100000E+01 /)
72   
73    cst_voigt = reshape( (/ &
74                  0.500E-02, 0.100E-01, 0.150E-01, 0.100E+00  &
75                         /) , (/ 2,2 /) )
76 
77  end subroutine ini_yomlw_h
78   
79end module yomlw_h
Note: See TracBrowser for help on using the repository browser.