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