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 | zzlay,zzlev,tsurf,fract,dist_star,dtau_aer,muvar, & |
---|
10 | dtlw,dtsw,fluxsurf_lw, & |
---|
11 | fluxsurf_sw,fluxsurfabs_sw,fluxtop_lw, & |
---|
12 | fluxabs_sw,fluxtop_dn, & |
---|
13 | OLR_nu,OSR_nu,GSR_nu, & |
---|
14 | int_dtaui,int_dtauv, & |
---|
15 | tau_col,firstcall,lastcall) |
---|
16 | |
---|
17 | use mod_phys_lmdz_para, only : is_master |
---|
18 | use radinc_h, only: L_NSPECTV, L_NSPECTI, naerkind, banddir, corrkdir,& |
---|
19 | L_LEVELS, L_NGAUSS, L_NLEVRAD, L_NLAYRAD, L_REFVAR |
---|
20 | use radcommon_h, only: wrefvar, Cmk, fzeroi, fzerov, gasi, gasv, & |
---|
21 | glat_ig, gweight, pfgasref, pgasmax, pgasmin, & |
---|
22 | pgasref, tgasmax, tgasmin, tgasref, scalep, & |
---|
23 | ubari, wnoi, stellarf, glat, dwnv, dwni, tauray |
---|
24 | use datafile_mod, only: datadir |
---|
25 | use ioipsl_getin_p_mod, only: getin_p |
---|
26 | use gases_h, only: ngasmx |
---|
27 | use radii_mod, only : su_aer_radii, haze_reffrad_fix |
---|
28 | use aerosol_mod, only : iaero_haze |
---|
29 | use aeropacity_mod, only: aeropacity |
---|
30 | use aeroptproperties_mod, only: aeroptproperties |
---|
31 | use tracer_h, only: igcm_ch4_gas,igcm_n2,mmol |
---|
32 | use comcstfi_mod, only: pi, mugaz, cpp, r, g |
---|
33 | use callkeys_mod, only: diurnal,tracer,varfixed, & |
---|
34 | kastprof,strictboundcorrk,specOLR, & |
---|
35 | tplanckmin,tplanckmax,global1d, & |
---|
36 | optichaze,haze_radproffix,& |
---|
37 | methane,carbox,cooling,nlte,strobel,& |
---|
38 | ch4fix,vmrch4_proffix,vmrch4fix,& |
---|
39 | callmufi |
---|
40 | use optcv_mod, only: optcv |
---|
41 | use optci_mod, only: optci |
---|
42 | use sfluxi_mod, only: sfluxi |
---|
43 | use sfluxv_mod, only: sfluxv |
---|
44 | use recombin_corrk_mod, only: corrk_recombin, call_recombin |
---|
45 | use planetwide_mod, only: planetwide_maxval, planetwide_minval |
---|
46 | use radcommon_h, only: wavev,wavei |
---|
47 | use mp2m_diagnostics |
---|
48 | implicit none |
---|
49 | |
---|
50 | !================================================================== |
---|
51 | ! |
---|
52 | ! Purpose |
---|
53 | ! ------- |
---|
54 | ! Solve the radiative transfer using the correlated-k method for |
---|
55 | ! the gaseous absorption and the Toon et al. (1989) method for |
---|
56 | ! scatttering due to aerosols. |
---|
57 | ! |
---|
58 | ! Authors |
---|
59 | ! ------- |
---|
60 | ! Emmanuel 01/2001, Forget 09/2001 |
---|
61 | ! Robin Wordsworth (2009) |
---|
62 | ! |
---|
63 | !================================================================== |
---|
64 | |
---|
65 | !----------------------------------------------------------------------- |
---|
66 | ! Declaration of the arguments (INPUT - OUTPUT) on the LMD GCM grid |
---|
67 | ! Layer #1 is the layer near the ground. |
---|
68 | ! Layer #nlayer is the layer at the top. |
---|
69 | !----------------------------------------------------------------------- |
---|
70 | |
---|
71 | |
---|
72 | ! INPUT |
---|
73 | INTEGER,INTENT(IN) :: ngrid ! Number of atmospheric columns. |
---|
74 | INTEGER,INTENT(IN) :: nlayer ! Number of atmospheric layers. |
---|
75 | REAL,INTENT(IN) :: pq(ngrid,nlayer,nq) ! Tracers (kg/kg_of_air). |
---|
76 | INTEGER,INTENT(IN) :: nq ! Number of tracers. |
---|
77 | REAL,INTENT(IN) :: qsurf(ngrid,nq) ! Tracers on surface (kg.m-2). |
---|
78 | REAL,INTENT(IN) :: albedo(ngrid,L_NSPECTV) ! Spectral Short Wavelengths Albedo. By MT2015 |
---|
79 | REAL,INTENT(IN) :: emis(ngrid) ! Long Wave emissivity. |
---|
80 | REAL,INTENT(IN) :: mu0(ngrid) ! Cosine of sun incident angle. |
---|
81 | REAL,INTENT(IN) :: pplev(ngrid,nlayer+1) ! Inter-layer pressure (Pa). |
---|
82 | REAL,INTENT(IN) :: pplay(ngrid,nlayer) ! Mid-layer pressure (Pa). |
---|
83 | REAL,INTENT(IN) :: pt(ngrid,nlayer) ! Air temperature (K). |
---|
84 | REAL,INTENT(IN) :: zzlay(ngrid,nlayer) ! Mid-layer altitude |
---|
85 | REAL,INTENT(IN) :: zzlev(ngrid,nlayer) ! Altitude at the layer boundaries. |
---|
86 | REAL,INTENT(IN) :: tsurf(ngrid) ! Surface temperature (K). |
---|
87 | REAL,INTENT(IN) :: fract(ngrid) ! Fraction of day. |
---|
88 | REAL,INTENT(IN) :: dist_star ! Distance star-planet (AU). |
---|
89 | REAL,INTENT(IN) :: muvar(ngrid,nlayer+1) |
---|
90 | logical,intent(in) :: firstcall ! Signals first call to physics. |
---|
91 | logical,intent(in) :: lastcall ! Signals last call to physics. |
---|
92 | |
---|
93 | ! OUTPUT |
---|
94 | REAL,INTENT(OUT) :: dtau_aer(ngrid,nlayer,naerkind) ! Aerosol tau at reference wavelenght. |
---|
95 | REAL,INTENT(OUT) :: dtlw(ngrid,nlayer) ! Heating rate (K/s) due to LW radiation. |
---|
96 | REAL,INTENT(OUT) :: dtsw(ngrid,nlayer) ! Heating rate (K/s) due to SW radiation. |
---|
97 | REAL,INTENT(OUT) :: fluxsurf_lw(ngrid) ! Incident LW flux to surf (W/m2). |
---|
98 | REAL,INTENT(OUT) :: fluxsurf_sw(ngrid) ! Incident SW flux to surf (W/m2) |
---|
99 | REAL,INTENT(OUT) :: fluxsurfabs_sw(ngrid) ! Absorbed SW flux by the surface (W/m2). By MT2015. |
---|
100 | REAL,INTENT(OUT) :: fluxtop_lw(ngrid) ! Outgoing LW flux to space (W/m2). |
---|
101 | REAL,INTENT(OUT) :: fluxabs_sw(ngrid) ! SW flux absorbed by the planet (W/m2). |
---|
102 | REAL,INTENT(OUT) :: fluxtop_dn(ngrid) ! Incident top of atmosphere SW flux (W/m2). |
---|
103 | REAL,INTENT(OUT) :: OLR_nu(ngrid,L_NSPECTI) ! Outgoing LW radiation in each band (Normalized to the band width (W/m2/cm-1). |
---|
104 | REAL,INTENT(OUT) :: OSR_nu(ngrid,L_NSPECTV) ! Outgoing SW radiation in each band (Normalized to the band width (W/m2/cm-1). |
---|
105 | REAL,INTENT(OUT) :: GSR_nu(ngrid,L_NSPECTV) ! Surface SW radiation in each band (Normalized to the band width (W/m2/cm-1). |
---|
106 | REAL,INTENT(OUT) :: int_dtaui(ngrid,nlayer,L_NSPECTI) ! VI optical thickness of layers within narrowbands for diags (). |
---|
107 | REAL,INTENT(OUT) :: int_dtauv(ngrid,nlayer,L_NSPECTV) ! IR optical thickness of layers within narrowbands for diags (). |
---|
108 | REAL,INTENT(OUT) :: tau_col(ngrid) ! Diagnostic from aeropacity. |
---|
109 | REAL,INTENT(OUT) :: albedo_equivalent(ngrid) ! Spectrally Integrated Albedo. For Diagnostic. By MT2015 |
---|
110 | |
---|
111 | |
---|
112 | ! Globally varying aerosol optical properties on GCM grid ; not needed everywhere so not in radcommon_h. |
---|
113 | ! made "save" variables so they are allocated once in for all, not because |
---|
114 | ! the values need be saved from a time step to the next |
---|
115 | REAL,SAVE,ALLOCATABLE :: QVISsQREF3d(:,:,:,:) |
---|
116 | REAL,SAVE,ALLOCATABLE :: omegaVIS3d(:,:,:,:) |
---|
117 | REAL,SAVE,ALLOCATABLE :: gVIS3d(:,:,:,:) |
---|
118 | !$OMP THREADPRIVATE(QVISsQREF3d,omegaVIS3d,gVIS3d) |
---|
119 | REAL,SAVE,ALLOCATABLE :: QIRsQREF3d(:,:,:,:) |
---|
120 | REAL,SAVE,ALLOCATABLE :: omegaIR3d(:,:,:,:) |
---|
121 | REAL,SAVE,ALLOCATABLE :: gIR3d(:,:,:,:) |
---|
122 | !$OMP THREADPRIVATE(QIRsQREF3d,omegaIR3d,gIR3d) |
---|
123 | |
---|
124 | ! REAL :: omegaREFvis3d(ngrid,nlayer,naerkind) |
---|
125 | ! REAL :: omegaREFir3d(ngrid,nlayer,naerkind) ! not sure of the point of these... |
---|
126 | |
---|
127 | REAL,ALLOCATABLE,SAVE :: reffrad(:,:,:) ! aerosol effective radius (m) |
---|
128 | REAL,ALLOCATABLE,SAVE :: nueffrad(:,:,:) ! aerosol effective variance |
---|
129 | !$OMP THREADPRIVATE(reffrad,nueffrad) |
---|
130 | |
---|
131 | !----------------------------------------------------------------------- |
---|
132 | ! Declaration of the variables required by correlated-k subroutines |
---|
133 | ! Numbered from top to bottom (unlike in the GCM) |
---|
134 | !----------------------------------------------------------------------- |
---|
135 | |
---|
136 | REAL*8 tmid(L_LEVELS),pmid(L_LEVELS) |
---|
137 | REAL*8 tlevrad(L_LEVELS),plevrad(L_LEVELS) |
---|
138 | |
---|
139 | ! Optical values for the optci/cv subroutines |
---|
140 | REAL*8 stel(L_NSPECTV),stel_fract(L_NSPECTV) |
---|
141 | ! NB: Arrays below are "save" to avoid reallocating them at every call |
---|
142 | ! not because their content needs be reused from call to the next |
---|
143 | REAL*8,allocatable,save :: dtaui(:,:,:) |
---|
144 | REAL*8,allocatable,save :: dtauv(:,:,:) |
---|
145 | REAL*8,allocatable,save :: cosbv(:,:,:) |
---|
146 | REAL*8,allocatable,save :: cosbi(:,:,:) |
---|
147 | REAL*8,allocatable,save :: wbari(:,:,:) |
---|
148 | REAL*8,allocatable,save :: wbarv(:,:,:) |
---|
149 | !$OMP THREADPRIVATE(dtaui,dtauv,cosbv,cosbi,wbari,wbarv) |
---|
150 | REAL*8,allocatable,save :: tauv(:,:,:) |
---|
151 | REAL*8,allocatable,save :: taucumv(:,:,:) |
---|
152 | REAL*8,allocatable,save :: taucumi(:,:,:) |
---|
153 | !$OMP THREADPRIVATE(tauv,taucumv,taucumi) |
---|
154 | REAL*8,allocatable,save :: tauaero(:,:) |
---|
155 | !$OMP THREADPRIVATE(tauaero) |
---|
156 | REAL*8 nfluxtopv,nfluxtopi,nfluxtop,fluxtopvdn |
---|
157 | real*8 nfluxtopv_nu(L_NSPECTV) |
---|
158 | REAL*8 nfluxoutv_nu(L_NSPECTV) ! Outgoing band-resolved VI flux at TOA (W/m2). |
---|
159 | REAL*8 nfluxtopi_nu(L_NSPECTI) ! Net band-resolved IR flux at TOA (W/m2). |
---|
160 | REAL*8 fluxupi_nu(L_NLAYRAD,L_NSPECTI) ! For 1D diagnostic. |
---|
161 | REAL*8 fmneti(L_NLAYRAD),fmnetv(L_NLAYRAD) |
---|
162 | real*8 fmneti_nu(L_NLAYRAD,L_NSPECTI) ! |
---|
163 | real*8 fmnetv_nu(L_NLAYRAD,L_NSPECTV) ! |
---|
164 | REAL*8 fluxupv(L_NLAYRAD),fluxupi(L_NLAYRAD) |
---|
165 | REAL*8 fluxdnv(L_NLAYRAD),fluxdni(L_NLAYRAD) |
---|
166 | REAL*8 albi,acosz |
---|
167 | REAL*8 albv(L_NSPECTV) ! Spectral Visible Albedo. |
---|
168 | |
---|
169 | INTEGER ig,l,k,nw,iaer,iq |
---|
170 | |
---|
171 | real*8,allocatable,save :: taugsurf(:,:) |
---|
172 | real*8,allocatable,save :: taugsurfi(:,:) |
---|
173 | !$OMP THREADPRIVATE(taugsurf,taugsurfi) |
---|
174 | real*8 qvar(L_LEVELS) ! Mixing ratio of variable component (mol/mol). index 1 is the top of the atmosphere, index L_LEVELS is the bottom |
---|
175 | |
---|
176 | ! Local aerosol optical properties for each column on RADIATIVE grid. |
---|
177 | real*8,save,allocatable :: QXVAER(:,:,:) ! Extinction coeff (QVISsQREF*QREFvis) |
---|
178 | real*8,save,allocatable :: QSVAER(:,:,:) |
---|
179 | real*8,save,allocatable :: GVAER(:,:,:) |
---|
180 | real*8,save,allocatable :: QXIAER(:,:,:) ! Extinction coeff (QIRsQREF*QREFir) |
---|
181 | real*8,save,allocatable :: QSIAER(:,:,:) |
---|
182 | real*8,save,allocatable :: GIAER(:,:,:) |
---|
183 | !$OMP THREADPRIVATE(QXVAER,QSVAER,GVAER,QXIAER,QSIAER,GIAER) |
---|
184 | real, dimension(:,:,:), save, allocatable :: QREFvis3d |
---|
185 | real, dimension(:,:,:), save, allocatable :: QREFir3d |
---|
186 | !$OMP THREADPRIVATE(QREFvis3d,QREFir3d) |
---|
187 | |
---|
188 | ! Miscellaneous : |
---|
189 | real*8 temp,temp1,temp2,pweight,sig |
---|
190 | character(len=10) :: tmp1 |
---|
191 | character(len=10) :: tmp2 |
---|
192 | character(len=100) :: message |
---|
193 | character(len=10),parameter :: subname="callcorrk" |
---|
194 | |
---|
195 | logical OLRz |
---|
196 | real*8 NFLUXGNDV_nu(L_NSPECTV) |
---|
197 | |
---|
198 | ! Included by RW for runaway greenhouse 1D study. |
---|
199 | real vtmp(nlayer) |
---|
200 | REAL*8 muvarrad(L_LEVELS) |
---|
201 | |
---|
202 | ! Included by MT for albedo calculations. |
---|
203 | REAL*8 albedo_temp(L_NSPECTV) ! For equivalent albedo calculation. |
---|
204 | REAL*8 surface_stellar_flux ! Stellar flux reaching the surface. Useful for equivalent albedo calculation. |
---|
205 | |
---|
206 | ! NLTE factor for CH4 |
---|
207 | real eps_nlte_sw23(ngrid,nlayer) ! CH4 NLTE efficiency factor for zdtsw |
---|
208 | real eps_nlte_sw33(ngrid,nlayer) ! CH4 NLTE efficiency factor for zdtsw |
---|
209 | real eps_nlte_lw(ngrid,nlayer) ! CH4 NLTE efficiency factor for zdtsw |
---|
210 | integer Nfine,ifine |
---|
211 | parameter(Nfine=701) |
---|
212 | real,save :: levdat(Nfine),vmrdat(Nfine) |
---|
213 | REAL dtlw_hcn_c2h2(ngrid, nlayer) ! cooling rate (K/s) due to C2H2/HCN (diagnostic) |
---|
214 | real :: vmrch4(ngrid,nlayer) ! vmr ch4 from vmrch4_proffix |
---|
215 | |
---|
216 | REAL dtlw_nu(nlayer,L_NSPECTI) ! heating rate (K/s) due to LW in spectral bands |
---|
217 | REAL dtsw_nu(nlayer,L_NSPECTV) ! heating rate (K/s) due to SW in spectral bands |
---|
218 | |
---|
219 | ! local variable |
---|
220 | REAL, save :: dpp ! intermediate |
---|
221 | !$OMP THREADPRIVATE(dpp) |
---|
222 | |
---|
223 | integer ok ! status (returned by NetCDF functions) |
---|
224 | |
---|
225 | REAL :: maxvalue,minvalue |
---|
226 | |
---|
227 | !=============================================================== |
---|
228 | ! I.a Initialization on first call |
---|
229 | !=============================================================== |
---|
230 | |
---|
231 | |
---|
232 | if(firstcall) then |
---|
233 | |
---|
234 | ! test on allocated necessary because of CLFvarying (two calls to callcorrk in physiq) |
---|
235 | if(.not.allocated(QVISsQREF3d)) then |
---|
236 | allocate(QVISsQREF3d(ngrid,nlayer,L_NSPECTV,naerkind)) |
---|
237 | endif |
---|
238 | if(.not.allocated(omegaVIS3d)) then |
---|
239 | allocate(omegaVIS3d(ngrid,nlayer,L_NSPECTV,naerkind)) |
---|
240 | endif |
---|
241 | if(.not.allocated(gVIS3d)) then |
---|
242 | allocate(gVIS3d(ngrid,nlayer,L_NSPECTV,naerkind)) |
---|
243 | endif |
---|
244 | if (.not.allocated(QIRsQREF3d)) then |
---|
245 | allocate(QIRsQREF3d(ngrid,nlayer,L_NSPECTI,naerkind)) |
---|
246 | endif |
---|
247 | if (.not.allocated(omegaIR3d)) then |
---|
248 | allocate(omegaIR3d(ngrid,nlayer,L_NSPECTI,naerkind)) |
---|
249 | endif |
---|
250 | if (.not.allocated(gIR3d)) then |
---|
251 | allocate(gIR3d(ngrid,nlayer,L_NSPECTI,naerkind)) |
---|
252 | endif |
---|
253 | if (.not.allocated(tauaero)) then |
---|
254 | allocate(tauaero(L_LEVELS,naerkind)) |
---|
255 | endif |
---|
256 | |
---|
257 | if(.not.allocated(QXVAER)) then |
---|
258 | allocate(QXVAER(L_LEVELS,L_NSPECTV,naerkind), stat=ok) |
---|
259 | if (ok /= 0) then |
---|
260 | write(*,*) "memory allocation failed for QXVAER!" |
---|
261 | call abort_physic(subname,'allocation failure for QXVAER',1) |
---|
262 | endif |
---|
263 | endif |
---|
264 | if(.not.allocated(QSVAER)) then |
---|
265 | allocate(QSVAER(L_LEVELS,L_NSPECTV,naerkind), stat=ok) |
---|
266 | if (ok /= 0) then |
---|
267 | write(*,*) "memory allocation failed for QSVAER!" |
---|
268 | call abort_physic(subname,'allocation failure for QSVAER',1) |
---|
269 | endif |
---|
270 | endif |
---|
271 | if(.not.allocated(GVAER)) then |
---|
272 | allocate(GVAER(L_LEVELS,L_NSPECTV,naerkind), stat=ok) |
---|
273 | if (ok /= 0) then |
---|
274 | write(*,*) "memory allocation failed for GVAER!" |
---|
275 | call abort_physic(subname,'allocation failure for GVAER',1) |
---|
276 | endif |
---|
277 | endif |
---|
278 | if(.not.allocated(QXIAER)) then |
---|
279 | allocate(QXIAER(L_LEVELS,L_NSPECTI,naerkind), stat=ok) |
---|
280 | if (ok /= 0) then |
---|
281 | write(*,*) "memory allocation failed for QXIAER!" |
---|
282 | call abort_physic(subname,'allocation failure for QXIAER',1) |
---|
283 | endif |
---|
284 | endif |
---|
285 | if(.not.allocated(QSIAER)) then |
---|
286 | allocate(QSIAER(L_LEVELS,L_NSPECTI,naerkind), stat=ok) |
---|
287 | if (ok /= 0) then |
---|
288 | write(*,*) "memory allocation failed for QSIAER!" |
---|
289 | call abort_physic(subname,'allocation failure for QSIAER',1) |
---|
290 | endif |
---|
291 | endif |
---|
292 | if(.not.allocated(GIAER)) then |
---|
293 | allocate(GIAER(L_LEVELS,L_NSPECTI,naerkind), stat=ok) |
---|
294 | if (ok /= 0) then |
---|
295 | write(*,*) "memory allocation failed for GIAER!" |
---|
296 | call abort_physic(subname,'allocation failure for GIAER',1) |
---|
297 | endif |
---|
298 | endif |
---|
299 | |
---|
300 | !!! ALLOCATED instances are necessary because of CLFvarying (strategy to call callcorrk twice in physiq...) |
---|
301 | IF(.not.ALLOCATED(QREFvis3d))THEN |
---|
302 | ALLOCATE(QREFvis3d(ngrid,nlayer,naerkind), stat=ok) |
---|
303 | IF (ok/=0) THEN |
---|
304 | write(*,*) "memory allocation failed for QREFvis3d!" |
---|
305 | call abort_physic(subname,'allocation failure for QREFvis3d',1) |
---|
306 | ENDIF |
---|
307 | ENDIF |
---|
308 | IF(.not.ALLOCATED(QREFir3d)) THEN |
---|
309 | ALLOCATE(QREFir3d(ngrid,nlayer,naerkind), stat=ok) |
---|
310 | IF (ok/=0) THEN |
---|
311 | write(*,*) "memory allocation failed for QREFir3d!" |
---|
312 | call abort_physic(subname,'allocation failure for QREFir3d',1) |
---|
313 | ENDIF |
---|
314 | ENDIF |
---|
315 | ! Effective radius and variance of the aerosols |
---|
316 | IF(.not.ALLOCATED(reffrad)) THEN |
---|
317 | allocate(reffrad(ngrid,nlayer,naerkind), stat=ok) |
---|
318 | IF (ok/=0) THEN |
---|
319 | write(*,*) "memory allocation failed for reffrad!" |
---|
320 | call abort_physic(subname,'allocation failure for reffrad',1) |
---|
321 | ENDIF |
---|
322 | ENDIF |
---|
323 | IF(.not.ALLOCATED(nueffrad)) THEN |
---|
324 | allocate(nueffrad(ngrid,nlayer,naerkind), stat=ok) |
---|
325 | IF (ok/=0) THEN |
---|
326 | write(*,*) "memory allocation failed for nueffrad!" |
---|
327 | call abort_physic(subname,'allocation failure for nueffrad',1) |
---|
328 | ENDIF |
---|
329 | ENDIF |
---|
330 | |
---|
331 | if (is_master) call system('rm -f surf_vals_long.out') |
---|
332 | |
---|
333 | !-------------------------------------------------- |
---|
334 | ! Set up correlated k |
---|
335 | !-------------------------------------------------- |
---|
336 | |
---|
337 | !this block is now done at firstcall of physiq_mod |
---|
338 | ! print*, "callcorrk: Correlated-k data base folder:",trim(datadir) |
---|
339 | ! call getin_p("corrkdir",corrkdir) |
---|
340 | ! print*, "corrkdir = ",corrkdir |
---|
341 | ! write( tmp1, '(i3)' ) L_NSPECTI |
---|
342 | ! write( tmp2, '(i3)' ) L_NSPECTV |
---|
343 | ! banddir=trim(adjustl(tmp1))//'x'//trim(adjustl(tmp2)) |
---|
344 | ! banddir=trim(adjustl(corrkdir))//'/'//trim(adjustl(banddir)) |
---|
345 | |
---|
346 | ! call setspi ! Basic infrared properties. |
---|
347 | ! call setspv ! Basic visible properties. |
---|
348 | ! call sugas_corrk ! Set up gaseous absorption properties. |
---|
349 | ! call suaer_corrk ! Set up aerosol optical properties. |
---|
350 | |
---|
351 | |
---|
352 | ! now that L_NGAUSS has been initialized (by sugas_corrk) |
---|
353 | ! allocate related arrays |
---|
354 | if(.not.allocated(dtaui)) then |
---|
355 | ALLOCATE(dtaui(L_NLAYRAD,L_NSPECTI,L_NGAUSS), stat=ok) |
---|
356 | if (ok/=0) then |
---|
357 | write(*,*) "memory allocation failed for dtaui!" |
---|
358 | call abort_physic(subname,'allocation failure for dtaui',1) |
---|
359 | endif |
---|
360 | endif |
---|
361 | if(.not.allocated(dtauv)) then |
---|
362 | ALLOCATE(dtauv(L_NLAYRAD,L_NSPECTV,L_NGAUSS), stat=ok) |
---|
363 | if (ok/=0) then |
---|
364 | write(*,*) "memory allocation failed for dtauv!" |
---|
365 | call abort_physic(subname,'allocation failure for dtauv',1) |
---|
366 | endif |
---|
367 | endif |
---|
368 | if(.not.allocated(cosbv)) then |
---|
369 | ALLOCATE(cosbv(L_NLAYRAD,L_NSPECTV,L_NGAUSS), stat=ok) |
---|
370 | if (ok/=0) then |
---|
371 | write(*,*) "memory allocation failed for cosbv!" |
---|
372 | call abort_physic(subname,'allocation failure for cobsv',1) |
---|
373 | endif |
---|
374 | endif |
---|
375 | if(.not.allocated(cosbi)) then |
---|
376 | ALLOCATE(cosbi(L_NLAYRAD,L_NSPECTI,L_NGAUSS), stat=ok) |
---|
377 | if (ok/=0) then |
---|
378 | write(*,*) "memory allocation failed for cosbi!" |
---|
379 | call abort_physic(subname,'allocation failure for cobsi',1) |
---|
380 | endif |
---|
381 | endif |
---|
382 | if(.not.allocated(wbari)) then |
---|
383 | ALLOCATE(wbari(L_NLAYRAD,L_NSPECTI,L_NGAUSS), stat=ok) |
---|
384 | if (ok/=0) then |
---|
385 | write(*,*) "memory allocation failed for wbari!" |
---|
386 | call abort_physic(subname,'allocation failure for wbari',1) |
---|
387 | endif |
---|
388 | endif |
---|
389 | if(.not.allocated(wbarv)) then |
---|
390 | ALLOCATE(wbarv(L_NLAYRAD,L_NSPECTV,L_NGAUSS), stat=ok) |
---|
391 | if (ok/=0) then |
---|
392 | write(*,*) "memory allocation failed for wbarv!" |
---|
393 | call abort_physic(subname,'allocation failure for wbarv',1) |
---|
394 | endif |
---|
395 | endif |
---|
396 | if(.not.allocated(tauv)) then |
---|
397 | ALLOCATE(tauv(L_NLEVRAD,L_NSPECTV,L_NGAUSS), stat=ok) |
---|
398 | if (ok/=0) then |
---|
399 | write(*,*) "memory allocation failed for tauv!" |
---|
400 | call abort_physic(subname,'allocation failure for tauv',1) |
---|
401 | endif |
---|
402 | endif |
---|
403 | if(.not.allocated(taucumv)) then |
---|
404 | ALLOCATE(taucumv(L_LEVELS,L_NSPECTV,L_NGAUSS), stat=ok) |
---|
405 | if (ok/=0) then |
---|
406 | write(*,*) "memory allocation failed for taucumv!" |
---|
407 | call abort_physic(subname,'allocation failure for taucumv',1) |
---|
408 | endif |
---|
409 | endif |
---|
410 | if(.not.allocated(taucumi)) then |
---|
411 | ALLOCATE(taucumi(L_LEVELS,L_NSPECTI,L_NGAUSS), stat=ok) |
---|
412 | if (ok/=0) then |
---|
413 | write(*,*) "memory allocation failed for taucumi!" |
---|
414 | call abort_physic(subname,'allocation failure for taucumi',1) |
---|
415 | endif |
---|
416 | endif |
---|
417 | if(.not.allocated(taugsurf)) then |
---|
418 | ALLOCATE(taugsurf(L_NSPECTV,L_NGAUSS-1), stat=ok) |
---|
419 | if (ok/=0) then |
---|
420 | write(*,*) "memory allocation failed for taugsurf!" |
---|
421 | call abort_physic(subname,'allocation failure for taugsurf',1) |
---|
422 | endif |
---|
423 | endif |
---|
424 | if(.not.allocated(taugsurfi)) then |
---|
425 | ALLOCATE(taugsurfi(L_NSPECTI,L_NGAUSS-1), stat=ok) |
---|
426 | if (ok/=0) then |
---|
427 | write(*,*) "memory allocation failed for taugsurfi!" |
---|
428 | call abort_physic(subname,'allocation failure for taugsurfi',1) |
---|
429 | endif |
---|
430 | endif |
---|
431 | |
---|
432 | end if ! of if (firstcall) |
---|
433 | |
---|
434 | !======================================================================= |
---|
435 | ! I.b Initialization on every call |
---|
436 | !======================================================================= |
---|
437 | |
---|
438 | qxvaer(:,:,:)=0.0 |
---|
439 | qsvaer(:,:,:)=0.0 |
---|
440 | gvaer(:,:,:) =0.0 |
---|
441 | |
---|
442 | qxiaer(:,:,:)=0.0 |
---|
443 | qsiaer(:,:,:)=0.0 |
---|
444 | giaer(:,:,:) =0.0 |
---|
445 | |
---|
446 | OLR_nu(:,:) = 0. |
---|
447 | OSR_nu(:,:) = 0. |
---|
448 | GSR_nu(:,:) = 0. |
---|
449 | |
---|
450 | !-------------------------------------------------- |
---|
451 | ! Effective radius and variance of the aerosols |
---|
452 | ! Madeleine's PhD (eq. 2.3-2.4): |
---|
453 | ! --> r_eff = <r^3> / <r^2> |
---|
454 | ! --> nu_eff = <r^4>*<r^2> / <r^3>^2 - 1 |
---|
455 | !-------------------------------------------------- |
---|
456 | ! Radiative Hazes |
---|
457 | if (optichaze) then |
---|
458 | if (callmufi) then |
---|
459 | ! Spherical aerosols |
---|
460 | sig = 0.2 |
---|
461 | where (mp2m_rc_sph(:,:) >= 5e-9) |
---|
462 | reffrad(:,:,1) = mp2m_rc_sph(:,:) * exp(5.*sig**2 / 2.) |
---|
463 | elsewhere |
---|
464 | reffrad(:,:,1) = 0d0 |
---|
465 | endwhere |
---|
466 | nueffrad(:,:,1) = exp(sig**2) - 1 |
---|
467 | !if (exp(sig**2) - 1 > 0.1) then |
---|
468 | ! nueffrad(:,:,1) = exp(sig**2) - 1 |
---|
469 | !else |
---|
470 | ! nueffrad(:,:,1) = 0.1 |
---|
471 | !endif |
---|
472 | ! Fractal aerosols |
---|
473 | sig = 0.35 |
---|
474 | where (mp2m_rc_fra(:,:) >= 1e-8) |
---|
475 | reffrad(:,:,2) = mp2m_rc_fra(:,:) * exp(5.*sig**2 / 2.) |
---|
476 | elsewhere |
---|
477 | reffrad(:,:,2) = 0d0 |
---|
478 | endwhere |
---|
479 | nueffrad(:,:,2) = exp(sig**2) - 1 |
---|
480 | !if (exp(sig**2) - 1 > 0.1) then |
---|
481 | ! nueffrad(:,:,2) = exp(sig**2) - 1 |
---|
482 | !else |
---|
483 | ! nueffrad(:,:,2) = 0.1 |
---|
484 | !endif |
---|
485 | |
---|
486 | else |
---|
487 | do iaer=1,naerkind |
---|
488 | if ((iaer.eq.iaero_haze)) then |
---|
489 | call su_aer_radii(ngrid,nlayer,reffrad(1,1,iaer),nueffrad(1,1,iaer)) |
---|
490 | endif |
---|
491 | end do ! iaer = 1, naerkind. |
---|
492 | if (haze_radproffix) then |
---|
493 | call haze_reffrad_fix(ngrid,nlayer,zzlay,reffrad,nueffrad) |
---|
494 | if (is_master) print*, 'haze_radproffix=T : fixed profile for haze rad' |
---|
495 | else |
---|
496 | if (is_master) print*,'reffrad haze:',reffrad(1,1,iaero_haze) |
---|
497 | if (is_master) print*,'nueff haze',nueffrad(1,1,iaero_haze) |
---|
498 | endif ! end haze_radproffix |
---|
499 | endif ! end callmufi |
---|
500 | endif ! end radiative haze |
---|
501 | |
---|
502 | ! How much light do we get ? |
---|
503 | do nw=1,L_NSPECTV |
---|
504 | stel(nw)=stellarf(nw)/(dist_star**2) |
---|
505 | end do |
---|
506 | |
---|
507 | if (optichaze) then |
---|
508 | |
---|
509 | ! Get 3D aerosol optical properties. |
---|
510 | call aeroptproperties(ngrid,nlayer,reffrad,nueffrad, & |
---|
511 | QVISsQREF3d,omegaVIS3d,gVIS3d, & |
---|
512 | QIRsQREF3d,omegaIR3d,gIR3d, & |
---|
513 | QREFvis3d,QREFir3d) |
---|
514 | |
---|
515 | ! Get aerosol optical depths dtau_aer and diagnostic tau_col. |
---|
516 | call aeropacity(ngrid,nlayer,nq,pplay,pplev,zzlev,pt,pq,dtau_aer, & |
---|
517 | reffrad,nueffrad,QREFvis3d,QREFir3d, & |
---|
518 | tau_col) |
---|
519 | endif |
---|
520 | |
---|
521 | |
---|
522 | !----------------------------------------------------------------------- |
---|
523 | ! Prepare CH4 mixing ratio for radiative transfer |
---|
524 | IF (methane) then |
---|
525 | vmrch4(:,:)=0. |
---|
526 | |
---|
527 | if (ch4fix) then |
---|
528 | if (vmrch4_proffix) then |
---|
529 | !! Interpolate on the model vertical grid |
---|
530 | do ig=1,ngrid |
---|
531 | CALL interp_line(levdat,vmrdat,Nfine, & |
---|
532 | zzlay(ig,:)/1000.,vmrch4(ig,:),nlayer) |
---|
533 | enddo |
---|
534 | else |
---|
535 | vmrch4(:,:)=vmrch4fix |
---|
536 | endif |
---|
537 | else |
---|
538 | vmrch4(:,:)=pq(:,:,igcm_ch4_gas)*100.* & |
---|
539 | mmol(igcm_n2)/mmol(igcm_ch4_gas) |
---|
540 | endif |
---|
541 | ENDIF |
---|
542 | |
---|
543 | ! Prepare NON LTE correction in Pluto atmosphere |
---|
544 | IF (nlte) then |
---|
545 | CALL nlte_ch4(ngrid,nlayer,nq,pplay,pplev,pt,vmrch4,& |
---|
546 | eps_nlte_sw23,eps_nlte_sw33,eps_nlte_lw) |
---|
547 | ENDIF |
---|
548 | ! Net atmospheric radiative cooling rate from C2H2 (K.s-1): |
---|
549 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
550 | dtlw_hcn_c2h2=0. |
---|
551 | if (cooling) then |
---|
552 | CALL cooling_hcn_c2h2(ngrid,nlayer,pplay,& |
---|
553 | pt,dtlw_hcn_c2h2) |
---|
554 | endif |
---|
555 | |
---|
556 | |
---|
557 | !----------------------------------------------------------------------- |
---|
558 | do ig=1,ngrid ! Starting Big Loop over every GCM column |
---|
559 | !----------------------------------------------------------------------- |
---|
560 | |
---|
561 | |
---|
562 | !======================================================================= |
---|
563 | ! II. Transformation of the GCM variables |
---|
564 | !======================================================================= |
---|
565 | |
---|
566 | |
---|
567 | !----------------------------------------------------------------------- |
---|
568 | ! Aerosol optical properties Qext, Qscat and g. |
---|
569 | ! The transformation in the vertical is the same as for temperature. |
---|
570 | !----------------------------------------------------------------------- |
---|
571 | |
---|
572 | |
---|
573 | ! AF24: for now only consider one aerosol (=haze) |
---|
574 | if (optichaze) then |
---|
575 | do iaer=1,naerkind |
---|
576 | ! Shortwave. |
---|
577 | do nw=1,L_NSPECTV |
---|
578 | |
---|
579 | do l=1,nlayer |
---|
580 | |
---|
581 | temp1=QVISsQREF3d(ig,nlayer+1-l,nw,iaer) & |
---|
582 | *QREFvis3d(ig,nlayer+1-l,iaer) |
---|
583 | |
---|
584 | temp2=QVISsQREF3d(ig,max(nlayer-l,1),nw,iaer) & |
---|
585 | *QREFvis3d(ig,max(nlayer-l,1),iaer) |
---|
586 | |
---|
587 | qxvaer(2*l,nw,iaer) = temp1 |
---|
588 | qxvaer(2*l+1,nw,iaer)=(temp1+temp2)/2 |
---|
589 | |
---|
590 | temp1=temp1*omegavis3d(ig,nlayer+1-l,nw,iaer) |
---|
591 | temp2=temp2*omegavis3d(ig,max(nlayer-l,1),nw,iaer) |
---|
592 | |
---|
593 | qsvaer(2*l,nw,iaer) = temp1 |
---|
594 | qsvaer(2*l+1,nw,iaer)=(temp1+temp2)/2 |
---|
595 | |
---|
596 | temp1=gvis3d(ig,nlayer+1-l,nw,iaer) |
---|
597 | temp2=gvis3d(ig,max(nlayer-l,1),nw,iaer) |
---|
598 | |
---|
599 | gvaer(2*l,nw,iaer) = temp1 |
---|
600 | gvaer(2*l+1,nw,iaer)=(temp1+temp2)/2 |
---|
601 | |
---|
602 | end do ! nlayer |
---|
603 | |
---|
604 | qxvaer(1,nw,iaer)=qxvaer(2,nw,iaer) |
---|
605 | qxvaer(2*nlayer+1,nw,iaer)=0. |
---|
606 | |
---|
607 | qsvaer(1,nw,iaer)=qsvaer(2,nw,iaer) |
---|
608 | qsvaer(2*nlayer+1,nw,iaer)=0. |
---|
609 | |
---|
610 | gvaer(1,nw,iaer)=gvaer(2,nw,iaer) |
---|
611 | gvaer(2*nlayer+1,nw,iaer)=0. |
---|
612 | |
---|
613 | end do ! L_NSPECTV |
---|
614 | |
---|
615 | do nw=1,L_NSPECTI |
---|
616 | ! Longwave |
---|
617 | do l=1,nlayer |
---|
618 | |
---|
619 | temp1=QIRsQREF3d(ig,nlayer+1-l,nw,iaer) & |
---|
620 | *QREFir3d(ig,nlayer+1-l,iaer) |
---|
621 | |
---|
622 | temp2=QIRsQREF3d(ig,max(nlayer-l,1),nw,iaer) & |
---|
623 | *QREFir3d(ig,max(nlayer-l,1),iaer) |
---|
624 | |
---|
625 | qxiaer(2*l,nw,iaer) = temp1 |
---|
626 | qxiaer(2*l+1,nw,iaer)=(temp1+temp2)/2 |
---|
627 | |
---|
628 | temp1=temp1*omegair3d(ig,nlayer+1-l,nw,iaer) |
---|
629 | temp2=temp2*omegair3d(ig,max(nlayer-l,1),nw,iaer) |
---|
630 | |
---|
631 | qsiaer(2*l,nw,iaer) = temp1 |
---|
632 | qsiaer(2*l+1,nw,iaer)=(temp1+temp2)/2 |
---|
633 | |
---|
634 | temp1=gir3d(ig,nlayer+1-l,nw,iaer) |
---|
635 | temp2=gir3d(ig,max(nlayer-l,1),nw,iaer) |
---|
636 | |
---|
637 | giaer(2*l,nw,iaer) = temp1 |
---|
638 | giaer(2*l+1,nw,iaer)=(temp1+temp2)/2 |
---|
639 | |
---|
640 | end do ! nlayer |
---|
641 | |
---|
642 | qxiaer(1,nw,iaer)=qxiaer(2,nw,iaer) |
---|
643 | qxiaer(2*nlayer+1,nw,iaer)=0. |
---|
644 | |
---|
645 | qsiaer(1,nw,iaer)=qsiaer(2,nw,iaer) |
---|
646 | qsiaer(2*nlayer+1,nw,iaer)=0. |
---|
647 | |
---|
648 | giaer(1,nw,iaer)=giaer(2,nw,iaer) |
---|
649 | giaer(2*nlayer+1,nw,iaer)=0. |
---|
650 | |
---|
651 | end do ! L_NSPECTI |
---|
652 | |
---|
653 | end do ! naerkind |
---|
654 | |
---|
655 | ! Test / Correct for freaky s. s. albedo values. |
---|
656 | do iaer=1,naerkind |
---|
657 | do k=1,L_LEVELS |
---|
658 | |
---|
659 | do nw=1,L_NSPECTV |
---|
660 | if(qsvaer(k,nw,iaer).gt.1.05*qxvaer(k,nw,iaer))then |
---|
661 | message='Serious problems with qsvaer values' |
---|
662 | call abort_physic(subname,message,1) |
---|
663 | endif |
---|
664 | if(qsvaer(k,nw,iaer).gt.qxvaer(k,nw,iaer))then |
---|
665 | qsvaer(k,nw,iaer)=qxvaer(k,nw,iaer) |
---|
666 | endif |
---|
667 | end do |
---|
668 | |
---|
669 | do nw=1,L_NSPECTI |
---|
670 | if(qsiaer(k,nw,iaer).gt.1.05*qxiaer(k,nw,iaer))then |
---|
671 | message='Serious problems with qsvaer values' |
---|
672 | call abort_physic(subname,message,1) |
---|
673 | endif |
---|
674 | if(qsiaer(k,nw,iaer).gt.qxiaer(k,nw,iaer))then |
---|
675 | qsiaer(k,nw,iaer)=qxiaer(k,nw,iaer) |
---|
676 | endif |
---|
677 | end do |
---|
678 | |
---|
679 | end do ! L_LEVELS |
---|
680 | end do ! naerkind |
---|
681 | end if ! optichaze |
---|
682 | |
---|
683 | !----------------------------------------------------------------------- |
---|
684 | ! Aerosol optical depths |
---|
685 | !----------------------------------------------------------------------- |
---|
686 | if (optichaze) then |
---|
687 | do iaer=1,naerkind ! a bug was here |
---|
688 | do k=0,nlayer-1 |
---|
689 | |
---|
690 | pweight=(pplay(ig,L_NLAYRAD-k)-pplev(ig,L_NLAYRAD-k+1))/ & |
---|
691 | (pplev(ig,L_NLAYRAD-k)-pplev(ig,L_NLAYRAD-k+1)) |
---|
692 | ! As 'dtau_aer' is at reference (visible) wavelenght we scale it as |
---|
693 | ! it will be multplied by qxi/v in optci/v |
---|
694 | temp=dtau_aer(ig,L_NLAYRAD-k,iaer)/QREFvis3d(ig,L_NLAYRAD-k,iaer) |
---|
695 | tauaero(2*k+2,iaer)=max(temp*pweight,0.d0) |
---|
696 | tauaero(2*k+3,iaer)=max(temp-tauaero(2*k+2,iaer),0.d0) |
---|
697 | |
---|
698 | end do |
---|
699 | ! boundary conditions |
---|
700 | tauaero(1,iaer) = tauaero(2,iaer) |
---|
701 | !tauaero(1,iaer) = 0. |
---|
702 | !JL18 at time of testing, the two above conditions gave the same results bit for bit. |
---|
703 | |
---|
704 | end do ! naerkind |
---|
705 | else |
---|
706 | tauaero(:,:)=0 |
---|
707 | end if ! optichaze |
---|
708 | |
---|
709 | ! Albedo and Emissivity. |
---|
710 | albi=1-emis(ig) ! Long Wave. |
---|
711 | DO nw=1,L_NSPECTV ! Short Wave loop. |
---|
712 | albv(nw)=albedo(ig,nw) |
---|
713 | ENDDO |
---|
714 | |
---|
715 | acosz=mu0(ig) ! Cosine of sun incident angle : 3D simulations or local 1D simulations using latitude. |
---|
716 | |
---|
717 | !----------------------------------------------------------------------- |
---|
718 | ! Methane Vapor |
---|
719 | !----------------------------------------------------------------------- |
---|
720 | if (methane) then |
---|
721 | do l=1,nlayer |
---|
722 | qvar(2*l) = vmrch4(ig,nlayer+1-l)/100.* & |
---|
723 | mmol(igcm_ch4_gas)/mmol(igcm_n2) |
---|
724 | qvar(2*l+1) = ((vmrch4(ig,nlayer+1-l)+vmrch4(ig, & |
---|
725 | max(nlayer-l,1)))/2.)/100.* & |
---|
726 | mmol(igcm_ch4_gas)/mmol(igcm_n2) |
---|
727 | end do |
---|
728 | qvar(1)=qvar(2) |
---|
729 | |
---|
730 | else |
---|
731 | do k=1,L_LEVELS |
---|
732 | qvar(k) = 1.0D-7 |
---|
733 | end do |
---|
734 | end if ! end methane |
---|
735 | |
---|
736 | !----------------------------------------------------------------------- |
---|
737 | ! kcm mode only ! |
---|
738 | !----------------------------------------------------------------------- |
---|
739 | |
---|
740 | DO l=1,nlayer |
---|
741 | muvarrad(2*l) = muvar(ig,nlayer+2-l) |
---|
742 | muvarrad(2*l+1) = (muvar(ig,nlayer+2-l)+muvar(ig,max(nlayer+1-l,1)))/2 |
---|
743 | END DO |
---|
744 | |
---|
745 | muvarrad(1) = muvarrad(2) |
---|
746 | muvarrad(2*nlayer+1)=muvar(ig,1) |
---|
747 | |
---|
748 | ! Keep values inside limits for which we have radiative transfer coefficients !!! |
---|
749 | if(L_REFVAR.gt.1)then ! (there was a bug here) |
---|
750 | do k=1,L_LEVELS |
---|
751 | if(qvar(k).lt.wrefvar(1))then |
---|
752 | qvar(k)=wrefvar(1)+1.0e-8 |
---|
753 | elseif(qvar(k).gt.wrefvar(L_REFVAR))then |
---|
754 | qvar(k)=wrefvar(L_REFVAR)-1.0e-8 |
---|
755 | endif |
---|
756 | end do |
---|
757 | endif |
---|
758 | |
---|
759 | !----------------------------------------------------------------------- |
---|
760 | ! Pressure and temperature |
---|
761 | !----------------------------------------------------------------------- |
---|
762 | |
---|
763 | DO l=1,nlayer |
---|
764 | plevrad(2*l) = pplay(ig,nlayer+1-l)/scalep |
---|
765 | plevrad(2*l+1) = pplev(ig,nlayer+1-l)/scalep |
---|
766 | tlevrad(2*l) = pt(ig,nlayer+1-l) |
---|
767 | tlevrad(2*l+1) = (pt(ig,nlayer+1-l)+pt(ig,max(nlayer-l,1)))/2 |
---|
768 | END DO |
---|
769 | |
---|
770 | plevrad(1) = 0. |
---|
771 | !!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. |
---|
772 | |
---|
773 | tlevrad(1) = tlevrad(2) |
---|
774 | tlevrad(2*nlayer+1)=tsurf(ig) |
---|
775 | |
---|
776 | pmid(1) = pplay(ig,nlayer)/scalep |
---|
777 | pmid(2) = pmid(1) |
---|
778 | |
---|
779 | tmid(1) = tlevrad(2) |
---|
780 | tmid(2) = tmid(1) |
---|
781 | |
---|
782 | DO l=1,L_NLAYRAD-1 |
---|
783 | tmid(2*l+1) = tlevrad(2*l+1) |
---|
784 | tmid(2*l+2) = tlevrad(2*l+1) |
---|
785 | pmid(2*l+1) = plevrad(2*l+1) |
---|
786 | pmid(2*l+2) = plevrad(2*l+1) |
---|
787 | END DO |
---|
788 | pmid(L_LEVELS) = plevrad(L_LEVELS) |
---|
789 | tmid(L_LEVELS) = tlevrad(L_LEVELS) |
---|
790 | |
---|
791 | !!Alternative interpolation: |
---|
792 | ! pmid(3) = pmid(1) |
---|
793 | ! pmid(4) = pmid(1) |
---|
794 | ! tmid(3) = tmid(1) |
---|
795 | ! tmid(4) = tmid(1) |
---|
796 | ! DO l=2,L_NLAYRAD-1 |
---|
797 | ! tmid(2*l+1) = tlevrad(2*l) |
---|
798 | ! tmid(2*l+2) = tlevrad(2*l) |
---|
799 | ! pmid(2*l+1) = plevrad(2*l) |
---|
800 | ! pmid(2*l+2) = plevrad(2*l) |
---|
801 | ! END DO |
---|
802 | ! pmid(L_LEVELS) = plevrad(L_LEVELS-1) |
---|
803 | ! tmid(L_LEVELS) = tlevrad(L_LEVELS-1) |
---|
804 | |
---|
805 | ! Test for out-of-bounds pressure. |
---|
806 | if(plevrad(3).lt.pgasmin)then |
---|
807 | print*,'Warning: minimum pressure is outside the radiative' |
---|
808 | print*,'transfer kmatrix bounds, exiting.' |
---|
809 | print*,'Pressure:', plevrad(3), 'Pa' |
---|
810 | message="Minimum pressure outside of kmatrix bounds" |
---|
811 | !call abort_physic(subname,message,1) |
---|
812 | elseif(plevrad(L_LEVELS).gt.pgasmax)then |
---|
813 | print*,'Maximum pressure is outside the radiative' |
---|
814 | print*,'transfer kmatrix bounds, exiting.' |
---|
815 | message="Minimum pressure outside of kmatrix bounds" |
---|
816 | call abort_physic(subname,message,1) |
---|
817 | endif |
---|
818 | |
---|
819 | ! Test for out-of-bounds temperature. |
---|
820 | ! -- JVO 20 : Also add a sanity test checking that tlevrad is |
---|
821 | ! within Planck function temperature boundaries, |
---|
822 | ! which would cause gfluxi/sfluxi to crash. |
---|
823 | do k=1,L_LEVELS |
---|
824 | |
---|
825 | if(tlevrad(k).lt.tgasmin)then |
---|
826 | print*,'Minimum temperature is outside the radiative' |
---|
827 | print*,'transfer kmatrix bounds' |
---|
828 | print*,"k=",k," tlevrad(k)=",tlevrad(k) |
---|
829 | print*,"tgasmin=",tgasmin |
---|
830 | if (strictboundcorrk) then |
---|
831 | message="Minimum temperature outside of kmatrix bounds" |
---|
832 | call abort_physic(subname,message,1) |
---|
833 | else |
---|
834 | print*,'***********************************************' |
---|
835 | print*,'we allow model to continue with tlevrad<tgasmin' |
---|
836 | print*,' ... we assume we know what you are doing ... ' |
---|
837 | print*,' ... but do not let this happen too often ... ' |
---|
838 | print*,'***********************************************' |
---|
839 | !tlevrad(k)=tgasmin ! Used in the source function ! |
---|
840 | endif |
---|
841 | elseif(tlevrad(k).gt.tgasmax)then |
---|
842 | print*,'Maximum temperature is outside the radiative' |
---|
843 | print*,'transfer kmatrix bounds, exiting.' |
---|
844 | print*,"k=",k," tlevrad(k)=",tlevrad(k) |
---|
845 | print*,"tgasmax=",tgasmax |
---|
846 | if (strictboundcorrk) then |
---|
847 | message="Maximum temperature outside of kmatrix bounds" |
---|
848 | call abort_physic(subname,message,1) |
---|
849 | else |
---|
850 | print*,'***********************************************' |
---|
851 | print*,'we allow model to continue with tlevrad>tgasmax' |
---|
852 | print*,' ... we assume we know what you are doing ... ' |
---|
853 | print*,' ... but do not let this happen too often ... ' |
---|
854 | print*,'***********************************************' |
---|
855 | !tlevrad(k)=tgasmax ! Used in the source function ! |
---|
856 | endif |
---|
857 | endif |
---|
858 | |
---|
859 | if (tlevrad(k).lt.tplanckmin) then |
---|
860 | print*,'Minimum temperature is outside the boundaries for' |
---|
861 | print*,'Planck function integration set in callphys.def, aborting.' |
---|
862 | print*,"k=",k," tlevrad(k)=",tlevrad(k) |
---|
863 | print*,"tplanckmin=",tplanckmin |
---|
864 | message="Minimum temperature outside Planck function bounds - Change tplanckmin in callphys.def" |
---|
865 | call abort_physic(subname,message,1) |
---|
866 | else if (tlevrad(k).gt.tplanckmax) then |
---|
867 | print*,'Maximum temperature is outside the boundaries for' |
---|
868 | print*,'Planck function integration set in callphys.def, aborting.' |
---|
869 | print*,"k=",k," tlevrad(k)=",tlevrad(k) |
---|
870 | print*,"tplanckmax=",tplanckmax |
---|
871 | message="Maximum temperature outside Planck function bounds - Change tplanckmax in callphys.def" |
---|
872 | call abort_physic(subname,message,1) |
---|
873 | endif |
---|
874 | |
---|
875 | enddo |
---|
876 | |
---|
877 | do k=1,L_NLAYRAD+1 |
---|
878 | if(tmid(k).lt.tgasmin)then |
---|
879 | print*,'Minimum temperature is outside the radiative' |
---|
880 | print*,'transfer kmatrix bounds, exiting.' |
---|
881 | print*,"k=",k," tmid(k)=",tmid(k) |
---|
882 | print*,"tgasmin=",tgasmin |
---|
883 | if (strictboundcorrk) then |
---|
884 | message="Minimum temperature outside of kmatrix bounds" |
---|
885 | call abort_physic(subname,message,1) |
---|
886 | else |
---|
887 | print*,'***********************************************' |
---|
888 | print*,'we allow model to continue but with tmid=tgasmin' |
---|
889 | print*,' ... we assume we know what you are doing ... ' |
---|
890 | print*,' ... but do not let this happen too often ... ' |
---|
891 | print*,'***********************************************' |
---|
892 | tmid(k)=tgasmin |
---|
893 | endif |
---|
894 | elseif(tmid(k).gt.tgasmax)then |
---|
895 | print*,'Maximum temperature is outside the radiative' |
---|
896 | print*,'transfer kmatrix bounds, exiting.' |
---|
897 | print*,"k=",k," tmid(k)=",tmid(k) |
---|
898 | print*,"tgasmax=",tgasmax |
---|
899 | if (strictboundcorrk) then |
---|
900 | message="Maximum temperature outside of kmatrix bounds" |
---|
901 | call abort_physic(subname,message,1) |
---|
902 | else |
---|
903 | print*,'***********************************************' |
---|
904 | print*,'we allow model to continue but with tmid=tgasmax' |
---|
905 | print*,' ... we assume we know what you are doing ... ' |
---|
906 | print*,' ... but do not let this happen too often ... ' |
---|
907 | print*,'***********************************************' |
---|
908 | tmid(k)=tgasmax |
---|
909 | endif |
---|
910 | endif |
---|
911 | enddo |
---|
912 | |
---|
913 | !======================================================================= |
---|
914 | ! III. Calling the main radiative transfer subroutines |
---|
915 | !======================================================================= |
---|
916 | |
---|
917 | ! ---------------------------------------------------------------- |
---|
918 | ! Recombine reference corrk tables if needed - Added by JVO, 2020. |
---|
919 | if (corrk_recombin) then |
---|
920 | call call_recombin(ig,nlayer,pq(ig,:,:),pplay(ig,:),pt(ig,:),qvar(:),tmid(:),pmid(:)) |
---|
921 | endif |
---|
922 | ! ---------------------------------------------------------------- |
---|
923 | |
---|
924 | Cmk= 0.01 * 1.0 / (glat(ig) * mugaz * 1.672621e-27) ! q_main=1.0 assumed. |
---|
925 | glat_ig=glat(ig) |
---|
926 | |
---|
927 | !----------------------------------------------------------------------- |
---|
928 | ! Short Wave Part |
---|
929 | !----------------------------------------------------------------------- |
---|
930 | |
---|
931 | ! Call everywhere for diagnostics. |
---|
932 | call optcv(dtauv,tauv,taucumv,plevrad, & |
---|
933 | qxvaer,qsvaer,gvaer,wbarv,cosbv,tauray,tauaero, & |
---|
934 | tmid,pmid,taugsurf,qvar,muvarrad) |
---|
935 | |
---|
936 | if(fract(ig) .ge. 1.0e-4) then ! Only during daylight. |
---|
937 | if((ngrid.eq.1).and.(global1d))then |
---|
938 | do nw=1,L_NSPECTV |
---|
939 | stel_fract(nw)= stel(nw)* 0.25 / acosz ! globally averaged = divide by 4, and we correct for solar zenith angle |
---|
940 | end do |
---|
941 | else |
---|
942 | do nw=1,L_NSPECTV |
---|
943 | stel_fract(nw)= stel(nw) * fract(ig) |
---|
944 | end do |
---|
945 | endif |
---|
946 | |
---|
947 | call sfluxv(dtauv,tauv,taucumv,albv,dwnv,wbarv,cosbv, & |
---|
948 | acosz,stel_fract,nfluxtopv,fluxtopvdn,nfluxoutv_nu,& |
---|
949 | nfluxgndv_nu,nfluxtopv_nu, & |
---|
950 | fmnetv,fmnetv_nu,fluxupv,fluxdnv,fzerov,taugsurf) |
---|
951 | |
---|
952 | else ! During the night, fluxes = 0. |
---|
953 | nfluxtopv = 0.0d0 |
---|
954 | fluxtopvdn = 0.0d0 |
---|
955 | nfluxtopv_nu(:) = 0.0d0 |
---|
956 | nfluxoutv_nu(:) = 0.0d0 |
---|
957 | nfluxgndv_nu(:) = 0.0d0 |
---|
958 | fmnetv_nu(:,:) = 0.0d0 |
---|
959 | fmnetv(:) = 0.0d0 |
---|
960 | fluxupv(:) = 0.0d0 |
---|
961 | fluxdnv(:) = 0.0d0 |
---|
962 | end if |
---|
963 | |
---|
964 | |
---|
965 | ! Equivalent Albedo Calculation (for OUTPUT). MT2015 |
---|
966 | if(fract(ig) .ge. 1.0e-4) then ! equivalent albedo makes sense only during daylight. |
---|
967 | surface_stellar_flux=sum(nfluxgndv_nu(1:L_NSPECTV)) |
---|
968 | if(surface_stellar_flux .gt. 1.0e-3) then ! equivalent albedo makes sense only if the stellar flux received by the surface is positive. |
---|
969 | DO nw=1,L_NSPECTV |
---|
970 | albedo_temp(nw)=albedo(ig,nw)*nfluxgndv_nu(nw) |
---|
971 | ENDDO |
---|
972 | albedo_temp(1:L_NSPECTV)=albedo_temp(1:L_NSPECTV)/surface_stellar_flux |
---|
973 | albedo_equivalent(ig)=sum(albedo_temp(1:L_NSPECTV)) |
---|
974 | else |
---|
975 | albedo_equivalent(ig)=0.0 ! Spectrally Integrated Albedo not defined for non-irradiated grid points. So we arbitrary set the equivalent albedo to 0. |
---|
976 | endif |
---|
977 | else |
---|
978 | albedo_equivalent(ig)=0.0 ! Spectrally Integrated Albedo not defined for non-irradiated grid points. So we arbitrary set the equivalent albedo to 0. |
---|
979 | endif |
---|
980 | |
---|
981 | |
---|
982 | !----------------------------------------------------------------------- |
---|
983 | ! Long Wave Part |
---|
984 | !----------------------------------------------------------------------- |
---|
985 | |
---|
986 | call optci(plevrad,tlevrad,dtaui,taucumi, & |
---|
987 | qxiaer,qsiaer,giaer,cosbi,wbari,tauaero,tmid,pmid, & |
---|
988 | taugsurfi,qvar,muvarrad) |
---|
989 | |
---|
990 | call sfluxi(plevrad,tlevrad,dtaui,taucumi,ubari,albi, & |
---|
991 | wnoi,dwni,cosbi,wbari,nfluxtopi,nfluxtopi_nu, & |
---|
992 | fmneti,fmneti_nu,fluxupi,fluxdni,fluxupi_nu,fzeroi,taugsurfi) |
---|
993 | |
---|
994 | !----------------------------------------------------------------------- |
---|
995 | ! Transformation of the correlated-k code outputs |
---|
996 | ! (into dtlw, dtsw, fluxsurf_lw, fluxsurf_sw, fluxtop_lw, fluxtop_sw) |
---|
997 | |
---|
998 | ! Flux incident at the top of the atmosphere |
---|
999 | fluxtop_dn(ig)=fluxtopvdn |
---|
1000 | |
---|
1001 | fluxtop_lw(ig) = real(nfluxtopi) |
---|
1002 | fluxabs_sw(ig) = real(-nfluxtopv) |
---|
1003 | fluxsurf_lw(ig) = real(fluxdni(L_NLAYRAD)) |
---|
1004 | fluxsurf_sw(ig) = real(fluxdnv(L_NLAYRAD)) |
---|
1005 | |
---|
1006 | ! Flux absorbed by the surface. By MT2015. |
---|
1007 | fluxsurfabs_sw(ig) = fluxsurf_sw(ig)*(1.-albedo_equivalent(ig)) |
---|
1008 | |
---|
1009 | if(fluxtop_dn(ig).lt.0.0)then |
---|
1010 | print*,'Achtung! fluxtop_dn has lost the plot!' |
---|
1011 | print*,'fluxtop_dn=',fluxtop_dn(ig) |
---|
1012 | print*,'acosz=',acosz |
---|
1013 | print*,'dtau_aer=',dtau_aer(ig,:,:) |
---|
1014 | print*,'temp= ',pt(ig,:) |
---|
1015 | print*,'pplay= ',pplay(ig,:) |
---|
1016 | message="Achtung! fluxtop_dn has lost the plot!" |
---|
1017 | call abort_physic(subname,message,1) |
---|
1018 | endif |
---|
1019 | |
---|
1020 | ! Spectral output, for exoplanet observational comparison |
---|
1021 | if(specOLR)then |
---|
1022 | do nw=1,L_NSPECTI |
---|
1023 | OLR_nu(ig,nw)=nfluxtopi_nu(nw)/DWNI(nw) !JL Normalize to the bandwidth |
---|
1024 | end do |
---|
1025 | do nw=1,L_NSPECTV |
---|
1026 | GSR_nu(ig,nw)=nfluxgndv_nu(nw)/DWNV(nw) |
---|
1027 | OSR_nu(ig,nw)=nfluxoutv_nu(nw)/DWNV(nw) !JL Normalize to the bandwidth |
---|
1028 | end do |
---|
1029 | endif |
---|
1030 | |
---|
1031 | ! Finally, the heating rates |
---|
1032 | DO l=2,L_NLAYRAD |
---|
1033 | ! dtsw(ig,L_NLAYRAD+1-l)=(fmnetv(l)-fmnetv(l-1)) & |
---|
1034 | ! *glat(ig)/(cpp*scalep*(plevrad(2*l+1)-plevrad(2*l-1))) |
---|
1035 | dpp = glat(ig)/(cpp*scalep*(plevrad(2*l+1)-plevrad(2*l-1))) |
---|
1036 | do nw=1,L_NSPECTV |
---|
1037 | dtsw_nu(L_NLAYRAD+1-l,nw)= & |
---|
1038 | (fmnetv_nu(l,nw)-fmnetv_nu(l-1,nw))*dpp |
---|
1039 | end do |
---|
1040 | |
---|
1041 | ! dtlw(ig,L_NLAYRAD+1-l)=(fmneti(l)-fmneti(l-1)) & |
---|
1042 | ! *glat(ig)/(cpp*scalep*(plevrad(2*l+1)-plevrad(2*l-1))) |
---|
1043 | do nw=1,L_NSPECTI |
---|
1044 | dtlw_nu(L_NLAYRAD+1-l,nw)= & |
---|
1045 | (fmneti_nu(l,nw)-fmneti_nu(l-1,nw))*dpp |
---|
1046 | end do |
---|
1047 | END DO |
---|
1048 | |
---|
1049 | ! These are values at top of atmosphere |
---|
1050 | ! dtsw(ig,L_NLAYRAD)=(fmnetv(1)-nfluxtopv) & |
---|
1051 | ! *glat(ig)/(cpp*scalep*(plevrad(3)-plevrad(2))) |
---|
1052 | ! dtlw(ig,L_NLAYRAD)=(fmneti(1)-nfluxtopi) & |
---|
1053 | ! *glat(ig)/(cpp*scalep*(plevrad(3)-plevrad(2))) |
---|
1054 | dpp = g/(cpp*scalep*(plevrad(3)-plevrad(1))) |
---|
1055 | do nw=1,L_NSPECTV |
---|
1056 | dtsw_nu(L_NLAYRAD,nw)= & |
---|
1057 | (fmnetv_nu(1,nw)-nfluxtopv_nu(nw))*dpp |
---|
1058 | end do |
---|
1059 | do nw=1,L_NSPECTI |
---|
1060 | dtlw_nu(L_NLAYRAD,nw)= & |
---|
1061 | (fmneti_nu(1,nw)-nfluxtopi_nu(nw))*dpp |
---|
1062 | end do |
---|
1063 | |
---|
1064 | ! Optical thickness diagnostics |
---|
1065 | ! Output exp(-tau) because gweight ponderates exp and not tau itself |
---|
1066 | int_dtauv(ig,:,:) = 0.0d0 |
---|
1067 | int_dtaui(ig,:,:) = 0.0d0 |
---|
1068 | do l=1,L_NLAYRAD |
---|
1069 | do nw=1,L_NSPECTV |
---|
1070 | do k=1,L_NGAUSS |
---|
1071 | int_dtauv(ig,l,nw) = int_dtauv(ig,l,nw) + exp(-dtauv(l,nw,k))*gweight(k) |
---|
1072 | enddo |
---|
1073 | enddo |
---|
1074 | do nw=1,L_NSPECTI |
---|
1075 | do k=1,L_NGAUSS |
---|
1076 | int_dtaui(ig,l,nw) = int_dtaui(ig,l,nw) + exp(-dtaui(l,nw,k))*gweight(k) |
---|
1077 | enddo |
---|
1078 | enddo |
---|
1079 | enddo |
---|
1080 | |
---|
1081 | ! ********************************************************** |
---|
1082 | ! NON NLTE correction in Pluto atmosphere |
---|
1083 | ! And conversion of LW spectral heating rates to total rates |
---|
1084 | ! ********************************************************** |
---|
1085 | |
---|
1086 | if (.not.nlte) then |
---|
1087 | eps_nlte_sw23(ig,:) =1. ! IF no NLTE |
---|
1088 | eps_nlte_sw33(ig,:) =1. ! IF no NLTE |
---|
1089 | eps_nlte_lw(ig,:) =1. ! IF no NLTE |
---|
1090 | endif |
---|
1091 | |
---|
1092 | do l=1,nlayer |
---|
1093 | |
---|
1094 | !LW |
---|
1095 | dtlw(ig,l) =0. |
---|
1096 | ! dtlw_co(ig,l) =0. ! only for diagnostic |
---|
1097 | do nw=1,L_NSPECTI |
---|
1098 | ! wewei : wavelength in micrometer |
---|
1099 | if ((wavei(nw).gt.6.).and.(wavei(nw).lt.9)) then |
---|
1100 | dtlw_nu(l,nw)=dtlw_nu(l,nw)*eps_nlte_lw(ig,l) |
---|
1101 | else |
---|
1102 | !dtlw_nu(l,nw)=1.*dtlw_nu(l,nw) ! no CO correction (Strobbel 1996) |
---|
1103 | dtlw_nu(l,nw)=0.33*dtlw_nu(l,nw) ! CO correction (Strobbel 1996) |
---|
1104 | ! dtlw_co(ig,l)=dtlw_co(ig,l)+ dtlw_nu(l,nw) ! diagnostic |
---|
1105 | end if |
---|
1106 | dtlw(ig,l)=dtlw(ig,l)+ dtlw_nu(l,nw) !average now on each wavelength |
---|
1107 | end do |
---|
1108 | ! adding c2h2 if cooling active |
---|
1109 | if (cooling) then |
---|
1110 | dtlw(ig,l)=dtlw(ig,l)+dtlw_hcn_c2h2(ig,l) |
---|
1111 | endif |
---|
1112 | |
---|
1113 | !SW |
---|
1114 | dtsw(ig,l) =0. |
---|
1115 | |
---|
1116 | if (strobel) then |
---|
1117 | |
---|
1118 | do nw=1,L_NSPECTV |
---|
1119 | if ((wavev(nw).gt.2).and.(wavev(nw).lt.2.6)) then |
---|
1120 | dtsw_nu(l,nw)=dtsw_nu(l,nw)*eps_nlte_sw23(ig,l) |
---|
1121 | elseif ((wavev(nw).gt.3).and.(wavev(nw).lt.3.6)) then |
---|
1122 | dtsw_nu(l,nw)=dtsw_nu(l,nw)*eps_nlte_sw33(ig,l) |
---|
1123 | else |
---|
1124 | dtsw_nu(l,nw)=dtsw_nu(l,nw) |
---|
1125 | end if |
---|
1126 | dtsw(ig,l)=dtsw(ig,l)+ dtsw_nu(l,nw) |
---|
1127 | end do |
---|
1128 | |
---|
1129 | else ! total heating rates multiplied by nlte coef |
---|
1130 | |
---|
1131 | do nw=1,L_NSPECTV |
---|
1132 | dtsw_nu(l,nw)=dtsw_nu(l,nw)*eps_nlte_sw23(ig,l) |
---|
1133 | dtsw(ig,l)=dtsw(ig,l)+ dtsw_nu(l,nw) |
---|
1134 | enddo |
---|
1135 | |
---|
1136 | endif |
---|
1137 | |
---|
1138 | |
---|
1139 | end do |
---|
1140 | ! ********************************************************** |
---|
1141 | |
---|
1142 | |
---|
1143 | !----------------------------------------------------------------------- |
---|
1144 | end do ! End of big loop over every GCM column. |
---|
1145 | !----------------------------------------------------------------------- |
---|
1146 | |
---|
1147 | |
---|
1148 | !----------------------------------------------------------------------- |
---|
1149 | ! Additional diagnostics |
---|
1150 | !----------------------------------------------------------------------- |
---|
1151 | |
---|
1152 | ! IR spectral output, for exoplanet observational comparison |
---|
1153 | if(lastcall.and.(ngrid.eq.1))then ! could disable the 1D output, they are in the diagfi and diagspec... JL12 |
---|
1154 | |
---|
1155 | print*,'Saving scalar quantities in surf_vals.out...' |
---|
1156 | print*,'psurf = ', pplev(1,1),' Pa' |
---|
1157 | open(116,file='surf_vals.out') |
---|
1158 | write(116,*) tsurf(1),pplev(1,1),fluxtop_dn(1), & |
---|
1159 | real(-nfluxtopv),real(nfluxtopi) |
---|
1160 | close(116) |
---|
1161 | |
---|
1162 | |
---|
1163 | ! USEFUL COMMENT - Do Not Remove. |
---|
1164 | ! |
---|
1165 | ! if(specOLR)then |
---|
1166 | ! open(117,file='OLRnu.out') |
---|
1167 | ! do nw=1,L_NSPECTI |
---|
1168 | ! write(117,*) OLR_nu(1,nw) |
---|
1169 | ! enddo |
---|
1170 | ! close(117) |
---|
1171 | ! |
---|
1172 | ! open(127,file='OSRnu.out') |
---|
1173 | ! do nw=1,L_NSPECTV |
---|
1174 | ! write(127,*) OSR_nu(1,nw) |
---|
1175 | ! enddo |
---|
1176 | ! close(127) |
---|
1177 | ! endif |
---|
1178 | |
---|
1179 | ! OLR vs altitude: do it as a .txt file. |
---|
1180 | OLRz=.false. |
---|
1181 | if(OLRz)then |
---|
1182 | print*,'saving IR vertical flux for OLRz...' |
---|
1183 | open(118,file='OLRz_plevs.out') |
---|
1184 | open(119,file='OLRz.out') |
---|
1185 | do l=1,L_NLAYRAD |
---|
1186 | write(118,*) plevrad(2*l) |
---|
1187 | do nw=1,L_NSPECTI |
---|
1188 | write(119,*) fluxupi_nu(l,nw) |
---|
1189 | enddo |
---|
1190 | enddo |
---|
1191 | close(118) |
---|
1192 | close(119) |
---|
1193 | endif |
---|
1194 | |
---|
1195 | endif |
---|
1196 | |
---|
1197 | ! See physiq.F for explanations about CLFvarying. This is temporary. |
---|
1198 | if (lastcall) then |
---|
1199 | IF( ALLOCATED( gasi ) ) DEALLOCATE( gasi ) |
---|
1200 | IF( ALLOCATED( gasv ) ) DEALLOCATE( gasv ) |
---|
1201 | !$OMP BARRIER |
---|
1202 | !$OMP MASTER |
---|
1203 | IF( ALLOCATED( pgasref ) ) DEALLOCATE( pgasref ) |
---|
1204 | IF( ALLOCATED( tgasref ) ) DEALLOCATE( tgasref ) |
---|
1205 | IF( ALLOCATED( wrefvar ) ) DEALLOCATE( wrefvar ) |
---|
1206 | IF( ALLOCATED( pfgasref ) ) DEALLOCATE( pfgasref ) |
---|
1207 | IF( ALLOCATED( gweight ) ) DEALLOCATE( gweight ) |
---|
1208 | !$OMP END MASTER |
---|
1209 | !$OMP BARRIER |
---|
1210 | IF ( ALLOCATED(reffrad)) DEALLOCATE(reffrad) |
---|
1211 | IF ( ALLOCATED(nueffrad)) DEALLOCATE(nueffrad) |
---|
1212 | endif |
---|
1213 | |
---|
1214 | |
---|
1215 | end subroutine callcorrk |
---|
1216 | |
---|
1217 | END MODULE callcorrk_mod |
---|