1 | MODULE lwmain_mod |
---|
2 | |
---|
3 | IMPLICIT NONE |
---|
4 | |
---|
5 | CONTAINS |
---|
6 | |
---|
7 | subroutine lwmain (ig0,icount,kdlon,kflev |
---|
8 | . ,dp,dt0,emis |
---|
9 | . ,plev,tlev,tlay,aerosol,coolrate |
---|
10 | . ,fluxground,fluxtop |
---|
11 | . ,netrad |
---|
12 | & ,QIRsQREF3d,omegaIR3d,gIR3d |
---|
13 | & ,co2ice) |
---|
14 | |
---|
15 | c---------------------------------------------------------------------- |
---|
16 | c LWMAIN organizes the LTE longwave calculations |
---|
17 | c for layer 1 to layer "nlaylte" (stored in "yomlw_h") |
---|
18 | c---------------------------------------------------------------------- |
---|
19 | |
---|
20 | use dimradmars_mod, only: ndlo2, nflev, nir, ndlon, nuco2 |
---|
21 | use dimradmars_mod, only: naerkind |
---|
22 | use yomlw_h, only: nlaylte, xi |
---|
23 | implicit none |
---|
24 | |
---|
25 | c declarations |
---|
26 | c ------------- |
---|
27 | include "callkeys.h" |
---|
28 | include "comg1d.h" |
---|
29 | |
---|
30 | c---------------------------------------------------------------------- |
---|
31 | c 0.1 arguments |
---|
32 | c --------- |
---|
33 | c inputs/outputs: |
---|
34 | c ------- |
---|
35 | integer, intent(in) :: ig0 |
---|
36 | integer, intent(in) :: icount |
---|
37 | integer, intent(in) :: kdlon ! part of ngrid |
---|
38 | integer, intent(in) :: kflev ! part of nlayer |
---|
39 | |
---|
40 | real, intent(in) :: dp (ndlo2,kflev) ! layer pressure thickness (Pa) |
---|
41 | real, intent(in) :: dt0 (ndlo2) ! surface temperature discontinuity (K) |
---|
42 | real, intent(in) :: emis (ndlo2) ! surface emissivity |
---|
43 | real, intent(in) :: plev (ndlo2,kflev+1) ! level pressure (Pa) |
---|
44 | real, intent(in) :: tlev (ndlo2,kflev+1) ! level temperature (K) |
---|
45 | real, intent(in) :: tlay (ndlo2,kflev) ! layer temperature (K) |
---|
46 | real, intent(in) :: aerosol(ndlo2,kflev,naerkind) ! aerosol extinction optical |
---|
47 | c depth at reference wavelength "longrefvis" set |
---|
48 | c in dimradmars_mod , in each layer, for one of |
---|
49 | c the "naerkind" kind of aerosol optical properties. |
---|
50 | |
---|
51 | |
---|
52 | c outputs: |
---|
53 | c -------- |
---|
54 | real, intent(out) :: coolrate(ndlo2,kflev) ! cooling rate (K/s) |
---|
55 | real, intent(out) :: fluxground(ndlo2) ! downward ground flux (W/m2) |
---|
56 | real, intent(out) :: fluxtop(ndlo2) ! outgoing upward flux (W/m2) ("OLR") |
---|
57 | real, intent(out) :: netrad (ndlo2,kflev) ! radiative budget (W/m2) |
---|
58 | c Aerosol optical properties |
---|
59 | real, intent(in) :: QIRsQREF3d(ndlo2,kflev,nir,naerkind) |
---|
60 | real, intent(in) :: omegaIR3d(ndlo2,kflev,nir,naerkind) |
---|
61 | real, intent(in) :: gIR3d(ndlo2,kflev,nir,naerkind) |
---|
62 | real, intent(in) :: co2ice(ndlo2) ! co2 ice surface layer (kg.m-2) |
---|
63 | c---------------------------------------------------------------------- |
---|
64 | c 0.2 local arrays |
---|
65 | c ------------ |
---|
66 | |
---|
67 | real aer_t (ndlon,nuco2,nflev+1) ! transmission (aer) |
---|
68 | real co2_u (ndlon,nuco2,nflev+1) ! absorber amounts (co2) |
---|
69 | real co2_up (ndlon,nuco2,nflev+1) ! idem scaled by the pressure (co2) |
---|
70 | |
---|
71 | real bsurf (ndlon,nir) ! surface spectral planck function |
---|
72 | real btop (ndlon,nir) ! top spectral planck function |
---|
73 | real blev (ndlon,nir,nflev+1) ! level spectral planck function |
---|
74 | real blay (ndlon,nir,nflev) ! layer spectral planck function |
---|
75 | real dblay (ndlon,nir,nflev) ! layer gradient spectral planck function |
---|
76 | real dbsublay (ndlon,nir,2*nflev) ! layer gradient spectral planck function |
---|
77 | ! in sub layers |
---|
78 | |
---|
79 | real tautotal(ndlon,nflev,nir) ! \ Total single scattering |
---|
80 | real omegtotal(ndlon,nflev,nir) ! > properties (Addition of the |
---|
81 | real gtotal(ndlon,nflev,nir) ! / NAERKIND aerosols prop.) |
---|
82 | |
---|
83 | real newcoolrate(ndlon,nflev) ! cooling rate (K/s) / with implicite scheme |
---|
84 | |
---|
85 | real emis_gaz(ndlo2) ! emissivity for gaz computations |
---|
86 | |
---|
87 | integer jk,jkk,ja,jl |
---|
88 | |
---|
89 | |
---|
90 | c---------------------------------------------------------------------- |
---|
91 | c 0.3 Initialisation |
---|
92 | c -------------- |
---|
93 | |
---|
94 | DO jl=1 , kdlon |
---|
95 | IF(co2ice(jl) .GT. 20.e-3) THEN |
---|
96 | emis_gaz(jl)=1. |
---|
97 | ELSE |
---|
98 | emis_gaz(jl)=emis(jl) |
---|
99 | ENDIF |
---|
100 | ENDDO |
---|
101 | |
---|
102 | c---------------------------------------------------------------------- |
---|
103 | c 1.0 planck function |
---|
104 | c --------------- |
---|
105 | |
---|
106 | call lwb ( kdlon, kflev, tlev, tlay, dt0 |
---|
107 | . , bsurf, btop, blay, blev, dblay, dbsublay) |
---|
108 | |
---|
109 | c---------------------------------------------------------------------- |
---|
110 | c 2.0 absorber amounts |
---|
111 | c ---------------- |
---|
112 | |
---|
113 | call lwu ( kdlon, kflev |
---|
114 | . , dp, plev, tlay, aerosol |
---|
115 | & , QIRsQREF3d,omegaIR3d,gIR3d |
---|
116 | . , aer_t, co2_u, co2_up |
---|
117 | . , tautotal,omegtotal,gtotal) |
---|
118 | |
---|
119 | c---------------------------------------------------------------------- |
---|
120 | c 3.0 transmission functions / exchange coefficiants |
---|
121 | c ---------------------------------------------- |
---|
122 | |
---|
123 | c distants |
---|
124 | c -------- |
---|
125 | if( mod(icount-1,ilwd).eq.0) then |
---|
126 | |
---|
127 | c print*, 'CALL of DISTANTS' |
---|
128 | call lwxd ( ig0, kdlon, kflev, emis_gaz |
---|
129 | . , aer_t, co2_u, co2_up) |
---|
130 | |
---|
131 | endif |
---|
132 | c neighbours |
---|
133 | c ---------- |
---|
134 | if( mod(icount-1,ilwn).eq.0) then |
---|
135 | |
---|
136 | c print*, 'CALL of NEIGHBOURS' |
---|
137 | call lwxn ( ig0, kdlon, kflev |
---|
138 | . , dp |
---|
139 | . , aer_t, co2_u, co2_up) |
---|
140 | |
---|
141 | endif |
---|
142 | c boundaries |
---|
143 | c ---------- |
---|
144 | if( mod(icount-1,ilwb).eq.0) then |
---|
145 | |
---|
146 | c print*, 'CALL of BOUNDARIES' |
---|
147 | call lwxb ( ig0, kdlon, kflev, emis_gaz |
---|
148 | . , aer_t, co2_u, co2_up) |
---|
149 | |
---|
150 | endif |
---|
151 | |
---|
152 | c---------------------------------------------------------------------- |
---|
153 | c 4.0 cooling rate |
---|
154 | c ------------ |
---|
155 | |
---|
156 | call lwflux ( ig0, kdlon, kflev, dp |
---|
157 | . , bsurf, btop, blev, blay, dbsublay |
---|
158 | . , tlay, tlev, dt0 ! pour sortie dans g2d uniquement |
---|
159 | . , emis |
---|
160 | . , tautotal,omegtotal,gtotal |
---|
161 | . , coolrate, fluxground, fluxtop |
---|
162 | . , netrad) |
---|
163 | |
---|
164 | c do jk = 1, nlaylte |
---|
165 | c print*,coolrate(1,jk) |
---|
166 | c enddo |
---|
167 | |
---|
168 | c do jkk = 0 , nlaylte+1 |
---|
169 | c do jk = 0 , nlaylte+1 |
---|
170 | c do ja = 1 , nuco2 |
---|
171 | c do jl = 1 , ngrid |
---|
172 | c if (xi (jl,ja,jk,jkk) .LT. 0 |
---|
173 | c . .OR. xi (jl,ja,jk,jkk) .GT. 1 ) then |
---|
174 | c print*,'xi bande',ja,jk,jkk,xi (jl,ja,jk,jkk) |
---|
175 | c endif |
---|
176 | c enddo |
---|
177 | c enddo |
---|
178 | c enddo |
---|
179 | c enddo |
---|
180 | |
---|
181 | c---------------------------------------------------------------------- |
---|
182 | c |
---|
183 | c 5. shema semi-implicite (lwi) |
---|
184 | c --------------------------- |
---|
185 | c |
---|
186 | c |
---|
187 | call lwi (ig0,kdlon,kflev,netrad,dblay,dp |
---|
188 | . , newcoolrate) |
---|
189 | c |
---|
190 | c Verif que (X sol,space) + somme(X i,sol) = 1 |
---|
191 | c |
---|
192 | do jkk = 1 , nlaylte |
---|
193 | do jl = 1 , kdlon |
---|
194 | c print*,'NEW et OLD coolrate :',jkk,newcoolrate(jl,jkk) |
---|
195 | c . ,coolrate(jl,jkk) |
---|
196 | coolrate(jl,jkk) = newcoolrate(jl,jkk) |
---|
197 | enddo |
---|
198 | enddo |
---|
199 | c |
---|
200 | c---------------------------------------------------------------------- |
---|
201 | |
---|
202 | END SUBROUTINE lwmain |
---|
203 | |
---|
204 | END MODULE lwmain_mod |
---|
205 | |
---|