1 | subroutine callcorrk(icount,ngrid,nlayer,pq,nq,qsurf, |
---|
2 | & albedo,emis,mu0,pplev,pplay,pt, |
---|
3 | & tsurf,fract,dist_star,igout,aerosol,cpp3D, |
---|
4 | & dtlw,dtsw,fluxsurf_lw, |
---|
5 | & fluxsurf_sw,fluxtop_lw,fluxtop_sw,fluxtop_dn, |
---|
6 | & reffrad,tau_col,ptime,pday,firstcall,lastcall) |
---|
7 | |
---|
8 | use radinc_h |
---|
9 | use radcommon_h |
---|
10 | use ioipsl_getincom |
---|
11 | |
---|
12 | implicit none |
---|
13 | |
---|
14 | !================================================================== |
---|
15 | ! |
---|
16 | ! Purpose |
---|
17 | ! ------- |
---|
18 | ! Solve the radiative transfer using the correlated-k method for |
---|
19 | ! the gaseous absorption and the Toon et al. (1989) method for |
---|
20 | ! scatttering due to aerosols. |
---|
21 | ! |
---|
22 | ! Authors |
---|
23 | ! ------- |
---|
24 | ! Emmanuel 01/2001, Forget 09/2001 |
---|
25 | ! Robin Wordsworth (2009) |
---|
26 | ! |
---|
27 | !================================================================== |
---|
28 | |
---|
29 | #include "dimphys.h" |
---|
30 | #include "datafile.h" |
---|
31 | #include "comcstfi.h" |
---|
32 | #include "callkeys.h" |
---|
33 | #include "tracer.h" |
---|
34 | |
---|
35 | !----------------------------------------------------------------------- |
---|
36 | ! Declaration of the arguments (INPUT - OUTPUT) on the LMD GCM grid |
---|
37 | ! Layer #1 is the layer near the ground. |
---|
38 | ! Layer #nlayermx is the layer at the top. |
---|
39 | |
---|
40 | ! INPUT |
---|
41 | INTEGER icount |
---|
42 | INTEGER ngrid,nlayer |
---|
43 | INTEGER igout |
---|
44 | REAL aerosol(ngrid,nlayermx,naerkind) ! aerosol tau (kg/kg) |
---|
45 | REAL albedo(ngrid) ! SW albedo |
---|
46 | REAL emis(ngrid) ! LW emissivity |
---|
47 | REAL pplay(ngrid,nlayermx) ! pres. level in GCM mid of layer |
---|
48 | REAL pplev(ngrid,nlayermx+1) ! pres. level at GCM layer boundaries |
---|
49 | |
---|
50 | REAL pt(ngrid,nlayermx) ! air temperature (K) |
---|
51 | REAL tsurf(ngrid) ! surface temperature (K) |
---|
52 | REAL dist_star,mu0(ngrid) ! distance star-planet (AU) |
---|
53 | REAL fract(ngrid) ! fraction of day |
---|
54 | |
---|
55 | ! Globally varying aerosol optical properties on GCM grid |
---|
56 | ! Not needed everywhere so not in radcommon_h |
---|
57 | REAL :: QVISsQREF3d(ngridmx,nlayermx,L_NSPECTV,naerkind) |
---|
58 | REAL :: omegaVIS3d(ngridmx,nlayermx,L_NSPECTV,naerkind) |
---|
59 | REAL :: gVIS3d(ngridmx,nlayermx,L_NSPECTV,naerkind) |
---|
60 | |
---|
61 | REAL :: QIRsQREF3d(ngridmx,nlayermx,L_NSPECTI,naerkind) |
---|
62 | REAL :: omegaIR3d(ngridmx,nlayermx,L_NSPECTI,naerkind) |
---|
63 | REAL :: gIR3d(ngridmx,nlayermx,L_NSPECTI,naerkind) |
---|
64 | |
---|
65 | REAL :: QREFvis3d(ngridmx,nlayermx,naerkind) |
---|
66 | REAL :: QREFir3d(ngridmx,nlayermx,naerkind) |
---|
67 | |
---|
68 | ! REAL :: omegaREFvis3d(ngridmx,nlayermx,naerkind) |
---|
69 | ! REAL :: omegaREFir3d(ngridmx,nlayermx,naerkind) ! not sure of the point of these... |
---|
70 | |
---|
71 | REAL reffrad(ngrid,nlayer,naerkind) |
---|
72 | REAL nueffrad(ngrid,nlayer,naerkind) |
---|
73 | |
---|
74 | ! OUTPUT |
---|
75 | REAL dtsw(ngridmx,nlayermx) ! heating rate (K/s) due to SW |
---|
76 | REAL dtlw(ngridmx,nlayermx) ! heating rate (K/s) due to LW |
---|
77 | REAL fluxsurf_lw(ngridmx) ! incident LW flux to surf (W/m2) |
---|
78 | REAL fluxtop_lw(ngridmx) ! outgoing LW flux to space (W/m2) |
---|
79 | REAL fluxsurf_sw(ngridmx) ! incident SW flux to surf (W/m2) |
---|
80 | REAL fluxtop_sw(ngridmx) ! outgoing LW flux to space (W/m2) |
---|
81 | REAL fluxtop_dn(ngridmx) ! incident top of atmosphere SW flux (W/m2) |
---|
82 | |
---|
83 | !----------------------------------------------------------------------- |
---|
84 | ! Declaration of the variables required by correlated-k subroutines |
---|
85 | ! Numbered from top to bottom unlike in the GCM! |
---|
86 | |
---|
87 | REAL*8 tmid(L_LEVELS),pmid(L_LEVELS) |
---|
88 | REAL*8 tlevrad(L_LEVELS),plevrad(L_LEVELS) |
---|
89 | |
---|
90 | ! Optical values for the optci/cv subroutines |
---|
91 | REAL*8 stel(L_NSPECTV),stel_fract(L_NSPECTV) |
---|
92 | REAL*8 dtaui(L_NLAYRAD,L_NSPECTI,L_NGAUSS) |
---|
93 | REAL*8 dtauv(L_NLAYRAD,L_NSPECTV,L_NGAUSS) |
---|
94 | REAL*8 cosbv(L_NLAYRAD,L_NSPECTV,L_NGAUSS) |
---|
95 | REAL*8 cosbi(L_NLAYRAD,L_NSPECTI,L_NGAUSS) |
---|
96 | REAL*8 wbari(L_NLAYRAD,L_NSPECTI,L_NGAUSS) |
---|
97 | REAL*8 wbarv(L_NLAYRAD,L_NSPECTV,L_NGAUSS) |
---|
98 | REAL*8 tauv(L_NLEVRAD,L_NSPECTV,L_NGAUSS) |
---|
99 | REAL*8 taucumv(L_LEVELS,L_NSPECTV,L_NGAUSS) |
---|
100 | REAL*8 taucumi(L_LEVELS,L_NSPECTI,L_NGAUSS) |
---|
101 | |
---|
102 | REAL*8 tauaero(L_LEVELS+1,naerkind) |
---|
103 | REAL*8 nfluxtopv,nfluxtopi,nfluxtop |
---|
104 | real*8 NFLUXTOPV_nu(L_NSPECTV) |
---|
105 | real*8 NFLUXTOPI_nu(L_NSPECTI) |
---|
106 | real*8 fluxupi_nu(L_NLAYRAD,L_NSPECTI) ! for 1D diagnostic |
---|
107 | REAL*8 fmneti(L_NLAYRAD),fmnetv(L_NLAYRAD) |
---|
108 | REAL*8 fluxupv(L_NLAYRAD),fluxupi(L_NLAYRAD) |
---|
109 | REAL*8 fluxdnv(L_NLAYRAD),fluxdni(L_NLAYRAD) |
---|
110 | REAL*8 albi,albv,acosz |
---|
111 | |
---|
112 | INTEGER ig,l,k,nw,iaer,irad |
---|
113 | |
---|
114 | real fluxtoplanet |
---|
115 | real*8 taugsurf(L_NSPECTV,L_NGAUSS-1) |
---|
116 | real*8 taugsurfi(L_NSPECTI,L_NGAUSS-1) |
---|
117 | |
---|
118 | real*8 qvar(L_LEVELS) ! mixing ratio of variable component |
---|
119 | REAL pq(ngridmx,nlayer,nq) |
---|
120 | REAL qsurf(ngridmx,nqmx) ! tracer on surface (e.g. kg.m-2) |
---|
121 | integer nq |
---|
122 | |
---|
123 | ! Local aerosol optical properties for each column on RADIATIVE grid |
---|
124 | real*8 QXVAER(L_LEVELS+1,L_NSPECTV,naerkind) |
---|
125 | real*8 QSVAER(L_LEVELS+1,L_NSPECTV,naerkind) |
---|
126 | real*8 GVAER(L_LEVELS+1,L_NSPECTV,naerkind) |
---|
127 | real*8 QXIAER(L_LEVELS+1,L_NSPECTI,naerkind) |
---|
128 | real*8 QSIAER(L_LEVELS+1,L_NSPECTI,naerkind) |
---|
129 | real*8 GIAER(L_LEVELS+1,L_NSPECTI,naerkind) |
---|
130 | |
---|
131 | save qxvaer, qsvaer, gvaer |
---|
132 | save qxiaer, qsiaer, giaer |
---|
133 | save QREFvis3d, QREFir3d |
---|
134 | |
---|
135 | REAL tau_col(ngrid) ! diagnostic from aeropacity |
---|
136 | |
---|
137 | ! Misc. |
---|
138 | logical firstcall, lastcall, nantest |
---|
139 | real*8 tempv(L_NSPECTV) |
---|
140 | real*8 tempi(L_NSPECTI) |
---|
141 | real*8 temp,temp1,temp2,pweight |
---|
142 | character(len=10) :: tmp1 |
---|
143 | character(len=10) :: tmp2 |
---|
144 | |
---|
145 | ! for fixed water vapour profiles |
---|
146 | integer i_var |
---|
147 | real RH |
---|
148 | real*8 pq_temp(nlayer) |
---|
149 | real ptemp, Ttemp, qsat |
---|
150 | |
---|
151 | ! real(KIND=r8) :: pq_temp(nlayer) ! better F90 way.. DOESNT PORT TO F77!!! |
---|
152 | |
---|
153 | ! for OLR spec |
---|
154 | integer OLRcount |
---|
155 | save OLRcount |
---|
156 | integer OLRcount2 |
---|
157 | save OLRcount2 |
---|
158 | character(len=2) :: tempOLR |
---|
159 | character(len=30) :: filenomOLR |
---|
160 | real ptime, pday |
---|
161 | logical OLRz |
---|
162 | real OLR_nu(ngrid,L_NSPECTI) |
---|
163 | |
---|
164 | ! Allow variations in cp with location |
---|
165 | REAL cpp3D(ngridmx,nlayermx) ! specific heat capacity at const. pressure |
---|
166 | |
---|
167 | ! for Dave Crisp LBL data comparison |
---|
168 | logical crisp |
---|
169 | crisp=.false. |
---|
170 | |
---|
171 | |
---|
172 | !======================================================================= |
---|
173 | ! Initialization on first call |
---|
174 | |
---|
175 | CALL zerophys((L_LEVELS+1)*L_NSPECTV*naerkind,qxvaer) |
---|
176 | CALL zerophys((L_LEVELS+1)*L_NSPECTV*naerkind,qsvaer) |
---|
177 | CALL zerophys((L_LEVELS+1)*L_NSPECTV*naerkind,gvaer) |
---|
178 | |
---|
179 | CALL zerophys((L_LEVELS+1)*L_NSPECTI*naerkind,qxiaer) |
---|
180 | CALL zerophys((L_LEVELS+1)*L_NSPECTI*naerkind,qsiaer) |
---|
181 | CALL zerophys((L_LEVELS+1)*L_NSPECTI*naerkind,giaer) |
---|
182 | |
---|
183 | if(firstcall) then |
---|
184 | |
---|
185 | |
---|
186 | !-------------------------------------------------- |
---|
187 | ! Effective radius and variance of the aerosols |
---|
188 | |
---|
189 | ! CO2 ice: |
---|
190 | DO l=1,nlayer |
---|
191 | DO ig=1,ngrid |
---|
192 | reffrad(ig,l,1) = 1.e-4 |
---|
193 | nueffrad(ig,l,1) = 0.1 |
---|
194 | ! these values will change once the microphysics gets to work |
---|
195 | ! UNLESS tracer=.false., in which case we should be working with |
---|
196 | ! a fixed aerosol layer, and be able to define reffrad in a |
---|
197 | ! .def file. To be improved! |
---|
198 | ENDDO |
---|
199 | ENDDO |
---|
200 | |
---|
201 | ! H2O ice: |
---|
202 | if(naerkind.eq.2)then |
---|
203 | DO l=1,nlayer |
---|
204 | DO ig=1,ngrid |
---|
205 | reffrad(ig,l,naerkind) = 1.e-5 |
---|
206 | nueffrad(ig,l,naerkind) = 0.1 |
---|
207 | ENDDO |
---|
208 | ENDDO |
---|
209 | endif |
---|
210 | |
---|
211 | print*, "Correlated-k data base folder:",datafile |
---|
212 | call getin("corrkdir",corrkdir) |
---|
213 | print*, "corrkdir = ",corrkdir |
---|
214 | write( tmp1, '(i3)' ) L_NSPECTI |
---|
215 | write( tmp2, '(i3)' ) L_NSPECTV |
---|
216 | banddir=trim(adjustl(tmp1))//'x'//trim(adjustl(tmp2)) |
---|
217 | banddir=trim(adjustl(corrkdir))//'/'//trim(adjustl(banddir)) |
---|
218 | |
---|
219 | print*,'starting sugas' |
---|
220 | call sugas_corrk ! set up gaseous absorption properties |
---|
221 | print*,'starting setspi' |
---|
222 | call setspi ! basic infrared properties |
---|
223 | print*,'starting setspv' |
---|
224 | call setspv ! basic visible properties |
---|
225 | print*,'starting suaer_corrk' |
---|
226 | call suaer_corrk ! set up aerosol optical properties |
---|
227 | |
---|
228 | Cmk= 0.01 * 1.0 / (g * mugaz * 1.672621e-27) ! q_main=1.0 assumed |
---|
229 | |
---|
230 | if((igcm_h2o_vap.eq.0) .and. varactive)then |
---|
231 | print*,'varactive in callcorrk but no h2o_vap tracer.' |
---|
232 | stop |
---|
233 | endif |
---|
234 | |
---|
235 | firstcall=.false. |
---|
236 | |
---|
237 | end if |
---|
238 | |
---|
239 | !======================================================================= |
---|
240 | ! Initialization on every call |
---|
241 | |
---|
242 | do l=1,nlayer |
---|
243 | do ig=1,ngrid |
---|
244 | do iaer=1,naerkind |
---|
245 | nueffrad(ig,l,iaer) = 0.1 ! this little bastard stays at 0.1 forever |
---|
246 | enddo |
---|
247 | enddo |
---|
248 | enddo |
---|
249 | |
---|
250 | ! if(aerofixed)then ! fixed particle radii if aerofixed=.true. |
---|
251 | do l=1,nlayer |
---|
252 | do ig=1,ngrid |
---|
253 | reffrad(ig,l,1) = 3.5e-5 |
---|
254 | enddo |
---|
255 | enddo |
---|
256 | if(naerkind.eq.2)then |
---|
257 | do l=1,nlayer |
---|
258 | do ig=1,ngrid |
---|
259 | reffrad(ig,l,2) = 5.e-6 |
---|
260 | enddo |
---|
261 | enddo |
---|
262 | endif |
---|
263 | ! endif |
---|
264 | |
---|
265 | ! how much light we get |
---|
266 | fluxtoplanet=0 |
---|
267 | DO nw=1,L_NSPECTV |
---|
268 | stel(nw)=stellarf(nw)/(dist_star**2) |
---|
269 | fluxtoplanet=fluxtoplanet + stel(nw) |
---|
270 | END DO |
---|
271 | |
---|
272 | ! print*,'nueffrad',nueffrad |
---|
273 | ! print*,'reffrad before',reffrad |
---|
274 | ! print*,'nueffrad before',nueffrad |
---|
275 | ! CALL aeroptproperties(ngrid,nlayer,reffrad,nueffrad, |
---|
276 | ! & QVISsQREF3d,omegaVIS3d,gVIS3d, |
---|
277 | ! & QIRsQREF3d,omegaIR3d,gIR3d, |
---|
278 | ! & QREFvis3d,QREFir3d) ! get 3D aerosol optical properties |
---|
279 | |
---|
280 | CALL aeroptpropertiesJBnew(ngrid,nlayer,reffrad,nueffrad, |
---|
281 | & QVISsQREF3d,omegaVIS3d,gVIS3d, |
---|
282 | & QIRsQREF3d,omegaIR3d,gIR3d, |
---|
283 | & QREFvis3d,QREFir3d) |
---|
284 | |
---|
285 | ! print*,'QVISsQREF3d=',QVISsQREF3d |
---|
286 | ! print*,'QIRsQREF3d=',QIRsQREF3d |
---|
287 | ! print*,'omegaVIS3d=',omegaVIS3d |
---|
288 | ! print*,'omegaIR3d=',omegaIR3d |
---|
289 | ! print*,'gVIS3d=',gVIS3d |
---|
290 | ! print*,'gIR3d=',gIR3d |
---|
291 | ! print*,'reffrad after',reffrad |
---|
292 | ! print*,'nueffrad after',nueffrad |
---|
293 | ! print*,'QREFvis3d=',QREFvis3d |
---|
294 | ! print*,'QREFir3d=',QREFir3d |
---|
295 | ! stop |
---|
296 | |
---|
297 | call aeropacity(ngrid,nlayer,nq,pplev,pq,aerosol, |
---|
298 | & reffrad,QREFvis3d,QREFir3d,tau_col) ! get aerosol optical depths |
---|
299 | |
---|
300 | |
---|
301 | !----------------------------------------------------------------------- |
---|
302 | ! Starting Big Loop over every GCM column |
---|
303 | do ig=1,ngridmx |
---|
304 | |
---|
305 | !======================================================================= |
---|
306 | ! Transformation of the GCM variables |
---|
307 | |
---|
308 | !----------------------------------------------------------------------- |
---|
309 | ! Aerosol optical properties Qext, Qscat and g |
---|
310 | ! The transformation in the vertical is the same as for temperature |
---|
311 | |
---|
312 | ! shortwave |
---|
313 | do iaer=1,naerkind |
---|
314 | DO nw=1,L_NSPECTV |
---|
315 | do l=1,nlayermx |
---|
316 | |
---|
317 | temp1=QVISsQREF3d(ig,nlayermx+1-l,nw,iaer) |
---|
318 | $ *QREFvis3d(ig,nlayermx+1-l,iaer) |
---|
319 | |
---|
320 | temp2=QVISsQREF3d(ig,max(nlayermx-l,1),nw,iaer) |
---|
321 | $ *QREFvis3d(ig,max(nlayermx-l,1),iaer) |
---|
322 | |
---|
323 | qxvaer(2*l,nw,iaer) = temp1 |
---|
324 | qxvaer(2*l+1,nw,iaer)=(temp1+temp2)/2 |
---|
325 | |
---|
326 | temp1=temp1*omegavis3d(ig,nlayermx+1-l,nw,iaer) |
---|
327 | temp2=temp2*omegavis3d(ig,max(nlayermx-l,1),nw,iaer) |
---|
328 | |
---|
329 | qsvaer(2*l,nw,iaer) = temp1 |
---|
330 | qsvaer(2*l+1,nw,iaer)=(temp1+temp2)/2 |
---|
331 | |
---|
332 | temp1=gvis3d(ig,nlayermx+1-l,nw,iaer) |
---|
333 | temp2=gvis3d(ig,max(nlayermx-l,1),nw,iaer) |
---|
334 | |
---|
335 | gvaer(2*l,nw,iaer) = temp1 |
---|
336 | gvaer(2*l+1,nw,iaer)=(temp1+temp2)/2 |
---|
337 | |
---|
338 | end do |
---|
339 | qxvaer(1,nw,iaer)=qxvaer(2,nw,iaer) |
---|
340 | qxvaer(2*nlayermx+1,nw,iaer)=0. |
---|
341 | |
---|
342 | qsvaer(1,nw,iaer)=qsvaer(2,nw,iaer) |
---|
343 | qsvaer(2*nlayermx+1,nw,iaer)=0. |
---|
344 | |
---|
345 | gvaer(1,nw,iaer)=gvaer(2,nw,iaer) |
---|
346 | gvaer(2*nlayermx+1,nw,iaer)=0. |
---|
347 | end do |
---|
348 | |
---|
349 | ! longwave |
---|
350 | DO nw=1,L_NSPECTI |
---|
351 | do l=1,nlayermx |
---|
352 | |
---|
353 | temp1=QIRsQREF3d(ig,nlayermx+1-l,nw,iaer) |
---|
354 | $ *QREFir3d(ig,nlayermx+1-l,iaer) |
---|
355 | |
---|
356 | temp2=QIRsQREF3d(ig,max(nlayermx-l,1),nw,iaer) |
---|
357 | $ *QREFir3d(ig,max(nlayermx-l,1),iaer) |
---|
358 | |
---|
359 | qxiaer(2*l,nw,iaer) = temp1 |
---|
360 | qxiaer(2*l+1,nw,iaer)=(temp1+temp2)/2 |
---|
361 | |
---|
362 | temp1=temp1*omegair3d(ig,nlayermx+1-l,nw,iaer) |
---|
363 | temp2=temp2*omegair3d(ig,max(nlayermx-l,1),nw,iaer) |
---|
364 | |
---|
365 | qsiaer(2*l,nw,iaer) = temp1 |
---|
366 | qsiaer(2*l+1,nw,iaer)=(temp1+temp2)/2 |
---|
367 | |
---|
368 | temp1=gir3d(ig,nlayermx+1-l,nw,iaer) |
---|
369 | temp2=gir3d(ig,max(nlayermx-l,1),nw,iaer) |
---|
370 | |
---|
371 | giaer(2*l,nw,iaer) = temp1 |
---|
372 | giaer(2*l+1,nw,iaer)=(temp1+temp2)/2 |
---|
373 | |
---|
374 | end do |
---|
375 | |
---|
376 | qxiaer(1,nw,iaer)=qxiaer(2,nw,iaer) |
---|
377 | qxiaer(2*nlayermx+1,nw,iaer)=0. |
---|
378 | |
---|
379 | qsiaer(1,nw,iaer)=qsiaer(2,nw,iaer) |
---|
380 | qsiaer(2*nlayermx+1,nw,iaer)=0. |
---|
381 | |
---|
382 | giaer(1,nw,iaer)=giaer(2,nw,iaer) |
---|
383 | giaer(2*nlayermx+1,nw,iaer)=0. |
---|
384 | end do |
---|
385 | end do |
---|
386 | |
---|
387 | !----------------------------------------------------------------------- |
---|
388 | ! Aerosol optical depths |
---|
389 | |
---|
390 | do iaer=1,naerkind ! a bug was here |
---|
391 | do k=0,nlayer-1 |
---|
392 | |
---|
393 | pweight= |
---|
394 | $ (pplay(ig,L_NLAYRAD-k)-pplev(ig,L_NLAYRAD-k+1))/ |
---|
395 | $ (pplev(ig,L_NLAYRAD-k)-pplev(ig,L_NLAYRAD-k+1)) |
---|
396 | |
---|
397 | temp=aerosol(ig,L_NLAYRAD-k,iaer)/ |
---|
398 | $ QREFvis3d(ig,L_NLAYRAD-k,iaer) |
---|
399 | |
---|
400 | tauaero(2*k+2,iaer)=max(temp*pweight,0.0) |
---|
401 | tauaero(2*k+3,iaer)=max(temp-tauaero(2*k+2,iaer),0.0) |
---|
402 | ! tauaero(2*k+2,iaer)= |
---|
403 | ! & real(max(temp*pweight,0.0),8) |
---|
404 | ! tauaero(2*k+3,iaer)= |
---|
405 | ! & real(max(temp-tauaero(2*k+2,iaer),0.0),8) |
---|
406 | |
---|
407 | end do |
---|
408 | tauaero(1,iaer)=0. |
---|
409 | end do |
---|
410 | |
---|
411 | ! Albedo and emissivity |
---|
412 | albi=1-emis(ig) ! longwave |
---|
413 | albv=albedo(ig) ! shortwave |
---|
414 | |
---|
415 | if(ngrid.eq.1) then ! fixed zenith angle in 1D |
---|
416 | acosz = 0.5 |
---|
417 | print*,'Solar zenith angle fixed at 60 deg in 1D (callcorrk.F)' |
---|
418 | else |
---|
419 | acosz=mu0(ig) ! cosine of sun incident angle |
---|
420 | endif |
---|
421 | |
---|
422 | |
---|
423 | !----------------------------------------------------------------------- |
---|
424 | ! Water vapour (to be generalised / ignored for other planets) |
---|
425 | |
---|
426 | if(varactive) then |
---|
427 | i_var=igcm_h2o_vap |
---|
428 | |
---|
429 | do l=1,nlayer |
---|
430 | qvar(2*l) = pq(ig,nlayer+1-l,i_var) |
---|
431 | qvar(2*l+1) = (pq(ig,nlayer+1-l,i_var)+pq(ig, |
---|
432 | $ max(nlayer-l,1),i_var))/2 ! Average approximation as for temperature... |
---|
433 | end do |
---|
434 | qvar(1)=qvar(2) |
---|
435 | qvar(2*nlayermx+1)=qsurf(ig,i_var) |
---|
436 | |
---|
437 | elseif(varfixed)then |
---|
438 | |
---|
439 | do l=1,nlayermx ! here we will assign fixed water vapour profiles globally |
---|
440 | RH = satval * ((pplay(ig,l)/pplev(ig,1) - 0.02) / 0.98) |
---|
441 | if(RH.lt.0.0) RH=0.0 |
---|
442 | |
---|
443 | ptemp=pplay(ig,l) |
---|
444 | Ttemp=pt(ig,l) |
---|
445 | call watersat_2(Ttemp,ptemp,qsat) |
---|
446 | |
---|
447 | !pq_temp(l) = qsat ! fully saturated everywhere |
---|
448 | pq_temp(l) = RH * qsat ! ~realistic profile (80% saturation at ground) |
---|
449 | end do |
---|
450 | |
---|
451 | do l=1,nlayer |
---|
452 | qvar(2*l) = pq_temp(nlayer+1-l) |
---|
453 | qvar(2*l+1) = (pq_temp(nlayer+1-l)+pq_temp( |
---|
454 | $ max(nlayer-l,1)))/2 |
---|
455 | end do |
---|
456 | qvar(1)=qvar(2) |
---|
457 | |
---|
458 | ! Lowest layer of atmosphere |
---|
459 | RH = satval * (1 - 0.02) / 0.98 |
---|
460 | if(RH.lt.0.0) RH=0.0 |
---|
461 | |
---|
462 | ptemp = pplev(ig,1) |
---|
463 | Ttemp = tsurf(ig) |
---|
464 | call watersat_2(Ttemp,ptemp,qsat) |
---|
465 | |
---|
466 | !qvar(2*nlayermx+1)=qsat ! fully saturated everywhere |
---|
467 | qvar(2*nlayermx+1)= RH * qsat ! ~realistic profile (80% saturation at ground) |
---|
468 | |
---|
469 | else |
---|
470 | do k=1,L_LEVELS |
---|
471 | qvar(k) = 1.0D-7 |
---|
472 | end do |
---|
473 | end if |
---|
474 | |
---|
475 | ! Keep values inside limits for which we have radiative transfer coefficients |
---|
476 | if(L_REFVAR.gt.1)then ! there was a bug here! |
---|
477 | do k=1,L_LEVELS |
---|
478 | if(qvar(k).lt.wrefvar(1))then |
---|
479 | qvar(k)=wrefvar(1)+1.0e-8 |
---|
480 | elseif(qvar(k).gt.wrefvar(L_REFVAR))then |
---|
481 | qvar(k)=wrefvar(L_REFVAR)-1.0e-8 |
---|
482 | endif |
---|
483 | end do |
---|
484 | endif |
---|
485 | |
---|
486 | !----------------------------------------------------------------------- |
---|
487 | ! Pressure and temperature |
---|
488 | |
---|
489 | DO l=1,nlayer |
---|
490 | plevrad(2*l) = pplay(ig,nlayer+1-l)/scalep |
---|
491 | plevrad(2*l+1) = pplev(ig,nlayer+1-l)/scalep |
---|
492 | tlevrad(2*l) = pt(ig,nlayer+1-l) |
---|
493 | tlevrad(2*l+1) = (pt(ig,nlayer+1-l)+pt(ig, |
---|
494 | $ max(nlayer-l,1)))/2 |
---|
495 | END DO |
---|
496 | |
---|
497 | plevrad(1) = 0 |
---|
498 | plevrad(2) = max(pgasmin,0.0001*plevrad(3)) |
---|
499 | |
---|
500 | tlevrad(1) = tlevrad(2) |
---|
501 | tlevrad(2*nlayermx+1)=tsurf(ig) |
---|
502 | |
---|
503 | |
---|
504 | if(crisp)then |
---|
505 | open(111,file='/u/rwlmd/CrispLBL/p.dat') |
---|
506 | open(112,file='/u/rwlmd/CrispLBL/T.dat') |
---|
507 | do k=1,L_LEVELS |
---|
508 | read(111,*) plevrad(k) |
---|
509 | read(112,*) tlevrad(k) |
---|
510 | plevrad(k)=plevrad(k)*1000.0 |
---|
511 | enddo |
---|
512 | close(111) |
---|
513 | close(112) |
---|
514 | endif |
---|
515 | |
---|
516 | tmid(1) = tlevrad(2) |
---|
517 | tmid(2) = tlevrad(2) |
---|
518 | pmid(1) = plevrad(2) |
---|
519 | pmid(2) = plevrad(2) |
---|
520 | |
---|
521 | DO l=1,L_NLAYRAD-1 |
---|
522 | tmid(2*l+1) = tlevrad(2*l+1) |
---|
523 | tmid(2*l+2) = tlevrad(2*l+1) |
---|
524 | pmid(2*l+1) = plevrad(2*l+1) |
---|
525 | pmid(2*l+2) = plevrad(2*l+1) |
---|
526 | END DO |
---|
527 | pmid(L_LEVELS) = plevrad(L_LEVELS) |
---|
528 | tmid(L_LEVELS) = tlevrad(L_LEVELS) |
---|
529 | |
---|
530 | ! test for out-of-bounds pressure |
---|
531 | if(plevrad(3).lt.pgasmin)then |
---|
532 | print*,'Minimum pressure is outside the radiative' |
---|
533 | print*,'transfer kmatrix bounds, exiting.' |
---|
534 | call abort |
---|
535 | elseif(plevrad(L_LEVELS).gt.pgasmax)then |
---|
536 | print*,'Maximum pressure is outside the radiative' |
---|
537 | print*,'transfer kmatrix bounds, exiting.' |
---|
538 | call abort |
---|
539 | endif |
---|
540 | |
---|
541 | ! test for out-of-bounds temperature |
---|
542 | do k=1,L_LEVELS |
---|
543 | if(tlevrad(k).lt.tgasmin)then |
---|
544 | print*,'Minimum temperature is outside the radiative' |
---|
545 | print*,'transfer kmatrix bounds, exiting.' |
---|
546 | call abort |
---|
547 | elseif(tlevrad(k).gt.tgasmax)then |
---|
548 | print*,'Maximum temperature is outside the radiative' |
---|
549 | print*,'transfer kmatrix bounds, exiting.' |
---|
550 | call abort |
---|
551 | endif |
---|
552 | enddo |
---|
553 | |
---|
554 | !======================================================================= |
---|
555 | ! Calling the main radiative transfer subroutines |
---|
556 | |
---|
557 | |
---|
558 | |
---|
559 | !----------------------------------------------------------------------- |
---|
560 | ! Shortwave |
---|
561 | |
---|
562 | IF(fract(ig) .GE. 1.0e-4) THEN ! only during daylight! |
---|
563 | |
---|
564 | fluxtoplanet=0. |
---|
565 | DO nw=1,L_NSPECTV |
---|
566 | stel_fract(nw)= stel(nw) * fract(ig) |
---|
567 | fluxtoplanet=fluxtoplanet + stel_fract(nw) |
---|
568 | END DO |
---|
569 | |
---|
570 | call optcv(dtauv,tauv,taucumv,plevrad, |
---|
571 | $ qxvaer,qsvaer,gvaer,wbarv,cosbv,tauray,tauaero, |
---|
572 | $ tmid,pmid,taugsurf,qvar) |
---|
573 | |
---|
574 | call sfluxv(dtauv,tauv,taucumv,albv,dwnv,wbarv,cosbv, |
---|
575 | $ acosz,stel_fract,gweight,nfluxtopv,nfluxtopv_nu, |
---|
576 | $ fmnetv,fluxupv,fluxdnv,fzerov,taugsurf) |
---|
577 | |
---|
578 | ELSE ! during the night, fluxes = 0 |
---|
579 | nfluxtopv=0.0 |
---|
580 | DO l=1,L_NLAYRAD |
---|
581 | fmnetv(l)=0.0 |
---|
582 | fluxupv(l)=0.0 |
---|
583 | fluxdnv(l)=0.0 |
---|
584 | END DO |
---|
585 | END IF |
---|
586 | |
---|
587 | !----------------------------------------------------------------------- |
---|
588 | ! Longwave |
---|
589 | |
---|
590 | call optci(plevrad,tlevrad,dtaui,taucumi, |
---|
591 | $ qxiaer,qsiaer,giaer,cosbi,wbari,tauaero,tmid,pmid, |
---|
592 | $ taugsurfi,qvar) |
---|
593 | |
---|
594 | call sfluxi(plevrad,tlevrad,dtaui,taucumi,ubari,albi, |
---|
595 | $ wnoi,dwni,cosbi,wbari,gweight,nfluxtopi,nfluxtopi_nu, |
---|
596 | $ fmneti,fluxupi,fluxdni,fluxupi_nu,fzeroi,taugsurfi) |
---|
597 | |
---|
598 | !----------------------------------------------------------------------- |
---|
599 | ! Transformation of the correlated-k code outputs |
---|
600 | ! (into dtlw, dtsw, fluxsurf_lw, fluxsurf_sw, fluxtop_lw, fluxtop_sw) |
---|
601 | |
---|
602 | fluxtop_lw(ig) = fluxupi(1) |
---|
603 | fluxsurf_lw(ig) = fluxdni(L_NLAYRAD) |
---|
604 | fluxtop_sw(ig) = fluxupv(1) |
---|
605 | fluxsurf_sw(ig) = fluxdnv(L_NLAYRAD) |
---|
606 | |
---|
607 | ! Flux incident at the top of the atmosphere |
---|
608 | fluxtop_dn(ig)=fluxdnv(1) |
---|
609 | |
---|
610 | !if(lastcall)then |
---|
611 | !print*,'albv=',albv |
---|
612 | !print*,'albi=',albi |
---|
613 | !print*,'fluxupi',fluxupi |
---|
614 | !print*,'fluxirnet',fluxdni-fluxupi |
---|
615 | !print*,'fmneti',fmneti |
---|
616 | !print*,'fluxdnv',fluxdnv |
---|
617 | !print*,'fluxtop_sw=',fluxtop_sw |
---|
618 | !print*,'fluxsurf_sw=',fluxsurf_sw |
---|
619 | !print*,'fluxtop_lw=',fluxtop_lw |
---|
620 | !print*,'fluxsurf_lw=',fluxsurf_lw |
---|
621 | !endif |
---|
622 | |
---|
623 | if(crisp)then |
---|
624 | open(111,file='/u/rwlmd/CrispLBL/GCMfluxdn.dat') |
---|
625 | open(112,file='/u/rwlmd/CrispLBL/GCMfluxup.dat') |
---|
626 | do k=1,L_NLAYRAD |
---|
627 | write(111,*) fluxdni(k) |
---|
628 | write(112,*) fluxupi(k) |
---|
629 | enddo |
---|
630 | close(111) |
---|
631 | close(112) |
---|
632 | print*,'fluxdni',fluxdni |
---|
633 | print*,'fluxupi',fluxupi |
---|
634 | call abort |
---|
635 | endif |
---|
636 | |
---|
637 | ! IR spectral output, for exoplanet observational comparison |
---|
638 | if(specOLR)then |
---|
639 | do nw=1,L_NSPECTI |
---|
640 | OLR_nu(ig,nw)=nfluxtopi_nu(nw) |
---|
641 | end do |
---|
642 | endif |
---|
643 | |
---|
644 | ! Finally, the heating rates |
---|
645 | |
---|
646 | if(nonideal)then |
---|
647 | |
---|
648 | DO l=2,L_NLAYRAD |
---|
649 | dtsw(ig,L_NLAYRAD+1-l)=(fmnetv(l)-fmnetv(l-1)) |
---|
650 | $ *g/(cpp3D(ig,L_NLAYRAD+1-l) |
---|
651 | $ *scalep*(plevrad(2*l+1)-plevrad(2*l-1))) |
---|
652 | dtlw(ig,L_NLAYRAD+1-l)=(fmneti(l)-fmneti(l-1)) |
---|
653 | $ *g/(cpp3D(ig,L_NLAYRAD+1-l) |
---|
654 | $ *scalep*(plevrad(2*l+1)-plevrad(2*l-1))) |
---|
655 | END DO |
---|
656 | |
---|
657 | ! These are values at top of atmosphere |
---|
658 | dtsw(ig,L_NLAYRAD)=(fmnetv(1)-nfluxtopv) |
---|
659 | $ *g/(cpp3D(ig,L_NLAYRAD)*scalep*(plevrad(3)-plevrad(1))) |
---|
660 | dtlw(ig,L_NLAYRAD)=(fmneti(1)-nfluxtopi) |
---|
661 | $ *g/(cpp3D(ig,L_NLAYRAD)*scalep*(plevrad(3)-plevrad(1))) |
---|
662 | |
---|
663 | else |
---|
664 | |
---|
665 | DO l=2,L_NLAYRAD |
---|
666 | dtsw(ig,L_NLAYRAD+1-l)=(fmnetv(l)-fmnetv(l-1)) |
---|
667 | $ *g/(cpp*scalep*(plevrad(2*l+1)-plevrad(2*l-1))) |
---|
668 | dtlw(ig,L_NLAYRAD+1-l)=(fmneti(l)-fmneti(l-1)) |
---|
669 | $ *g/(cpp*scalep*(plevrad(2*l+1)-plevrad(2*l-1))) |
---|
670 | END DO |
---|
671 | |
---|
672 | ! These are values at top of atmosphere |
---|
673 | dtsw(ig,L_NLAYRAD)=(fmnetv(1)-nfluxtopv) |
---|
674 | $ *g/(cpp*scalep*(plevrad(3)-plevrad(1))) |
---|
675 | dtlw(ig,L_NLAYRAD)=(fmneti(1)-nfluxtopi) |
---|
676 | $ *g/(cpp*scalep*(plevrad(3)-plevrad(1))) |
---|
677 | |
---|
678 | endif |
---|
679 | |
---|
680 | !if(lastcall)then |
---|
681 | ! ! print*,'tlevrad=',tlevrad |
---|
682 | ! print*,'tsurf=',tsurf |
---|
683 | ! print*,'pt=',pt |
---|
684 | ! print*,'fmnetv',shape(fmnetv) |
---|
685 | ! print*,'fmneti',shape(fmneti) |
---|
686 | !!endif |
---|
687 | |
---|
688 | ! --------------------------------------------------------------- |
---|
689 | end do ! end of big loop over every GCM column (ig = 1:ngrid) |
---|
690 | |
---|
691 | |
---|
692 | !----------------------------------------------------------------------- |
---|
693 | ! Additional diagnostics |
---|
694 | |
---|
695 | ! IR spectral output, for exoplanet observational comparison |
---|
696 | if(specOLR)then |
---|
697 | if(ngrid.ne.1)then |
---|
698 | call writediagspec(ngrid,"OLR3D", |
---|
699 | & "OLR(lon,lat,band)","W m^-2",3,OLR_nu) |
---|
700 | endif |
---|
701 | endif |
---|
702 | |
---|
703 | if(lastcall)then ! for kasthop |
---|
704 | |
---|
705 | if(ngrid.eq.1)then |
---|
706 | |
---|
707 | print*,'Saving tsurf,psurf in surf_vals.out...' |
---|
708 | open(116,file='surf_vals.out') |
---|
709 | write(116,*) tsurf(1),pplev(1,1), |
---|
710 | & fluxtop_dn(1) - fluxtop_sw(1),fluxtop_lw(1) |
---|
711 | close(116) |
---|
712 | |
---|
713 | if(specOLR)then |
---|
714 | open(117,file='OLRnu.out') |
---|
715 | do nw=1,L_NSPECTI |
---|
716 | write(117,*) OLR_nu(1,nw) |
---|
717 | enddo |
---|
718 | close(117) |
---|
719 | endif |
---|
720 | ! OLR vs altitude: do it as a .txt file |
---|
721 | OLRz=.false. |
---|
722 | if(OLRz)then |
---|
723 | print*,'saving IR vertical flux for OLRz...' |
---|
724 | open(118,file='OLRz_plevs.out') |
---|
725 | open(119,file='OLRz.out') |
---|
726 | do l=1,L_NLAYRAD |
---|
727 | write(118,*) plevrad(2*l) |
---|
728 | do nw=1,L_NSPECTI |
---|
729 | write(119,*) fluxupi_nu(l,nw) |
---|
730 | enddo |
---|
731 | enddo |
---|
732 | close(118) |
---|
733 | close(119) |
---|
734 | endif |
---|
735 | endif |
---|
736 | |
---|
737 | endif |
---|
738 | |
---|
739 | end |
---|