1 | MODULE callcorrk_mod |
---|
2 | |
---|
3 | IMPLICIT NONE |
---|
4 | |
---|
5 | CONTAINS |
---|
6 | |
---|
7 | subroutine callcorrk(ngrid,nlayer,pq,nq,qsurf, & |
---|
8 | albedo,albedo_equivalent,emis,mu0,pplev,pplay,pt, & |
---|
9 | tsurf,fract,dist_star,aerosol,muvar, & |
---|
10 | dtlw,dtsw,fluxsurf_lw, & |
---|
11 | fluxsurf_sw,fluxsurfabs_sw,fluxtop_lw, & |
---|
12 | fluxabs_sw,fluxtop_dn, & |
---|
13 | OLR_nu,OSR_nu, & |
---|
14 | int_dtaui,int_dtauv, & |
---|
15 | tau_col,cloudfrac,totcloudfrac, & |
---|
16 | clearsky,firstcall,lastcall) |
---|
17 | |
---|
18 | use mod_phys_lmdz_para, only : is_master |
---|
19 | use radinc_h, only: L_NSPECTV, L_NSPECTI, naerkind, banddir, corrkdir,& |
---|
20 | L_LEVELS, L_NGAUSS, L_NLEVRAD, L_NLAYRAD, L_REFVAR |
---|
21 | use radcommon_h, only: wrefvar, Cmk, fzeroi, fzerov, gasi, gasv, & |
---|
22 | glat_ig, gweight, pfgasref, pgasmax, pgasmin, & |
---|
23 | pgasref, tgasmax, tgasmin, tgasref, scalep, & |
---|
24 | ubari, wnoi, stellarf, glat, dwnv, dwni, tauray |
---|
25 | use watercommon_h, only: psat_water, epsi |
---|
26 | use datafile_mod, only: datadir |
---|
27 | use ioipsl_getin_p_mod, only: getin_p |
---|
28 | use gases_h, only: ngasmx |
---|
29 | use radii_mod, only : su_aer_radii,co2_reffrad,h2o_reffrad,dust_reffrad,h2so4_reffrad,back2lay_reffrad |
---|
30 | use aerosol_mod, only : iaero_co2,iaero_h2o,iaero_dust,iaero_h2so4, iaero_back2lay, iaero_nh3, iaero_aurora |
---|
31 | use tracer_h, only: igcm_h2o_vap, igcm_h2o_ice, igcm_co2_ice |
---|
32 | use comcstfi_mod, only: pi, mugaz, cpp |
---|
33 | use callkeys_mod, only: varactive,diurnal,tracer,water,varfixed,satval,diagdtau, & |
---|
34 | kastprof,strictboundcorrk,specOLR,CLFvarying |
---|
35 | use optcv_mod, only: optcv |
---|
36 | use optci_mod, only: optci |
---|
37 | implicit none |
---|
38 | |
---|
39 | !================================================================== |
---|
40 | ! |
---|
41 | ! Purpose |
---|
42 | ! ------- |
---|
43 | ! Solve the radiative transfer using the correlated-k method for |
---|
44 | ! the gaseous absorption and the Toon et al. (1989) method for |
---|
45 | ! scatttering due to aerosols. |
---|
46 | ! |
---|
47 | ! Authors |
---|
48 | ! ------- |
---|
49 | ! Emmanuel 01/2001, Forget 09/2001 |
---|
50 | ! Robin Wordsworth (2009) |
---|
51 | ! |
---|
52 | !================================================================== |
---|
53 | |
---|
54 | !----------------------------------------------------------------------- |
---|
55 | ! Declaration of the arguments (INPUT - OUTPUT) on the LMD GCM grid |
---|
56 | ! Layer #1 is the layer near the ground. |
---|
57 | ! Layer #nlayer is the layer at the top. |
---|
58 | !----------------------------------------------------------------------- |
---|
59 | |
---|
60 | |
---|
61 | ! INPUT |
---|
62 | INTEGER,INTENT(IN) :: ngrid ! Number of atmospheric columns. |
---|
63 | INTEGER,INTENT(IN) :: nlayer ! Number of atmospheric layers. |
---|
64 | REAL,INTENT(IN) :: pq(ngrid,nlayer,nq) ! Tracers (kg/kg_of_air). |
---|
65 | INTEGER,INTENT(IN) :: nq ! Number of tracers. |
---|
66 | REAL,INTENT(IN) :: qsurf(ngrid,nq) ! Tracers on surface (kg.m-2). |
---|
67 | REAL,INTENT(IN) :: albedo(ngrid,L_NSPECTV) ! Spectral Short Wavelengths Albedo. By MT2015 |
---|
68 | REAL,INTENT(IN) :: emis(ngrid) ! Long Wave emissivity. |
---|
69 | REAL,INTENT(IN) :: mu0(ngrid) ! Cosine of sun incident angle. |
---|
70 | REAL,INTENT(IN) :: pplev(ngrid,nlayer+1) ! Inter-layer pressure (Pa). |
---|
71 | REAL,INTENT(IN) :: pplay(ngrid,nlayer) ! Mid-layer pressure (Pa). |
---|
72 | REAL,INTENT(IN) :: pt(ngrid,nlayer) ! Air temperature (K). |
---|
73 | REAL,INTENT(IN) :: tsurf(ngrid) ! Surface temperature (K). |
---|
74 | REAL,INTENT(IN) :: fract(ngrid) ! Fraction of day. |
---|
75 | REAL,INTENT(IN) :: dist_star ! Distance star-planet (AU). |
---|
76 | REAL,INTENT(IN) :: muvar(ngrid,nlayer+1) |
---|
77 | REAL,INTENT(IN) :: cloudfrac(ngrid,nlayer) ! Fraction of clouds (%). |
---|
78 | logical,intent(in) :: clearsky |
---|
79 | logical,intent(in) :: firstcall ! Signals first call to physics. |
---|
80 | logical,intent(in) :: lastcall ! Signals last call to physics. |
---|
81 | |
---|
82 | ! OUTPUT |
---|
83 | REAL,INTENT(OUT) :: aerosol(ngrid,nlayer,naerkind) ! Aerosol tau (kg/kg). |
---|
84 | REAL,INTENT(OUT) :: dtlw(ngrid,nlayer) ! Heating rate (K/s) due to LW radiation. |
---|
85 | REAL,INTENT(OUT) :: dtsw(ngrid,nlayer) ! Heating rate (K/s) due to SW radiation. |
---|
86 | REAL,INTENT(OUT) :: fluxsurf_lw(ngrid) ! Incident LW flux to surf (W/m2). |
---|
87 | REAL,INTENT(OUT) :: fluxsurf_sw(ngrid) ! Incident SW flux to surf (W/m2) |
---|
88 | REAL,INTENT(OUT) :: fluxsurfabs_sw(ngrid) ! Absorbed SW flux by the surface (W/m2). By MT2015. |
---|
89 | REAL,INTENT(OUT) :: fluxtop_lw(ngrid) ! Outgoing LW flux to space (W/m2). |
---|
90 | REAL,INTENT(OUT) :: fluxabs_sw(ngrid) ! SW flux absorbed by the planet (W/m2). |
---|
91 | REAL,INTENT(OUT) :: fluxtop_dn(ngrid) ! Incident top of atmosphere SW flux (W/m2). |
---|
92 | REAL,INTENT(OUT) :: OLR_nu(ngrid,L_NSPECTI) ! Outgoing LW radition in each band (Normalized to the band width (W/m2/cm-1). |
---|
93 | REAL,INTENT(OUT) :: OSR_nu(ngrid,L_NSPECTV) ! Outgoing SW radition in each band (Normalized to the band width (W/m2/cm-1). |
---|
94 | REAL,INTENT(OUT) :: tau_col(ngrid) ! Diagnostic from aeropacity. |
---|
95 | REAL,INTENT(OUT) :: albedo_equivalent(ngrid) ! Spectrally Integrated Albedo. For Diagnostic. By MT2015 |
---|
96 | REAL,INTENT(OUT) :: totcloudfrac(ngrid) ! Column Fraction of clouds (%). |
---|
97 | REAL,INTENT(OUT) :: int_dtaui(ngrid,nlayer,L_NSPECTI) ! VI optical thickness of layers within narrowbands for diags (). |
---|
98 | REAL,INTENT(OUT) :: int_dtauv(ngrid,nlayer,L_NSPECTV) ! IR optical thickness of layers within narrowbands for diags (). |
---|
99 | |
---|
100 | |
---|
101 | |
---|
102 | |
---|
103 | |
---|
104 | ! Globally varying aerosol optical properties on GCM grid ; not needed everywhere so not in radcommon_h. |
---|
105 | REAL :: QVISsQREF3d(ngrid,nlayer,L_NSPECTV,naerkind) |
---|
106 | REAL :: omegaVIS3d(ngrid,nlayer,L_NSPECTV,naerkind) |
---|
107 | REAL :: gVIS3d(ngrid,nlayer,L_NSPECTV,naerkind) |
---|
108 | REAL :: QIRsQREF3d(ngrid,nlayer,L_NSPECTI,naerkind) |
---|
109 | REAL :: omegaIR3d(ngrid,nlayer,L_NSPECTI,naerkind) |
---|
110 | REAL :: gIR3d(ngrid,nlayer,L_NSPECTI,naerkind) |
---|
111 | |
---|
112 | ! REAL :: omegaREFvis3d(ngrid,nlayer,naerkind) |
---|
113 | ! REAL :: omegaREFir3d(ngrid,nlayer,naerkind) ! not sure of the point of these... |
---|
114 | |
---|
115 | REAL,ALLOCATABLE,SAVE :: reffrad(:,:,:) ! aerosol effective radius (m) |
---|
116 | REAL,ALLOCATABLE,SAVE :: nueffrad(:,:,:) ! aerosol effective variance |
---|
117 | !$OMP THREADPRIVATE(reffrad,nueffrad) |
---|
118 | |
---|
119 | !----------------------------------------------------------------------- |
---|
120 | ! Declaration of the variables required by correlated-k subroutines |
---|
121 | ! Numbered from top to bottom (unlike in the GCM) |
---|
122 | !----------------------------------------------------------------------- |
---|
123 | |
---|
124 | REAL*8 tmid(L_LEVELS),pmid(L_LEVELS) |
---|
125 | REAL*8 tlevrad(L_LEVELS),plevrad(L_LEVELS) |
---|
126 | |
---|
127 | ! Optical values for the optci/cv subroutines |
---|
128 | REAL*8 stel(L_NSPECTV),stel_fract(L_NSPECTV) |
---|
129 | ! NB: Arrays below are "save" to avoid reallocating them at every call |
---|
130 | ! not because their content needs be reused from call to the next |
---|
131 | REAL*8,allocatable,save :: dtaui(:,:,:) |
---|
132 | REAL*8,allocatable,save :: dtauv(:,:,:) |
---|
133 | REAL*8,allocatable,save :: cosbv(:,:,:) |
---|
134 | REAL*8,allocatable,save :: cosbi(:,:,:) |
---|
135 | REAL*8,allocatable,save :: wbari(:,:,:) |
---|
136 | REAL*8,allocatable,save :: wbarv(:,:,:) |
---|
137 | REAL*8,allocatable,save :: tauv(:,:,:) |
---|
138 | REAL*8,allocatable,save :: taucumv(:,:,:) |
---|
139 | REAL*8,allocatable,save :: taucumi(:,:,:) |
---|
140 | |
---|
141 | REAL*8 tauaero(L_LEVELS,naerkind) |
---|
142 | REAL*8 nfluxtopv,nfluxtopi,nfluxtop,fluxtopvdn |
---|
143 | REAL*8 nfluxoutv_nu(L_NSPECTV) ! Outgoing band-resolved VI flux at TOA (W/m2). |
---|
144 | REAL*8 nfluxtopi_nu(L_NSPECTI) ! Net band-resolved IR flux at TOA (W/m2). |
---|
145 | REAL*8 fluxupi_nu(L_NLAYRAD,L_NSPECTI) ! For 1D diagnostic. |
---|
146 | REAL*8 fmneti(L_NLAYRAD),fmnetv(L_NLAYRAD) |
---|
147 | REAL*8 fluxupv(L_NLAYRAD),fluxupi(L_NLAYRAD) |
---|
148 | REAL*8 fluxdnv(L_NLAYRAD),fluxdni(L_NLAYRAD) |
---|
149 | REAL*8 albi,acosz |
---|
150 | REAL*8 albv(L_NSPECTV) ! Spectral Visible Albedo. |
---|
151 | |
---|
152 | INTEGER ig,l,k,nw,iaer |
---|
153 | |
---|
154 | real,save :: szangle |
---|
155 | logical,save :: global1d |
---|
156 | !$OMP THREADPRIVATE(szangle,global1d) |
---|
157 | |
---|
158 | real*8,allocatable,save :: taugsurf(:,:) |
---|
159 | real*8,allocatable,save :: taugsurfi(:,:) |
---|
160 | real*8 qvar(L_LEVELS) ! Mixing ratio of variable component (mol/mol). |
---|
161 | |
---|
162 | ! Local aerosol optical properties for each column on RADIATIVE grid. |
---|
163 | real*8,save,allocatable :: QXVAER(:,:,:) |
---|
164 | real*8,save,allocatable :: QSVAER(:,:,:) |
---|
165 | real*8,save,allocatable :: GVAER(:,:,:) |
---|
166 | real*8,save,allocatable :: QXIAER(:,:,:) |
---|
167 | real*8,save,allocatable :: QSIAER(:,:,:) |
---|
168 | real*8,save,allocatable :: GIAER(:,:,:) |
---|
169 | |
---|
170 | real, dimension(:,:,:), save, allocatable :: QREFvis3d |
---|
171 | real, dimension(:,:,:), save, allocatable :: QREFir3d |
---|
172 | !$OMP THREADPRIVATE(QXVAER,QSVAER,GVAER,QXIAER,QSIAER,GIAER,QREFvis3d,QREFir3d) |
---|
173 | |
---|
174 | |
---|
175 | ! Miscellaneous : |
---|
176 | real*8 temp,temp1,temp2,pweight |
---|
177 | character(len=10) :: tmp1 |
---|
178 | character(len=10) :: tmp2 |
---|
179 | |
---|
180 | ! For fixed water vapour profiles. |
---|
181 | integer i_var |
---|
182 | real RH |
---|
183 | real*8 pq_temp(nlayer) |
---|
184 | ! real(KIND=r8) :: pq_temp(nlayer) ! better F90 way.. DOESNT PORT TO F77!!! |
---|
185 | real psat,qsat |
---|
186 | |
---|
187 | logical OLRz |
---|
188 | real*8 NFLUXGNDV_nu(L_NSPECTV) |
---|
189 | |
---|
190 | ! Included by RW for runaway greenhouse 1D study. |
---|
191 | real vtmp(nlayer) |
---|
192 | REAL*8 muvarrad(L_LEVELS) |
---|
193 | |
---|
194 | ! Included by MT for albedo calculations. |
---|
195 | REAL*8 albedo_temp(L_NSPECTV) ! For equivalent albedo calculation. |
---|
196 | REAL*8 surface_stellar_flux ! Stellar flux reaching the surface. Useful for equivalent albedo calculation. |
---|
197 | |
---|
198 | |
---|
199 | !=============================================================== |
---|
200 | ! I.a Initialization on first call |
---|
201 | !=============================================================== |
---|
202 | |
---|
203 | |
---|
204 | if(firstcall) then |
---|
205 | |
---|
206 | ! test on allocated necessary because of CLFvarying (two calls to callcorrk in physiq) |
---|
207 | if(.not.allocated(QXVAER)) allocate(QXVAER(L_LEVELS,L_NSPECTV,naerkind)) |
---|
208 | if(.not.allocated(QSVAER)) allocate(QSVAER(L_LEVELS,L_NSPECTV,naerkind)) |
---|
209 | if(.not.allocated(GVAER)) allocate(GVAER(L_LEVELS,L_NSPECTV,naerkind)) |
---|
210 | if(.not.allocated(QXIAER)) allocate(QXIAER(L_LEVELS,L_NSPECTI,naerkind)) |
---|
211 | if(.not.allocated(QSIAER)) allocate(QSIAER(L_LEVELS,L_NSPECTI,naerkind)) |
---|
212 | if(.not.allocated(GIAER)) allocate(GIAER(L_LEVELS,L_NSPECTI,naerkind)) |
---|
213 | |
---|
214 | !!! ALLOCATED instances are necessary because of CLFvarying (strategy to call callcorrk twice in physiq...) |
---|
215 | IF(.not.ALLOCATED(QREFvis3d)) ALLOCATE(QREFvis3d(ngrid,nlayer,naerkind)) |
---|
216 | IF(.not.ALLOCATED(QREFir3d)) ALLOCATE(QREFir3d(ngrid,nlayer,naerkind)) |
---|
217 | ! Effective radius and variance of the aerosols |
---|
218 | IF(.not.ALLOCATED(reffrad)) allocate(reffrad(ngrid,nlayer,naerkind)) |
---|
219 | IF(.not.ALLOCATED(nueffrad)) allocate(nueffrad(ngrid,nlayer,naerkind)) |
---|
220 | |
---|
221 | #ifndef MESOSCALE |
---|
222 | call system('rm -f surf_vals_long.out') |
---|
223 | #endif |
---|
224 | |
---|
225 | if(naerkind.gt.4)then |
---|
226 | print*,'Code not general enough to deal with naerkind > 4 yet.' |
---|
227 | call abort |
---|
228 | endif |
---|
229 | call su_aer_radii(ngrid,nlayer,reffrad,nueffrad) |
---|
230 | |
---|
231 | |
---|
232 | !-------------------------------------------------- |
---|
233 | ! Set up correlated k |
---|
234 | !-------------------------------------------------- |
---|
235 | |
---|
236 | |
---|
237 | print*, "callcorrk: Correlated-k data base folder:",trim(datadir) |
---|
238 | call getin_p("corrkdir",corrkdir) |
---|
239 | print*, "corrkdir = ",corrkdir |
---|
240 | write( tmp1, '(i3)' ) L_NSPECTI |
---|
241 | write( tmp2, '(i3)' ) L_NSPECTV |
---|
242 | banddir=trim(adjustl(tmp1))//'x'//trim(adjustl(tmp2)) |
---|
243 | banddir=trim(adjustl(corrkdir))//'/'//trim(adjustl(banddir)) |
---|
244 | |
---|
245 | call setspi ! Basic infrared properties. |
---|
246 | call setspv ! Basic visible properties. |
---|
247 | call sugas_corrk ! Set up gaseous absorption properties. |
---|
248 | call suaer_corrk ! Set up aerosol optical properties. |
---|
249 | |
---|
250 | |
---|
251 | ! now that L_NGAUSS has been initialized (by sugas_corrk) |
---|
252 | ! allocate related arrays |
---|
253 | if(.not.allocated(dtaui)) ALLOCATE(dtaui(L_NLAYRAD,L_NSPECTI,L_NGAUSS)) |
---|
254 | if(.not.allocated(dtauv)) ALLOCATE(dtauv(L_NLAYRAD,L_NSPECTV,L_NGAUSS)) |
---|
255 | if(.not.allocated(cosbv)) ALLOCATE(cosbv(L_NLAYRAD,L_NSPECTV,L_NGAUSS)) |
---|
256 | if(.not.allocated(cosbi)) ALLOCATE(cosbi(L_NLAYRAD,L_NSPECTI,L_NGAUSS)) |
---|
257 | if(.not.allocated(wbari)) ALLOCATE(wbari(L_NLAYRAD,L_NSPECTI,L_NGAUSS)) |
---|
258 | if(.not.allocated(wbarv)) ALLOCATE(wbarv(L_NLAYRAD,L_NSPECTV,L_NGAUSS)) |
---|
259 | if(.not.allocated(tauv)) ALLOCATE(tauv(L_NLEVRAD,L_NSPECTV,L_NGAUSS)) |
---|
260 | if(.not.allocated(taucumv)) ALLOCATE(taucumv(L_LEVELS,L_NSPECTV,L_NGAUSS)) |
---|
261 | if(.not.allocated(taucumi)) ALLOCATE(taucumi(L_LEVELS,L_NSPECTI,L_NGAUSS)) |
---|
262 | if(.not.allocated(taugsurf)) ALLOCATE(taugsurf(L_NSPECTV,L_NGAUSS-1)) |
---|
263 | if(.not.allocated(taugsurfi)) ALLOCATE(taugsurfi(L_NSPECTI,L_NGAUSS-1)) |
---|
264 | |
---|
265 | if((igcm_h2o_vap.eq.0) .and. varactive)then |
---|
266 | print*,'varactive in callcorrk but no h2o_vap tracer.' |
---|
267 | stop |
---|
268 | endif |
---|
269 | |
---|
270 | OLR_nu(:,:) = 0. |
---|
271 | OSR_nu(:,:) = 0. |
---|
272 | |
---|
273 | if (ngrid.eq.1) then |
---|
274 | PRINT*, 'Simulate global averaged conditions ?' |
---|
275 | global1d = .false. ! default value |
---|
276 | call getin_p("global1d",global1d) |
---|
277 | write(*,*) "global1d = ",global1d |
---|
278 | |
---|
279 | ! Test of incompatibility : if global1d is true, there should not be any diurnal cycle. |
---|
280 | if (global1d.and.diurnal) then |
---|
281 | print*,'if global1d is true, diurnal must be set to false' |
---|
282 | stop |
---|
283 | endif |
---|
284 | |
---|
285 | if (global1d) then |
---|
286 | PRINT *,'Solar Zenith angle (deg.) ?' |
---|
287 | PRINT *,'(assumed for averaged solar flux S/4)' |
---|
288 | szangle=60.0 ! default value |
---|
289 | call getin_p("szangle",szangle) |
---|
290 | write(*,*) "szangle = ",szangle |
---|
291 | endif |
---|
292 | endif |
---|
293 | |
---|
294 | end if ! of if (firstcall) |
---|
295 | |
---|
296 | !======================================================================= |
---|
297 | ! I.b Initialization on every call |
---|
298 | !======================================================================= |
---|
299 | |
---|
300 | qxvaer(:,:,:)=0.0 |
---|
301 | qsvaer(:,:,:)=0.0 |
---|
302 | gvaer(:,:,:) =0.0 |
---|
303 | |
---|
304 | qxiaer(:,:,:)=0.0 |
---|
305 | qsiaer(:,:,:)=0.0 |
---|
306 | giaer(:,:,:) =0.0 |
---|
307 | |
---|
308 | !-------------------------------------------------- |
---|
309 | ! Effective radius and variance of the aerosols |
---|
310 | !-------------------------------------------------- |
---|
311 | |
---|
312 | do iaer=1,naerkind |
---|
313 | |
---|
314 | if ((iaer.eq.iaero_co2).and.tracer.and.(igcm_co2_ice.gt.0)) then ! Treat condensed co2 particles. |
---|
315 | call co2_reffrad(ngrid,nlayer,nq,pq,reffrad(1,1,iaero_co2)) |
---|
316 | if (is_master) then |
---|
317 | print*,'Max. CO2 ice particle size = ',maxval(reffrad(1:ngrid,1:nlayer,iaer))/1.e-6,' um' |
---|
318 | print*,'Min. CO2 ice particle size = ',minval(reffrad(1:ngrid,1:nlayer,iaer))/1.e-6,' um' |
---|
319 | end if |
---|
320 | end if |
---|
321 | |
---|
322 | if ((iaer.eq.iaero_h2o).and.water) then ! Treat condensed water particles. To be generalized for other aerosols ... |
---|
323 | call h2o_reffrad(ngrid,nlayer,pq(1,1,igcm_h2o_ice),pt, & |
---|
324 | reffrad(1,1,iaero_h2o),nueffrad(1,1,iaero_h2o)) |
---|
325 | if (is_master) then |
---|
326 | print*,'Max. H2O cloud particle size = ',maxval(reffrad(1:ngrid,1:nlayer,iaer))/1.e-6,' um' |
---|
327 | print*,'Min. H2O cloud particle size = ',minval(reffrad(1:ngrid,1:nlayer,iaer))/1.e-6,' um' |
---|
328 | end if |
---|
329 | endif |
---|
330 | |
---|
331 | if(iaer.eq.iaero_dust)then |
---|
332 | call dust_reffrad(ngrid,nlayer,reffrad(1,1,iaero_dust)) |
---|
333 | if (is_master) then |
---|
334 | print*,'Dust particle size = ',reffrad(1,1,iaer)/1.e-6,' um' |
---|
335 | end if |
---|
336 | endif |
---|
337 | |
---|
338 | if(iaer.eq.iaero_h2so4)then |
---|
339 | call h2so4_reffrad(ngrid,nlayer,reffrad(1,1,iaero_h2so4)) |
---|
340 | if (is_master) then |
---|
341 | print*,'H2SO4 particle size =',reffrad(1,1,iaer)/1.e-6,' um' |
---|
342 | end if |
---|
343 | endif |
---|
344 | |
---|
345 | if(iaer.eq.iaero_back2lay)then |
---|
346 | call back2lay_reffrad(ngrid,reffrad(1,1,iaero_back2lay),nlayer,pplev) |
---|
347 | endif |
---|
348 | ! if(iaer.eq.iaero_nh3)then |
---|
349 | ! call nh3_reffrad(ngrid,nlayer,reffrad(1,1,iaero_nh3)) |
---|
350 | ! endif |
---|
351 | ! if(iaer.eq.iaero_aurora)then |
---|
352 | ! call aurora_reffrad(ngrid,nlayer,reffrad(1,1,iaero_aurora)) |
---|
353 | ! endif |
---|
354 | |
---|
355 | end do !iaer=1,naerkind. |
---|
356 | |
---|
357 | |
---|
358 | ! How much light do we get ? |
---|
359 | do nw=1,L_NSPECTV |
---|
360 | stel(nw)=stellarf(nw)/(dist_star**2) |
---|
361 | end do |
---|
362 | |
---|
363 | ! Get 3D aerosol optical properties. |
---|
364 | call aeroptproperties(ngrid,nlayer,reffrad,nueffrad, & |
---|
365 | QVISsQREF3d,omegaVIS3d,gVIS3d, & |
---|
366 | QIRsQREF3d,omegaIR3d,gIR3d, & |
---|
367 | QREFvis3d,QREFir3d) |
---|
368 | |
---|
369 | ! Get aerosol optical depths. |
---|
370 | call aeropacity(ngrid,nlayer,nq,pplay,pplev,pq,aerosol, & |
---|
371 | reffrad,QREFvis3d,QREFir3d, & |
---|
372 | tau_col,cloudfrac,totcloudfrac,clearsky) |
---|
373 | |
---|
374 | |
---|
375 | |
---|
376 | !----------------------------------------------------------------------- |
---|
377 | do ig=1,ngrid ! Starting Big Loop over every GCM column |
---|
378 | !----------------------------------------------------------------------- |
---|
379 | |
---|
380 | |
---|
381 | !======================================================================= |
---|
382 | ! II. Transformation of the GCM variables |
---|
383 | !======================================================================= |
---|
384 | |
---|
385 | |
---|
386 | !----------------------------------------------------------------------- |
---|
387 | ! Aerosol optical properties Qext, Qscat and g. |
---|
388 | ! The transformation in the vertical is the same as for temperature. |
---|
389 | !----------------------------------------------------------------------- |
---|
390 | |
---|
391 | |
---|
392 | do iaer=1,naerkind |
---|
393 | ! Shortwave. |
---|
394 | do nw=1,L_NSPECTV |
---|
395 | |
---|
396 | do l=1,nlayer |
---|
397 | |
---|
398 | temp1=QVISsQREF3d(ig,nlayer+1-l,nw,iaer) & |
---|
399 | *QREFvis3d(ig,nlayer+1-l,iaer) |
---|
400 | |
---|
401 | temp2=QVISsQREF3d(ig,max(nlayer-l,1),nw,iaer) & |
---|
402 | *QREFvis3d(ig,max(nlayer-l,1),iaer) |
---|
403 | |
---|
404 | qxvaer(2*l,nw,iaer) = temp1 |
---|
405 | qxvaer(2*l+1,nw,iaer)=(temp1+temp2)/2 |
---|
406 | |
---|
407 | temp1=temp1*omegavis3d(ig,nlayer+1-l,nw,iaer) |
---|
408 | temp2=temp2*omegavis3d(ig,max(nlayer-l,1),nw,iaer) |
---|
409 | |
---|
410 | qsvaer(2*l,nw,iaer) = temp1 |
---|
411 | qsvaer(2*l+1,nw,iaer)=(temp1+temp2)/2 |
---|
412 | |
---|
413 | temp1=gvis3d(ig,nlayer+1-l,nw,iaer) |
---|
414 | temp2=gvis3d(ig,max(nlayer-l,1),nw,iaer) |
---|
415 | |
---|
416 | gvaer(2*l,nw,iaer) = temp1 |
---|
417 | gvaer(2*l+1,nw,iaer)=(temp1+temp2)/2 |
---|
418 | |
---|
419 | end do ! nlayer |
---|
420 | |
---|
421 | qxvaer(1,nw,iaer)=qxvaer(2,nw,iaer) |
---|
422 | qxvaer(2*nlayer+1,nw,iaer)=0. |
---|
423 | |
---|
424 | qsvaer(1,nw,iaer)=qsvaer(2,nw,iaer) |
---|
425 | qsvaer(2*nlayer+1,nw,iaer)=0. |
---|
426 | |
---|
427 | gvaer(1,nw,iaer)=gvaer(2,nw,iaer) |
---|
428 | gvaer(2*nlayer+1,nw,iaer)=0. |
---|
429 | |
---|
430 | end do ! L_NSPECTV |
---|
431 | |
---|
432 | do nw=1,L_NSPECTI |
---|
433 | ! Longwave |
---|
434 | do l=1,nlayer |
---|
435 | |
---|
436 | temp1=QIRsQREF3d(ig,nlayer+1-l,nw,iaer) & |
---|
437 | *QREFir3d(ig,nlayer+1-l,iaer) |
---|
438 | |
---|
439 | temp2=QIRsQREF3d(ig,max(nlayer-l,1),nw,iaer) & |
---|
440 | *QREFir3d(ig,max(nlayer-l,1),iaer) |
---|
441 | |
---|
442 | qxiaer(2*l,nw,iaer) = temp1 |
---|
443 | qxiaer(2*l+1,nw,iaer)=(temp1+temp2)/2 |
---|
444 | |
---|
445 | temp1=temp1*omegair3d(ig,nlayer+1-l,nw,iaer) |
---|
446 | temp2=temp2*omegair3d(ig,max(nlayer-l,1),nw,iaer) |
---|
447 | |
---|
448 | qsiaer(2*l,nw,iaer) = temp1 |
---|
449 | qsiaer(2*l+1,nw,iaer)=(temp1+temp2)/2 |
---|
450 | |
---|
451 | temp1=gir3d(ig,nlayer+1-l,nw,iaer) |
---|
452 | temp2=gir3d(ig,max(nlayer-l,1),nw,iaer) |
---|
453 | |
---|
454 | giaer(2*l,nw,iaer) = temp1 |
---|
455 | giaer(2*l+1,nw,iaer)=(temp1+temp2)/2 |
---|
456 | |
---|
457 | end do ! nlayer |
---|
458 | |
---|
459 | qxiaer(1,nw,iaer)=qxiaer(2,nw,iaer) |
---|
460 | qxiaer(2*nlayer+1,nw,iaer)=0. |
---|
461 | |
---|
462 | qsiaer(1,nw,iaer)=qsiaer(2,nw,iaer) |
---|
463 | qsiaer(2*nlayer+1,nw,iaer)=0. |
---|
464 | |
---|
465 | giaer(1,nw,iaer)=giaer(2,nw,iaer) |
---|
466 | giaer(2*nlayer+1,nw,iaer)=0. |
---|
467 | |
---|
468 | end do ! L_NSPECTI |
---|
469 | |
---|
470 | end do ! naerkind |
---|
471 | |
---|
472 | ! Test / Correct for freaky s. s. albedo values. |
---|
473 | do iaer=1,naerkind |
---|
474 | do k=1,L_LEVELS |
---|
475 | |
---|
476 | do nw=1,L_NSPECTV |
---|
477 | if(qsvaer(k,nw,iaer).gt.1.05*qxvaer(k,nw,iaer))then |
---|
478 | print*,'Serious problems with qsvaer values' |
---|
479 | print*,'in callcorrk' |
---|
480 | call abort |
---|
481 | endif |
---|
482 | if(qsvaer(k,nw,iaer).gt.qxvaer(k,nw,iaer))then |
---|
483 | qsvaer(k,nw,iaer)=qxvaer(k,nw,iaer) |
---|
484 | endif |
---|
485 | end do |
---|
486 | |
---|
487 | do nw=1,L_NSPECTI |
---|
488 | if(qsiaer(k,nw,iaer).gt.1.05*qxiaer(k,nw,iaer))then |
---|
489 | print*,'Serious problems with qsiaer values' |
---|
490 | print*,'in callcorrk' |
---|
491 | call abort |
---|
492 | endif |
---|
493 | if(qsiaer(k,nw,iaer).gt.qxiaer(k,nw,iaer))then |
---|
494 | qsiaer(k,nw,iaer)=qxiaer(k,nw,iaer) |
---|
495 | endif |
---|
496 | end do |
---|
497 | |
---|
498 | end do ! L_LEVELS |
---|
499 | end do ! naerkind |
---|
500 | |
---|
501 | !----------------------------------------------------------------------- |
---|
502 | ! Aerosol optical depths |
---|
503 | !----------------------------------------------------------------------- |
---|
504 | |
---|
505 | do iaer=1,naerkind ! a bug was here |
---|
506 | do k=0,nlayer-1 |
---|
507 | |
---|
508 | pweight=(pplay(ig,L_NLAYRAD-k)-pplev(ig,L_NLAYRAD-k+1))/ & |
---|
509 | (pplev(ig,L_NLAYRAD-k)-pplev(ig,L_NLAYRAD-k+1)) |
---|
510 | temp=aerosol(ig,L_NLAYRAD-k,iaer)/QREFvis3d(ig,L_NLAYRAD-k,iaer) |
---|
511 | tauaero(2*k+2,iaer)=max(temp*pweight,0.d0) |
---|
512 | tauaero(2*k+3,iaer)=max(temp-tauaero(2*k+2,iaer),0.d0) |
---|
513 | |
---|
514 | end do |
---|
515 | ! boundary conditions |
---|
516 | tauaero(1,iaer) = tauaero(2,iaer) |
---|
517 | !tauaero(1,iaer) = 0. |
---|
518 | !JL18 at time of testing, the two above conditions gave the same results bit for bit. |
---|
519 | |
---|
520 | end do ! naerkind |
---|
521 | |
---|
522 | ! Albedo and Emissivity. |
---|
523 | albi=1-emis(ig) ! Long Wave. |
---|
524 | DO nw=1,L_NSPECTV ! Short Wave loop. |
---|
525 | albv(nw)=albedo(ig,nw) |
---|
526 | ENDDO |
---|
527 | |
---|
528 | if ((ngrid.eq.1).and.(global1d)) then ! Fixed zenith angle 'szangle' in 1D simulations w/ globally-averaged sunlight. |
---|
529 | acosz = cos(pi*szangle/180.0) |
---|
530 | print*,'acosz=',acosz,', szangle=',szangle |
---|
531 | else |
---|
532 | acosz=mu0(ig) ! Cosine of sun incident angle : 3D simulations or local 1D simulations using latitude. |
---|
533 | endif |
---|
534 | |
---|
535 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
536 | !!! Note by JL13 : In the following, some indices were changed in the interpolations, |
---|
537 | !!! so that the model results are less dependent on the number of layers ! |
---|
538 | !!! |
---|
539 | !!! --- The older versions are commented with the comment !JL13index --- |
---|
540 | !!! |
---|
541 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
542 | |
---|
543 | |
---|
544 | !----------------------------------------------------------------------- |
---|
545 | ! Water vapour (to be generalised for other gases eventually ...) |
---|
546 | !----------------------------------------------------------------------- |
---|
547 | |
---|
548 | if(varactive)then |
---|
549 | |
---|
550 | i_var=igcm_h2o_vap |
---|
551 | do l=1,nlayer |
---|
552 | qvar(2*l) = pq(ig,nlayer+1-l,i_var) |
---|
553 | qvar(2*l+1) = pq(ig,nlayer+1-l,i_var) |
---|
554 | !JL13index qvar(2*l+1) = (pq(ig,nlayer+1-l,i_var)+pq(ig,max(nlayer-l,1),i_var))/2 |
---|
555 | !JL13index ! Average approximation as for temperature... |
---|
556 | end do |
---|
557 | qvar(1)=qvar(2) |
---|
558 | |
---|
559 | elseif(varfixed)then |
---|
560 | |
---|
561 | do l=1,nlayer ! Here we will assign fixed water vapour profiles globally. |
---|
562 | RH = satval * ((pplay(ig,l)/pplev(ig,1) - 0.02) / 0.98) |
---|
563 | if(RH.lt.0.0) RH=0.0 |
---|
564 | |
---|
565 | call Psat_water(pt(ig,l),pplay(ig,l),psat,qsat) |
---|
566 | |
---|
567 | !pq_temp(l) = qsat ! fully saturated everywhere |
---|
568 | pq_temp(l) = RH * qsat ! ~realistic profile (e.g. 80% saturation at ground) |
---|
569 | end do |
---|
570 | |
---|
571 | do l=1,nlayer |
---|
572 | qvar(2*l) = pq_temp(nlayer+1-l) |
---|
573 | qvar(2*l+1) = (pq_temp(nlayer+1-l)+pq_temp(max(nlayer-l,1)))/2 |
---|
574 | end do |
---|
575 | |
---|
576 | qvar(1)=qvar(2) |
---|
577 | |
---|
578 | ! Lowest layer of atmosphere |
---|
579 | RH = satval * (1 - 0.02) / 0.98 |
---|
580 | if(RH.lt.0.0) RH=0.0 |
---|
581 | |
---|
582 | qvar(2*nlayer+1)= RH * qsat ! ~realistic profile (e.g. 80% saturation at ground) |
---|
583 | |
---|
584 | else |
---|
585 | do k=1,L_LEVELS |
---|
586 | qvar(k) = 1.0D-7 |
---|
587 | end do |
---|
588 | end if ! varactive/varfixed |
---|
589 | |
---|
590 | if(.not.kastprof)then |
---|
591 | ! IMPORTANT: Now convert from kg/kg to mol/mol. |
---|
592 | do k=1,L_LEVELS |
---|
593 | qvar(k) = qvar(k)/(epsi+qvar(k)*(1.-epsi)) |
---|
594 | end do |
---|
595 | end if |
---|
596 | |
---|
597 | !----------------------------------------------------------------------- |
---|
598 | ! kcm mode only ! |
---|
599 | !----------------------------------------------------------------------- |
---|
600 | |
---|
601 | if(kastprof)then |
---|
602 | |
---|
603 | if(.not.global1d)then ! garde-fou/safeguard added by MT (to be removed in the future) |
---|
604 | write(*,*) 'You have to fix mu0, ' |
---|
605 | write(*,*) 'the cosinus of the solar angle' |
---|
606 | stop |
---|
607 | endif |
---|
608 | |
---|
609 | ! Initial values equivalent to mugaz. |
---|
610 | DO l=1,nlayer |
---|
611 | muvarrad(2*l) = mugaz |
---|
612 | muvarrad(2*l+1) = mugaz |
---|
613 | END DO |
---|
614 | |
---|
615 | if(ngasmx.gt.1)then |
---|
616 | |
---|
617 | DO l=1,nlayer |
---|
618 | muvarrad(2*l) = muvar(ig,nlayer+2-l) |
---|
619 | muvarrad(2*l+1) = (muvar(ig,nlayer+2-l) + & |
---|
620 | muvar(ig,max(nlayer+1-l,1)))/2 |
---|
621 | END DO |
---|
622 | |
---|
623 | muvarrad(1) = muvarrad(2) |
---|
624 | muvarrad(2*nlayer+1) = muvar(ig,1) |
---|
625 | |
---|
626 | print*,'Recalculating qvar with VARIABLE epsi for kastprof' |
---|
627 | print*,'Assumes that the variable gas is H2O!!!' |
---|
628 | print*,'Assumes that there is only one tracer' |
---|
629 | |
---|
630 | !i_var=igcm_h2o_vap |
---|
631 | i_var=1 |
---|
632 | |
---|
633 | if(nq.gt.1)then |
---|
634 | print*,'Need 1 tracer only to run kcm1d.e' |
---|
635 | stop |
---|
636 | endif |
---|
637 | |
---|
638 | do l=1,nlayer |
---|
639 | vtmp(l)=pq(ig,l,i_var)/(epsi+pq(ig,l,i_var)*(1.-epsi)) |
---|
640 | !vtmp(l)=pq(ig,l,i_var)*muvar(ig,l+1)/mH2O !JL to be changed |
---|
641 | end do |
---|
642 | |
---|
643 | do l=1,nlayer |
---|
644 | qvar(2*l) = vtmp(nlayer+1-l) |
---|
645 | qvar(2*l+1) = vtmp(nlayer+1-l) |
---|
646 | ! qvar(2*l+1) = ( vtmp(nlayer+1-l) + vtmp(max(nlayer-l,1)) )/2 |
---|
647 | end do |
---|
648 | qvar(1)=qvar(2) |
---|
649 | |
---|
650 | print*,'Warning: reducing qvar in callcorrk.' |
---|
651 | print*,'Temperature profile no longer consistent ', & |
---|
652 | 'with saturated H2O. qsat=',satval |
---|
653 | |
---|
654 | do k=1,L_LEVELS |
---|
655 | qvar(k) = qvar(k)*satval |
---|
656 | end do |
---|
657 | |
---|
658 | endif |
---|
659 | else ! if kastprof |
---|
660 | DO l=1,nlayer |
---|
661 | muvarrad(2*l) = muvar(ig,nlayer+2-l) |
---|
662 | muvarrad(2*l+1) = (muvar(ig,nlayer+2-l)+muvar(ig,max(nlayer+1-l,1)))/2 |
---|
663 | END DO |
---|
664 | |
---|
665 | muvarrad(1) = muvarrad(2) |
---|
666 | muvarrad(2*nlayer+1)=muvar(ig,1) |
---|
667 | endif ! if kastprof |
---|
668 | |
---|
669 | ! Keep values inside limits for which we have radiative transfer coefficients !!! |
---|
670 | if(L_REFVAR.gt.1)then ! (there was a bug here) |
---|
671 | do k=1,L_LEVELS |
---|
672 | if(qvar(k).lt.wrefvar(1))then |
---|
673 | qvar(k)=wrefvar(1)+1.0e-8 |
---|
674 | elseif(qvar(k).gt.wrefvar(L_REFVAR))then |
---|
675 | qvar(k)=wrefvar(L_REFVAR)-1.0e-8 |
---|
676 | endif |
---|
677 | end do |
---|
678 | endif |
---|
679 | |
---|
680 | !----------------------------------------------------------------------- |
---|
681 | ! Pressure and temperature |
---|
682 | !----------------------------------------------------------------------- |
---|
683 | |
---|
684 | DO l=1,nlayer |
---|
685 | plevrad(2*l) = pplay(ig,nlayer+1-l)/scalep |
---|
686 | plevrad(2*l+1) = pplev(ig,nlayer+1-l)/scalep |
---|
687 | tlevrad(2*l) = pt(ig,nlayer+1-l) |
---|
688 | tlevrad(2*l+1) = (pt(ig,nlayer+1-l)+pt(ig,max(nlayer-l,1)))/2 |
---|
689 | END DO |
---|
690 | |
---|
691 | plevrad(1) = 0. |
---|
692 | ! plevrad(2) = 0. !! JL18 enabling this line puts the radiative top at p=0 which was the idea before, but does not seem to perform best after all. |
---|
693 | |
---|
694 | tlevrad(1) = tlevrad(2) |
---|
695 | tlevrad(2*nlayer+1)=tsurf(ig) |
---|
696 | |
---|
697 | pmid(1) = pplay(ig,nlayer)/scalep |
---|
698 | pmid(2) = pmid(1) |
---|
699 | |
---|
700 | tmid(1) = tlevrad(2) |
---|
701 | tmid(2) = tmid(1) |
---|
702 | |
---|
703 | DO l=1,L_NLAYRAD-1 |
---|
704 | tmid(2*l+1) = tlevrad(2*l+1) |
---|
705 | tmid(2*l+2) = tlevrad(2*l+1) |
---|
706 | pmid(2*l+1) = plevrad(2*l+1) |
---|
707 | pmid(2*l+2) = plevrad(2*l+1) |
---|
708 | END DO |
---|
709 | pmid(L_LEVELS) = plevrad(L_LEVELS) |
---|
710 | tmid(L_LEVELS) = tlevrad(L_LEVELS) |
---|
711 | |
---|
712 | !!Alternative interpolation: |
---|
713 | ! pmid(3) = pmid(1) |
---|
714 | ! pmid(4) = pmid(1) |
---|
715 | ! tmid(3) = tmid(1) |
---|
716 | ! tmid(4) = tmid(1) |
---|
717 | ! DO l=2,L_NLAYRAD-1 |
---|
718 | ! tmid(2*l+1) = tlevrad(2*l) |
---|
719 | ! tmid(2*l+2) = tlevrad(2*l) |
---|
720 | ! pmid(2*l+1) = plevrad(2*l) |
---|
721 | ! pmid(2*l+2) = plevrad(2*l) |
---|
722 | ! END DO |
---|
723 | ! pmid(L_LEVELS) = plevrad(L_LEVELS-1) |
---|
724 | ! tmid(L_LEVELS) = tlevrad(L_LEVELS-1) |
---|
725 | |
---|
726 | ! Test for out-of-bounds pressure. |
---|
727 | if(plevrad(3).lt.pgasmin)then |
---|
728 | print*,'Minimum pressure is outside the radiative' |
---|
729 | print*,'transfer kmatrix bounds, exiting.' |
---|
730 | call abort |
---|
731 | elseif(plevrad(L_LEVELS).gt.pgasmax)then |
---|
732 | print*,'Maximum pressure is outside the radiative' |
---|
733 | print*,'transfer kmatrix bounds, exiting.' |
---|
734 | call abort |
---|
735 | endif |
---|
736 | |
---|
737 | ! Test for out-of-bounds temperature. |
---|
738 | do k=1,L_LEVELS |
---|
739 | if(tlevrad(k).lt.tgasmin)then |
---|
740 | print*,'Minimum temperature is outside the radiative' |
---|
741 | print*,'transfer kmatrix bounds' |
---|
742 | print*,"k=",k," tlevrad(k)=",tlevrad(k) |
---|
743 | print*,"tgasmin=",tgasmin |
---|
744 | if (strictboundcorrk) then |
---|
745 | call abort |
---|
746 | else |
---|
747 | print*,'***********************************************' |
---|
748 | print*,'we allow model to continue with tlevrad<tgasmin' |
---|
749 | print*,' ... we assume we know what you are doing ... ' |
---|
750 | print*,' ... but do not let this happen too often ... ' |
---|
751 | print*,'***********************************************' |
---|
752 | !tlevrad(k)=tgasmin ! Used in the source function ! |
---|
753 | endif |
---|
754 | elseif(tlevrad(k).gt.tgasmax)then |
---|
755 | print*,'Maximum temperature is outside the radiative' |
---|
756 | print*,'transfer kmatrix bounds, exiting.' |
---|
757 | print*,"k=",k," tlevrad(k)=",tlevrad(k) |
---|
758 | print*,"tgasmax=",tgasmax |
---|
759 | if (strictboundcorrk) then |
---|
760 | call abort |
---|
761 | else |
---|
762 | print*,'***********************************************' |
---|
763 | print*,'we allow model to continue with tlevrad<tgasmax' |
---|
764 | print*,' ... we assume we know what you are doing ... ' |
---|
765 | print*,' ... but do not let this happen too often ... ' |
---|
766 | print*,'***********************************************' |
---|
767 | !tlevrad(k)=tgasmax ! Used in the source function ! |
---|
768 | endif |
---|
769 | endif |
---|
770 | enddo |
---|
771 | do k=1,L_NLAYRAD+1 |
---|
772 | if(tmid(k).lt.tgasmin)then |
---|
773 | print*,'Minimum temperature is outside the radiative' |
---|
774 | print*,'transfer kmatrix bounds, exiting.' |
---|
775 | print*,"k=",k," tmid(k)=",tmid(k) |
---|
776 | print*,"tgasmin=",tgasmin |
---|
777 | if (strictboundcorrk) then |
---|
778 | call abort |
---|
779 | else |
---|
780 | print*,'***********************************************' |
---|
781 | print*,'we allow model to continue but with tmid=tgasmin' |
---|
782 | print*,' ... we assume we know what you are doing ... ' |
---|
783 | print*,' ... but do not let this happen too often ... ' |
---|
784 | print*,'***********************************************' |
---|
785 | tmid(k)=tgasmin |
---|
786 | endif |
---|
787 | elseif(tmid(k).gt.tgasmax)then |
---|
788 | print*,'Maximum temperature is outside the radiative' |
---|
789 | print*,'transfer kmatrix bounds, exiting.' |
---|
790 | print*,"k=",k," tmid(k)=",tmid(k) |
---|
791 | print*,"tgasmax=",tgasmax |
---|
792 | if (strictboundcorrk) then |
---|
793 | call abort |
---|
794 | else |
---|
795 | print*,'***********************************************' |
---|
796 | print*,'we allow model to continue but with tmid=tgasmin' |
---|
797 | print*,' ... we assume we know what you are doing ... ' |
---|
798 | print*,' ... but do not let this happen too often ... ' |
---|
799 | print*,'***********************************************' |
---|
800 | tmid(k)=tgasmax |
---|
801 | endif |
---|
802 | endif |
---|
803 | enddo |
---|
804 | |
---|
805 | !======================================================================= |
---|
806 | ! III. Calling the main radiative transfer subroutines |
---|
807 | !======================================================================= |
---|
808 | |
---|
809 | |
---|
810 | Cmk= 0.01 * 1.0 / (glat(ig) * mugaz * 1.672621e-27) ! q_main=1.0 assumed. |
---|
811 | glat_ig=glat(ig) |
---|
812 | |
---|
813 | !----------------------------------------------------------------------- |
---|
814 | ! Short Wave Part |
---|
815 | !----------------------------------------------------------------------- |
---|
816 | |
---|
817 | if(fract(ig) .ge. 1.0e-4) then ! Only during daylight. |
---|
818 | if((ngrid.eq.1).and.(global1d))then |
---|
819 | do nw=1,L_NSPECTV |
---|
820 | stel_fract(nw)= stel(nw)* 0.25 / acosz ! globally averaged = divide by 4, and we correct for solar zenith angle |
---|
821 | end do |
---|
822 | else |
---|
823 | do nw=1,L_NSPECTV |
---|
824 | stel_fract(nw)= stel(nw) * fract(ig) |
---|
825 | end do |
---|
826 | endif |
---|
827 | |
---|
828 | call optcv(dtauv,tauv,taucumv,plevrad, & |
---|
829 | qxvaer,qsvaer,gvaer,wbarv,cosbv,tauray,tauaero, & |
---|
830 | tmid,pmid,taugsurf,qvar,muvarrad) |
---|
831 | |
---|
832 | call sfluxv(dtauv,tauv,taucumv,albv,dwnv,wbarv,cosbv, & |
---|
833 | acosz,stel_fract, & |
---|
834 | nfluxtopv,fluxtopvdn,nfluxoutv_nu,nfluxgndv_nu, & |
---|
835 | fmnetv,fluxupv,fluxdnv,fzerov,taugsurf) |
---|
836 | |
---|
837 | else ! During the night, fluxes = 0. |
---|
838 | nfluxtopv = 0.0d0 |
---|
839 | fluxtopvdn = 0.0d0 |
---|
840 | nfluxoutv_nu(:) = 0.0d0 |
---|
841 | nfluxgndv_nu(:) = 0.0d0 |
---|
842 | do l=1,L_NLAYRAD |
---|
843 | fmnetv(l)=0.0d0 |
---|
844 | fluxupv(l)=0.0d0 |
---|
845 | fluxdnv(l)=0.0d0 |
---|
846 | end do |
---|
847 | end if |
---|
848 | |
---|
849 | |
---|
850 | ! Equivalent Albedo Calculation (for OUTPUT). MT2015 |
---|
851 | if(fract(ig) .ge. 1.0e-4) then ! equivalent albedo makes sense only during daylight. |
---|
852 | surface_stellar_flux=sum(nfluxgndv_nu(1:L_NSPECTV)) |
---|
853 | if(surface_stellar_flux .gt. 1.0e-3) then ! equivalent albedo makes sense only if the stellar flux received by the surface is positive. |
---|
854 | DO nw=1,L_NSPECTV |
---|
855 | albedo_temp(nw)=albedo(ig,nw)*nfluxgndv_nu(nw) |
---|
856 | ENDDO |
---|
857 | albedo_temp(1:L_NSPECTV)=albedo_temp(1:L_NSPECTV)/surface_stellar_flux |
---|
858 | albedo_equivalent(ig)=sum(albedo_temp(1:L_NSPECTV)) |
---|
859 | else |
---|
860 | albedo_equivalent(ig)=0.0 ! Spectrally Integrated Albedo not defined for non-irradiated grid points. So we arbitrary set the equivalent albedo to 0. |
---|
861 | endif |
---|
862 | else |
---|
863 | albedo_equivalent(ig)=0.0 ! Spectrally Integrated Albedo not defined for non-irradiated grid points. So we arbitrary set the equivalent albedo to 0. |
---|
864 | endif |
---|
865 | |
---|
866 | |
---|
867 | !----------------------------------------------------------------------- |
---|
868 | ! Long Wave Part |
---|
869 | !----------------------------------------------------------------------- |
---|
870 | |
---|
871 | call optci(plevrad,tlevrad,dtaui,taucumi, & |
---|
872 | qxiaer,qsiaer,giaer,cosbi,wbari,tauaero,tmid,pmid, & |
---|
873 | taugsurfi,qvar,muvarrad) |
---|
874 | |
---|
875 | call sfluxi(plevrad,tlevrad,dtaui,taucumi,ubari,albi, & |
---|
876 | wnoi,dwni,cosbi,wbari,nfluxtopi,nfluxtopi_nu, & |
---|
877 | fmneti,fluxupi,fluxdni,fluxupi_nu,fzeroi,taugsurfi) |
---|
878 | |
---|
879 | !----------------------------------------------------------------------- |
---|
880 | ! Transformation of the correlated-k code outputs |
---|
881 | ! (into dtlw, dtsw, fluxsurf_lw, fluxsurf_sw, fluxtop_lw, fluxtop_sw) |
---|
882 | |
---|
883 | ! Flux incident at the top of the atmosphere |
---|
884 | fluxtop_dn(ig)=fluxtopvdn |
---|
885 | |
---|
886 | fluxtop_lw(ig) = real(nfluxtopi) |
---|
887 | fluxabs_sw(ig) = real(-nfluxtopv) |
---|
888 | fluxsurf_lw(ig) = real(fluxdni(L_NLAYRAD)) |
---|
889 | fluxsurf_sw(ig) = real(fluxdnv(L_NLAYRAD)) |
---|
890 | |
---|
891 | ! Flux absorbed by the surface. By MT2015. |
---|
892 | fluxsurfabs_sw(ig) = fluxsurf_sw(ig)*(1.-albedo_equivalent(ig)) |
---|
893 | |
---|
894 | if(fluxtop_dn(ig).lt.0.0)then |
---|
895 | print*,'Achtung! fluxtop_dn has lost the plot!' |
---|
896 | print*,'fluxtop_dn=',fluxtop_dn(ig) |
---|
897 | print*,'acosz=',acosz |
---|
898 | print*,'aerosol=',aerosol(ig,:,:) |
---|
899 | print*,'temp= ',pt(ig,:) |
---|
900 | print*,'pplay= ',pplay(ig,:) |
---|
901 | call abort |
---|
902 | endif |
---|
903 | |
---|
904 | ! Spectral output, for exoplanet observational comparison |
---|
905 | if(specOLR)then |
---|
906 | do nw=1,L_NSPECTI |
---|
907 | OLR_nu(ig,nw)=nfluxtopi_nu(nw)/DWNI(nw) !JL Normalize to the bandwidth |
---|
908 | end do |
---|
909 | do nw=1,L_NSPECTV |
---|
910 | !GSR_nu(ig,nw)=nfluxgndv_nu(nw) |
---|
911 | OSR_nu(ig,nw)=nfluxoutv_nu(nw)/DWNV(nw) !JL Normalize to the bandwidth |
---|
912 | end do |
---|
913 | endif |
---|
914 | |
---|
915 | ! Finally, the heating rates |
---|
916 | |
---|
917 | DO l=2,L_NLAYRAD |
---|
918 | dtsw(ig,L_NLAYRAD+1-l)=(fmnetv(l)-fmnetv(l-1)) & |
---|
919 | *glat(ig)/(cpp*scalep*(plevrad(2*l+1)-plevrad(2*l-1))) |
---|
920 | dtlw(ig,L_NLAYRAD+1-l)=(fmneti(l)-fmneti(l-1)) & |
---|
921 | *glat(ig)/(cpp*scalep*(plevrad(2*l+1)-plevrad(2*l-1))) |
---|
922 | END DO |
---|
923 | |
---|
924 | ! These are values at top of atmosphere |
---|
925 | dtsw(ig,L_NLAYRAD)=(fmnetv(1)-nfluxtopv) & |
---|
926 | *glat(ig)/(cpp*scalep*(plevrad(3)-plevrad(2))) |
---|
927 | dtlw(ig,L_NLAYRAD)=(fmneti(1)-nfluxtopi) & |
---|
928 | *glat(ig)/(cpp*scalep*(plevrad(3)-plevrad(2))) |
---|
929 | |
---|
930 | ! Optical thickness diagnostics (added by JVO) |
---|
931 | if (diagdtau) then |
---|
932 | do l=1,L_NLAYRAD |
---|
933 | do nw=1,L_NSPECTV |
---|
934 | int_dtauv(ig,l,nw) = 0.0d0 |
---|
935 | DO k=1,L_NGAUSS |
---|
936 | ! Output exp(-tau) because gweight ponderates exp and not tau itself |
---|
937 | int_dtauv(ig,l,nw)= int_dtauv(ig,l,nw) + exp(-dtauv(l,nw,k))*gweight(k) |
---|
938 | ENDDO |
---|
939 | enddo |
---|
940 | do nw=1,L_NSPECTI |
---|
941 | int_dtaui(ig,l,nw) = 0.0d0 |
---|
942 | DO k=1,L_NGAUSS |
---|
943 | ! Output exp(-tau) because gweight ponderates exp and not tau itself |
---|
944 | int_dtaui(ig,l,nw)= int_dtaui(ig,l,nw) + exp(-dtaui(l,nw,k))*gweight(k) |
---|
945 | ENDDO |
---|
946 | enddo |
---|
947 | enddo |
---|
948 | endif |
---|
949 | |
---|
950 | |
---|
951 | !----------------------------------------------------------------------- |
---|
952 | end do ! End of big loop over every GCM column. |
---|
953 | !----------------------------------------------------------------------- |
---|
954 | |
---|
955 | |
---|
956 | |
---|
957 | !----------------------------------------------------------------------- |
---|
958 | ! Additional diagnostics |
---|
959 | !----------------------------------------------------------------------- |
---|
960 | |
---|
961 | ! IR spectral output, for exoplanet observational comparison |
---|
962 | if(lastcall.and.(ngrid.eq.1))then ! could disable the 1D output, they are in the diagfi and diagspec... JL12 |
---|
963 | |
---|
964 | print*,'Saving scalar quantities in surf_vals.out...' |
---|
965 | print*,'psurf = ', pplev(1,1),' Pa' |
---|
966 | open(116,file='surf_vals.out') |
---|
967 | write(116,*) tsurf(1),pplev(1,1),fluxtop_dn(1), & |
---|
968 | real(-nfluxtopv),real(nfluxtopi) |
---|
969 | close(116) |
---|
970 | |
---|
971 | |
---|
972 | ! USEFUL COMMENT - Do Not Remove. |
---|
973 | ! |
---|
974 | ! if(specOLR)then |
---|
975 | ! open(117,file='OLRnu.out') |
---|
976 | ! do nw=1,L_NSPECTI |
---|
977 | ! write(117,*) OLR_nu(1,nw) |
---|
978 | ! enddo |
---|
979 | ! close(117) |
---|
980 | ! |
---|
981 | ! open(127,file='OSRnu.out') |
---|
982 | ! do nw=1,L_NSPECTV |
---|
983 | ! write(127,*) OSR_nu(1,nw) |
---|
984 | ! enddo |
---|
985 | ! close(127) |
---|
986 | ! endif |
---|
987 | |
---|
988 | ! OLR vs altitude: do it as a .txt file. |
---|
989 | OLRz=.false. |
---|
990 | if(OLRz)then |
---|
991 | print*,'saving IR vertical flux for OLRz...' |
---|
992 | open(118,file='OLRz_plevs.out') |
---|
993 | open(119,file='OLRz.out') |
---|
994 | do l=1,L_NLAYRAD |
---|
995 | write(118,*) plevrad(2*l) |
---|
996 | do nw=1,L_NSPECTI |
---|
997 | write(119,*) fluxupi_nu(l,nw) |
---|
998 | enddo |
---|
999 | enddo |
---|
1000 | close(118) |
---|
1001 | close(119) |
---|
1002 | endif |
---|
1003 | |
---|
1004 | endif |
---|
1005 | |
---|
1006 | ! See physiq.F for explanations about CLFvarying. This is temporary. |
---|
1007 | if (lastcall .and. .not.CLFvarying) then |
---|
1008 | IF( ALLOCATED( gasi ) ) DEALLOCATE( gasi ) |
---|
1009 | IF( ALLOCATED( gasv ) ) DEALLOCATE( gasv ) |
---|
1010 | !$OMP BARRIER |
---|
1011 | !$OMP MASTER |
---|
1012 | IF( ALLOCATED( pgasref ) ) DEALLOCATE( pgasref ) |
---|
1013 | IF( ALLOCATED( tgasref ) ) DEALLOCATE( tgasref ) |
---|
1014 | IF( ALLOCATED( wrefvar ) ) DEALLOCATE( wrefvar ) |
---|
1015 | IF( ALLOCATED( pfgasref ) ) DEALLOCATE( pfgasref ) |
---|
1016 | IF( ALLOCATED( gweight ) ) DEALLOCATE( gweight ) |
---|
1017 | !$OMP END MASTER |
---|
1018 | !$OMP BARRIER |
---|
1019 | IF ( ALLOCATED(reffrad)) DEALLOCATE(reffrad) |
---|
1020 | IF ( ALLOCATED(nueffrad)) DEALLOCATE(nueffrad) |
---|
1021 | endif |
---|
1022 | |
---|
1023 | |
---|
1024 | end subroutine callcorrk |
---|
1025 | |
---|
1026 | END MODULE callcorrk_mod |
---|