source: trunk/LMDZ.PLUTO/libf/phypluto/callcorrk_pluto_mod.F90 @ 3585

Last change on this file since 3585 was 3585, checked in by debatzbr, 9 hours ago

Connecting microphysics to radiative transfer + miscellaneous cleans

File size: 33.0 KB
Line 
1MODULE callcorrk_pluto_mod
2
3IMPLICIT NONE
4
5CONTAINS
6
7    subroutine callcorrk_pluto(icount,ngrid,nlayer,pq,nq,qsurf,   &
8          albedo,emis,mu0,pplev,pplay,pt,   &
9          zzlay,zzlev,tsurf,fract,dist_star,dtau_aer,    &
10          dtlw,dtsw,fluxsurf_lw,    &
11          fluxsurf_sw,fluxtop_lw,fluxtop_sw,fluxtop_dn, &
12          reffrad,tau_col,ptime,pday,firstcall,lastcall)
13
14      use radinc_h
15      use radcommon_h
16      use ioipsl_getincom
17      use suaer_corrk_mod, only: suaer_corrk
18      use radii_mod, only: su_aer_radii,haze_reffrad_fix
19      use aeropacity_mod, only: aeropacity
20      use aeroptproperties_mod, only: aeroptproperties
21      use aerosol_mod, only : iaero_haze
22      use datafile_mod, only: datadir
23      use comcstfi_mod, only: pi,g,cpp,mugaz
24      use tracer_h, only: igcm_n2,igcm_ch4_gas,igcm_ch4_ice,rho_ch4_ice,lw_ch4,&
25                          mmol
26      use callkeys_mod, only: optichaze,ch4fix,cooling,methane,nlte,&
27                              strobel,vmrch4_proffix,specOLR,vmrch4fix,&
28                              haze_radproffix,callmufi
29      use optcv_pluto_mod, only: optcv_pluto
30      use optci_pluto_mod, only: optci_pluto
31      use sfluxi_pluto_mod, only: sfluxi_pluto
32      use sfluxv_pluto_mod, only: sfluxv_pluto
33      use mod_phys_lmdz_para, only : is_master
34      use mp2m_diagnostics
35
36      implicit none
37
38!==================================================================
39!
40!     Purpose
41!     -------
42!     Solve the radiative transfer using the correlated-k method for
43!     the gaseous absorption and the Toon et al. (1989) method for
44!     scatttering due to aerosols.
45!
46!     Authors
47!     -------
48!     Emmanuel 01/2001, Forget 09/2001
49!     Robin Wordsworth (2009)
50!     Modif Pluton Tanguy Bertrand 2017
51!==================================================================
52
53include "dimensions.h"
54
55!-----------------------------------------------------------------------
56!     Declaration of the arguments (INPUT - OUTPUT) on the LMD GCM grid
57!     Layer #1 is the layer near the ground.
58!     Layer #nlayer is the layer at the top.
59
60!     INPUT
61      INTEGER icount
62      INTEGER ngrid,nlayer
63      REAL dtau_aer(ngrid,nlayer,naerkind) ! aerosol opacity tau
64      REAL albedo(ngrid)                   ! SW albedo
65      REAL emis(ngrid)                     ! LW emissivity
66      REAL pplay(ngrid,nlayer)             ! pres. level in GCM mid of layer
67      REAL zzlay(ngrid,nlayer)             ! Altitude at the middle of the layers
68      REAL zzlev(ngrid,nlayer)             ! Altitude at the layer boundaries.
69      REAL pplev(ngrid,nlayer+1)           ! pres. level at GCM layer boundaries
70 
71      REAL pt(ngrid,nlayer)               ! air temperature (K)
72      REAL tsurf(ngrid)                     ! surface temperature (K)
73      REAL dist_star,mu0(ngrid)             ! distance star-planet (AU)
74      REAL fract(ngrid)                     ! fraction of day
75
76!     Globally varying aerosol optical properties on GCM grid
77!     Not needed everywhere so not in radcommon_h
78      REAL :: QVISsQREF3d(ngrid,nlayer,L_NSPECTV,naerkind)
79      REAL :: omegaVIS3d(ngrid,nlayer,L_NSPECTV,naerkind)
80      REAL :: gVIS3d(ngrid,nlayer,L_NSPECTV,naerkind)
81
82      REAL :: QIRsQREF3d(ngrid,nlayer,L_NSPECTI,naerkind)
83      REAL :: omegaIR3d(ngrid,nlayer,L_NSPECTI,naerkind)
84      REAL :: gIR3d(ngrid,nlayer,L_NSPECTI,naerkind)
85
86      REAL :: QREFvis3d(ngrid,nlayer,naerkind)
87      REAL :: QREFir3d(ngrid,nlayer,naerkind)
88
89!      REAL :: omegaREFvis3d(ngrid,nlayer,naerkind)
90!      REAL :: omegaREFir3d(ngrid,nlayer,naerkind) ! not sure of the point of these...
91
92      REAL reffrad(ngrid,nlayer,naerkind)
93      REAL nueffrad(ngrid,nlayer,naerkind)
94      REAL profhaze(ngrid,nlayer) ! TB17 fixed profile of haze mmr
95
96!     OUTPUT
97      REAL dtsw(ngrid,nlayer) ! heating rate (K/s) due to SW
98      REAL dtlw(ngrid,nlayer) ! heating rate (K/s) due to LW
99      REAL fluxsurf_lw(ngrid)   ! incident LW flux to surf (W/m2)
100      REAL fluxtop_lw(ngrid)    ! outgoing LW flux to space (W/m2)
101      REAL fluxsurf_sw(ngrid)   ! incident SW flux to surf (W/m2)
102      REAL fluxtop_sw(ngrid)    ! outgoing LW flux to space (W/m2)
103      REAL fluxtop_dn(ngrid)    ! incident top of atmosphere SW flux (W/m2)
104
105!-----------------------------------------------------------------------
106!     Declaration of the variables required by correlated-k subroutines
107!     Numbered from top to bottom unlike in the GCM!
108
109      REAL*8 tmid(L_LEVELS),pmid(L_LEVELS)
110      REAL*8 tlevrad(L_LEVELS),plevrad(L_LEVELS)
111
112!     Optical values for the optci/cv subroutines
113      REAL*8 stel(L_NSPECTV),stel_fract(L_NSPECTV)
114      REAL*8 dtaui(L_NLAYRAD,L_NSPECTI,L_NGAUSS)
115      REAL*8 dtauv(L_NLAYRAD,L_NSPECTV,L_NGAUSS)
116      REAL*8 cosbv(L_NLAYRAD,L_NSPECTV,L_NGAUSS)
117      REAL*8 cosbi(L_NLAYRAD,L_NSPECTI,L_NGAUSS)
118      REAL*8 wbari(L_NLAYRAD,L_NSPECTI,L_NGAUSS)
119      REAL*8 wbarv(L_NLAYRAD,L_NSPECTV,L_NGAUSS)
120      REAL*8 tauv(L_NLEVRAD,L_NSPECTV,L_NGAUSS)
121      REAL*8 taucumv(L_LEVELS,L_NSPECTV,L_NGAUSS)
122      REAL*8 taucumi(L_LEVELS,L_NSPECTI,L_NGAUSS)
123
124      REAL*8 tauaero(L_LEVELS+1,naerkind)
125      REAL*8 nfluxtopv,nfluxtopi,nfluxtop
126      real*8 NFLUXTOPV_nu(L_NSPECTV)
127      real*8 NFLUXTOPI_nu(L_NSPECTI)
128      real*8 fluxupi_nu(L_NLAYRAD,L_NSPECTI) ! for 1D diagnostic
129      REAL*8 fmneti(L_NLAYRAD),fmnetv(L_NLAYRAD)
130      real*8 fmneti_nu(L_NLAYRAD,L_NSPECTI) !
131      real*8 fmnetv_nu(L_NLAYRAD,L_NSPECTV) !
132      REAL*8 fluxupv(L_NLAYRAD),fluxupi(L_NLAYRAD)
133      REAL*8 fluxdnv(L_NLAYRAD),fluxdni(L_NLAYRAD)
134      REAL*8 albi,albv,acosz
135
136      INTEGER ig,l,k,nw,iaer,irad
137
138      real*8 sig
139
140      real fluxtoplanet
141      real*8 taugsurf(L_NSPECTV,L_NGAUSS-1)
142      real*8 taugsurfi(L_NSPECTI,L_NGAUSS-1)
143
144      real*8 qvar(L_LEVELS)          ! mixing ratio of variable component
145      REAL pq(ngrid,nlayer,nq)
146      REAL qsurf(ngrid,nq)       ! tracer on surface (e.g. kg.m-2)
147      integer nq
148
149!     Local aerosol optical properties for each column on RADIATIVE grid
150      real*8  QXVAER(L_LEVELS+1,L_NSPECTV,naerkind)
151      real*8  QSVAER(L_LEVELS+1,L_NSPECTV,naerkind)
152      real*8  GVAER(L_LEVELS+1,L_NSPECTV,naerkind)
153      real*8  QXIAER(L_LEVELS+1,L_NSPECTI,naerkind)
154      real*8  QSIAER(L_LEVELS+1,L_NSPECTI,naerkind)
155      real*8  GIAER(L_LEVELS+1,L_NSPECTI,naerkind)
156
157    !   save qxvaer, qsvaer, gvaer
158    !   save qxiaer, qsiaer, giaer
159    !   save QREFvis3d, QREFir3d
160
161      REAL tau_col(ngrid) ! diagnostic from aeropacity
162
163!     Misc.
164      logical firstcall, lastcall, nantest
165      real*8  tempv(L_NSPECTV)
166      real*8  tempi(L_NSPECTI)
167      real*8  temp,temp1,temp2,pweight
168      character(len=10) :: tmp1
169      character(len=10) :: tmp2
170
171!     for fixed vapour profiles
172      real RH
173      real*8 pq_temp(nlayer)
174      real ptemp, Ttemp, qsat
175
176!     for OLR spec
177      integer OLRcount
178      save OLRcount
179      integer OLRcount2
180      save OLRcount2
181      character(len=2)  :: tempOLR
182      character(len=30) :: filenomOLR
183      real ptime, pday
184
185      REAL epsi_ch4
186      SAVE epsi_ch4
187
188      logical diagrad_OLRz,diagrad_OLR,diagrad_surf,diagrad_rates
189      real OLR_nu(ngrid,L_NSPECTI)
190
191!     NLTE factor for CH4
192      real eps_nlte_sw23(ngrid,nlayer) ! CH4 NLTE efficiency factor for zdtsw
193      real eps_nlte_sw33(ngrid,nlayer) ! CH4 NLTE efficiency factor for zdtsw
194      real eps_nlte_lw(ngrid,nlayer) ! CH4 NLTE efficiency factor for zdtsw
195      REAL dtlw_nu(nlayer,L_NSPECTI) ! heating rate (K/s) due to LW in spectral bands
196      REAL dtsw_nu(nlayer,L_NSPECTV) ! heating rate (K/s) due to SW in spectral bands
197
198      REAL dpp  ! intermediate
199
200      REAL dtlw_co(ngrid, nlayer) ! cooling rate (K/s) due to CO (diagnostic)
201      REAL dtlw_hcn_c2h2(ngrid, nlayer) ! cooling rate (K/s) due to C2H2/HCN (diagnostic)
202
203      !!read altitudes and vmrch4
204      integer Nfine,ifine
205      parameter(Nfine=701)
206      character(len=100) :: file_path
207      real,save :: levdat(Nfine),vmrdat(Nfine)
208      real :: vmrch4(ngrid,nlayer)              ! vmr ch4 from vmrch4_proffix
209
210!=======================================================================
211!     Initialization on first call
212
213      qxvaer(:,:,:) = 0
214      qsvaer(:,:,:) = 0
215      gvaer(:,:,:) = 0
216
217      qxiaer(:,:,:) = 0
218      qsiaer(:,:,:) = 0
219      giaer(:,:,:) = 0
220
221
222      if(firstcall) then
223
224         if (is_master) print*, "callcorrk: Correlated-k data folder:",trim(datadir)
225         call getin("corrkdir",corrkdir)
226         print*, "corrkdir = ",corrkdir
227         write( tmp1, '(i3)' ) L_NSPECTI
228         write( tmp2, '(i3)' ) L_NSPECTV
229         banddir=trim(adjustl(tmp1))//'x'//trim(adjustl(tmp2))
230         banddir=trim(adjustl(corrkdir))//'/'//trim(adjustl(banddir))
231
232         if (is_master) print*,'starting sugas'
233         call sugas_corrk       ! set up gaseous absorption properties
234         if (is_master) print*,'starting setspi'
235         call setspi            ! basic infrared properties
236         if (is_master) print*,'starting setspv'
237         call setspv            ! basic visible properties
238
239         !--------------------------------------------------
240         !     Effective radius and variance of the aerosols
241         !--------------------------------------------------
242         ! Radiative Hazes
243         if (optichaze) then
244           if (callmufi) then
245            ! Spherical aerosols
246            sig = 0.2
247            WHERE(mp2m_rc_sph(:,:) > 1e-10)
248               reffrad(:,:,1) = mp2m_rc_sph(:,:) * exp(5.*sig**2 / 2.)
249            ELSEWHERE
250               reffrad(:,:,1) = 0d0         
251            ENDWHERE
252            nueffrad(:,:,1) = exp(sig**2) - 1
253            ! Fractal aerosols
254            sig = 0.35
255            WHERE(mp2m_rc_fra(:,:) > 1e-10)
256               reffrad(:,:,2) = mp2m_rc_fra(:,:) * exp(5.*sig**2 / 2.)
257            ELSEWHERE
258               reffrad(:,:,2) = 0d0
259            ENDWHERE
260            nueffrad(:,:,2) = exp(sig**2) - 1
261           
262           else
263            do iaer=1,naerkind
264              if ((iaer.eq.iaero_haze)) then
265              call su_aer_radii(ngrid,nlayer,reffrad(1,1,iaer),nueffrad(1,1,iaer))
266              endif
267            end do ! iaer = 1, naerkind.
268            if (haze_radproffix) then
269              call haze_reffrad_fix(ngrid,nlayer,zzlay,reffrad,nueffrad)
270              if (is_master) print*, 'haze_radproffix=T : fixed profile for haze rad'
271            else
272              if (is_master) print*,'reffrad haze:',reffrad(1,1,iaero_haze)
273              if (is_master) print*,'nueff haze',nueffrad(1,1,iaero_haze)
274            endif ! end haze_radproffix
275           endif ! end callmufi
276         endif ! end radiative haze
277
278         Cmk= 0.01 * 1.0 / (g * mugaz * 1.672621e-27) ! q_main=1.0 assumed
279
280         if (methane) then
281           epsi_ch4=mmol(igcm_ch4_gas)/mmol(igcm_n2)
282         endif
283
284         ! If fixed profile of CH4 gas
285         IF (vmrch4_proffix) then
286            file_path=trim(datadir)//'/gas_prop/vmr_ch4.txt'
287            open(115,file=file_path,form='formatted')
288            do ifine=1,Nfine
289              read(115,*) levdat(ifine), vmrdat(ifine)
290            enddo
291            close(115)
292         ENDIF
293
294      end if ! firstcall
295
296!=======================================================================
297!     L_NSPECTV is the number of Visual(or Solar) spectral intervals
298!     how much light we get
299      fluxtoplanet=0
300      DO nw=1,L_NSPECTV
301         stel(nw)=stellarf(nw)/(dist_star**2)  !flux
302         fluxtoplanet=fluxtoplanet + stel(nw)
303      END DO
304
305!-----------------------------------------------------------------------
306!     Get 3D aerosol optical properties.
307      ! ici on selectionne les proprietes opt correspondant a reffrad
308      if (optichaze) then
309        !--------------------------------------------------
310        ! Effective radius and variance of the aerosols if profil non
311        ! uniform. Called only if optichaze=true.
312        !--------------------------------------------------
313
314        call aeroptproperties(ngrid,nlayer,reffrad,nueffrad,         &
315            QVISsQREF3d,omegaVIS3d,gVIS3d,                          &
316            QIRsQREF3d,omegaIR3d,gIR3d,                             &
317            QREFvis3d,QREFir3d)
318
319        ! Get aerosol optical depths.
320        call aeropacity(ngrid,nlayer,nq,pplay,pplev,zzlev,pt,pq,dtau_aer,      &
321            reffrad,nueffrad,QREFvis3d,QREFir3d,                             &
322            tau_col)
323      endif
324
325!-----------------------------------------------------------------------
326!     Prepare CH4 mixing ratio for radiative transfer
327      IF (methane) then
328        vmrch4(:,:)=0.
329
330        if (ch4fix) then
331           if (vmrch4_proffix) then
332            !! Interpolate on the model vertical grid
333             do ig=1,ngrid
334               CALL interp_line(levdat,vmrdat,Nfine,    &
335                       zzlay(ig,:)/1000.,vmrch4(ig,:),nlayer)
336             enddo
337           else
338            vmrch4(:,:)=vmrch4fix
339           endif
340        else
341            vmrch4(:,:)=pq(:,:,igcm_ch4_gas)*100.*  &
342                         mmol(igcm_n2)/mmol(igcm_ch4_gas)
343        endif
344      ENDIF
345
346!     Prepare NON LTE correction in Pluto atmosphere
347      IF (nlte) then
348        CALL nlte_ch4(ngrid,nlayer,nq,pplay,pplev,pt,vmrch4,    &
349                  eps_nlte_sw23,eps_nlte_sw33,eps_nlte_lw)
350      ENDIF
351!     Net atmospheric radiative cooling rate from C2H2 (K.s-1):
352!     ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
353      dtlw_hcn_c2h2=0.
354      if (cooling) then
355         CALL cooling_hcn_c2h2(ngrid,nlayer,pplay,  &
356                                  pt,dtlw_hcn_c2h2)
357      endif
358
359!-----------------------------------------------------------------------
360!=======================================================================
361!-----------------------------------------------------------------------
362!     starting big loop over every GCM column
363      do ig=1,ngrid
364
365!=======================================================================
366!     Transformation of the GCM variables
367
368!-----------------------------------------------------------------------
369!     Aerosol optical properties Qext, Qscat and g on each band
370!     The transformation in the vertical is the same as for temperature
371         if (optichaze) then
372!     shortwave
373            do iaer=1,naerkind
374               DO nw=1,L_NSPECTV
375                  do l=1,nlayer
376
377                     temp1=QVISsQREF3d(ig,nlayer+1-l,nw,iaer) &
378                         *QREFvis3d(ig,nlayer+1-l,iaer)
379
380                     temp2=QVISsQREF3d(ig,max(nlayer-l,1),nw,iaer)    &
381                         *QREFvis3d(ig,max(nlayer-l,1),iaer)
382                     qxvaer(2*l,nw,iaer)  = temp1
383                     qxvaer(2*l+1,nw,iaer)=(temp1+temp2)/2
384
385                    temp1=temp1*omegavis3d(ig,nlayer+1-l,nw,iaer)
386                    temp2=temp2*omegavis3d(ig,max(nlayer-l,1),nw,iaer)
387
388                     qsvaer(2*l,nw,iaer)  = temp1
389                     qsvaer(2*l+1,nw,iaer)=(temp1+temp2)/2
390
391                     temp1=gvis3d(ig,nlayer+1-l,nw,iaer)
392                     temp2=gvis3d(ig,max(nlayer-l,1),nw,iaer)
393
394                     gvaer(2*l,nw,iaer)  = temp1
395                     gvaer(2*l+1,nw,iaer)=(temp1+temp2)/2
396
397                  end do
398                  qxvaer(1,nw,iaer)=qxvaer(2,nw,iaer)
399                  qxvaer(2*nlayer+1,nw,iaer)=0.
400
401                  qsvaer(1,nw,iaer)=qsvaer(2,nw,iaer)
402                  qsvaer(2*nlayer+1,nw,iaer)=0.
403
404                  gvaer(1,nw,iaer)=gvaer(2,nw,iaer)
405                  gvaer(2*nlayer+1,nw,iaer)=0.
406               end do
407
408!     longwave
409               DO nw=1,L_NSPECTI
410                  do l=1,nlayer
411
412                     temp1=QIRsQREF3d(ig,nlayer+1-l,nw,iaer)  &
413                         *QREFir3d(ig,nlayer+1-l,iaer)
414
415                     temp2=QIRsQREF3d(ig,max(nlayer-l,1),nw,iaer) &
416                         *QREFir3d(ig,max(nlayer-l,1),iaer)
417
418                     qxiaer(2*l,nw,iaer)  = temp1
419                     qxiaer(2*l+1,nw,iaer)=(temp1+temp2)/2
420
421                     temp1=temp1*omegair3d(ig,nlayer+1-l,nw,iaer)
422                     temp2=temp2*omegair3d(ig,max(nlayer-l,1),nw,iaer)
423
424                     qsiaer(2*l,nw,iaer)  = temp1
425                     qsiaer(2*l+1,nw,iaer)=(temp1+temp2)/2
426
427                     temp1=gir3d(ig,nlayer+1-l,nw,iaer)
428                     temp2=gir3d(ig,max(nlayer-l,1),nw,iaer)
429
430                     giaer(2*l,nw,iaer)  = temp1
431                     giaer(2*l+1,nw,iaer)=(temp1+temp2)/2
432
433                  end do
434
435                  qxiaer(1,nw,iaer)=qxiaer(2,nw,iaer)
436                  qxiaer(2*nlayer+1,nw,iaer)=0.
437
438                  qsiaer(1,nw,iaer)=qsiaer(2,nw,iaer)
439                  qsiaer(2*nlayer+1,nw,iaer)=0.
440
441                  giaer(1,nw,iaer)=giaer(2,nw,iaer)
442                  giaer(2*nlayer+1,nw,iaer)=0.
443               end do
444            end do   ! naerkind
445
446            ! Test / Correct for freaky s. s. albedo values.
447            do iaer=1,naerkind
448               do k=1,L_LEVELS
449
450                  do nw=1,L_NSPECTV
451                     if(qsvaer(k,nw,iaer).gt.1.05*qxvaer(k,nw,iaer))then
452                        print*,'Serious problems with qsvaer values'
453                        print*,'in callcorrk'
454                        call abort
455                     endif
456                     if(qsvaer(k,nw,iaer).gt.qxvaer(k,nw,iaer))then
457                        qsvaer(k,nw,iaer)=qxvaer(k,nw,iaer)
458                     endif
459                  end do
460
461                  do nw=1,L_NSPECTI
462                     if(qsiaer(k,nw,iaer).gt.1.05*qxiaer(k,nw,iaer))then
463                        print*,'Serious problems with qsiaer values'
464                        print*,'in callcorrk'
465                        call abort
466                     endif
467                     if(qsiaer(k,nw,iaer).gt.qxiaer(k,nw,iaer))then
468                        qsiaer(k,nw,iaer)=qxiaer(k,nw,iaer)
469                     endif
470                  end do
471               end do ! levels
472
473            end do ! naerkind
474
475         endif ! optichaze
476
477!-----------------------------------------------------------------------
478!     Aerosol optical depths
479         IF (optichaze) THEN
480          do iaer=1,naerkind     ! heritage generic
481            do k=0,nlayer-1
482               pweight= &
483                   (pplay(ig,L_NLAYRAD-k)-pplev(ig,L_NLAYRAD-k+1))/ &
484                   (pplev(ig,L_NLAYRAD-k)-pplev(ig,L_NLAYRAD-k+1))
485               if (QREFvis3d(ig,L_NLAYRAD-k,iaer).ne.0) then
486                 temp=dtau_aer(ig,L_NLAYRAD-k,iaer)/ &
487                   QREFvis3d(ig,L_NLAYRAD-k,iaer)
488               else
489                 print*, 'stop corrk',k,QREFvis3d(ig,L_NLAYRAD-k,iaer)
490                 stop
491               end if
492               tauaero(2*k+2,iaer)=max(temp*pweight,0.d0)
493               tauaero(2*k+3,iaer)=max(temp-tauaero(2*k+2,iaer),0.d0) ! tauaero en L_LEVELS soit deux fois plus que nlayer
494            end do
495
496            ! generic New boundary conditions
497            tauaero(1,iaer)          = tauaero(2,iaer)
498            !tauaero(L_LEVELS+1,iaer) = tauaero(L_LEVELS,iaer)
499            !tauaero(1,iaer)          = 0.
500            !tauaero(L_LEVELS+1,iaer) = 0.
501
502          end do   ! naerkind
503         ELSE
504           tauaero(:,:)=0
505         ENDIF
506!-----------------------------------------------------------------------
507
508!     Albedo and emissivity
509         albi=1-emis(ig)        ! longwave
510         albv=albedo(ig)        ! shortwave
511         acosz=mu0(ig)          ! cosine of sun incident angle
512
513!-----------------------------------------------------------------------
514!     Methane vapour
515
516!     qvar = mixing ratio
517!     L_LEVELS (51) different de GCM levels (25) . L_LEVELS = 2*(llm-1)+3=2*(ngrid-1)+3
518!     L_REFVAR  The number of different mixing ratio values in
519!     datagcm/composition.in for the k-coefficients.
520         qvar(:)=0.
521         IF (methane) then
522
523           do l=1,nlayer
524               qvar(2*l)   = vmrch4(ig,nlayer+1-l)/100.*    &
525                                    mmol(igcm_ch4_gas)/mmol(igcm_n2)
526               qvar(2*l+1) = ((vmrch4(ig,nlayer+1-l)+vmrch4(ig, &
527                            max(nlayer-l,1)))/2.)/100.* &
528                            mmol(igcm_ch4_gas)/mmol(igcm_n2)
529           end do
530           qvar(1)=qvar(2)
531
532
533      ! Keep values inside limits for which we have radiative transfer coefficients
534           if(L_REFVAR.gt.1)then ! there was a bug here!
535             do k=1,L_LEVELS
536               if(qvar(k).lt.wrefvar(1))then
537                 qvar(k)=wrefvar(1)+1.0e-8
538               elseif(qvar(k).gt.wrefvar(L_REFVAR))then
539                 qvar(k)=wrefvar(L_REFVAR)-1.0e-8
540               endif
541             end do
542           endif
543
544      ! IMPORTANT: Now convert from kg/kg to mol/mol
545           do k=1,L_LEVELS
546            qvar(k) = qvar(k)/(epsi_ch4+qvar(k)*(1.-epsi_ch4))
547           end do
548         ENDIF ! methane
549
550!-----------------------------------------------------------------------
551!     Pressure and temperature
552
553! generic updated:
554         DO l=1,nlayer
555           plevrad(2*l)   = pplay(ig,nlayer+1-l)/scalep
556           plevrad(2*l+1) = pplev(ig,nlayer+1-l)/scalep
557           tlevrad(2*l)   = pt(ig,nlayer+1-l)
558           tlevrad(2*l+1) = (pt(ig,nlayer+1-l)+pt(ig,max(nlayer-l,1)))/2
559         END DO
560
561         plevrad(1) = 0.
562         plevrad(2) = 0.   !! Trick to have correct calculations of fluxes
563                        ! in gflux(i/v).F, but the pmid levels are not impacted by this change.
564
565         tlevrad(1) = tlevrad(2)
566         tlevrad(2*nlayer+1)=tsurf(ig)
567
568         pmid(1) = max(pgasmin,0.0001*plevrad(3))
569         pmid(2) =  pmid(1)
570
571         tmid(1) = tlevrad(2)
572         tmid(2) = tmid(1)
573
574! INI
575!        DO l=1,nlayer
576!           plevrad(2*l)   = pplay(ig,nlayer+1-l)/scalep
577!           plevrad(2*l+1) = pplev(ig,nlayer+1-l)/scalep
578!           tlevrad(2*l)   = pt(ig,nlayer+1-l)
579!           tlevrad(2*l+1) = (pt(ig,nlayer+1-l)+pt(ig,
580!     $        max(nlayer-l,1)))/2
581!        END DO
582
583!! following lines changed in 03/2015 to solve upper atmosphere bug
584!        plevrad(1) = 0.
585!        plevrad(2) = max(pgasmin,0.0001*plevrad(3))
586!
587!        tlevrad(1) = tlevrad(2)
588!        tlevrad(2*nlayer+1)=tsurf(ig)
589!
590!        tmid(1) = tlevrad(2)
591!        tmid(2) = tlevrad(2)
592!
593!        pmid(1) = plevrad(2)
594!        pmid(2) = plevrad(2)
595
596         DO l=1,L_NLAYRAD-1
597           tmid(2*l+1) = tlevrad(2*l+1)
598           tmid(2*l+2) = tlevrad(2*l+1)
599           pmid(2*l+1) = plevrad(2*l+1)
600           pmid(2*l+2) = plevrad(2*l+1)
601         END DO
602! end of changes
603         pmid(L_LEVELS) = plevrad(L_LEVELS)
604         tmid(L_LEVELS) = tlevrad(L_LEVELS)
605
606      !TB
607         if ((PMID(2).le.1.e-5).and.(ig.eq.1)) then
608          if ((TMID(2).le.30.).and.(ig.eq.1)) then
609           write(*,*) 'Caution! Pres/temp of upper levels lower than'
610           write(*,*) 'ref pressure/temp: kcoef fixed for upper levels'
611!!     cf tpindex.F
612          endif
613         endif
614
615      !  test for out-of-bounds pressure
616         if(plevrad(3).lt.pgasmin)then
617           print*,'Minimum pressure is outside the radiative'
618           print*,'transfer kmatrix bounds, exiting.'
619      !     call abort
620         elseif(plevrad(L_LEVELS).gt.pgasmax)then
621           print*,'Maximum pressure is outside the radiative'
622           print*,'transfer kmatrix bounds, exiting.'
623      !     call abort
624         endif
625
626      !  test for out-of-bounds temperature
627         do k=1,L_LEVELS
628          if(tlevrad(k).lt.tgasmin)then
629            print*,'Minimum temperature is outside the radiative'
630            print*,'transfer kmatrix bounds, exiting.'
631            print*,'t(',k,')=',tlevrad(k),' < ',tgasmin
632            ! call abort
633          elseif(tlevrad(k).gt.tgasmax)then
634            print*,'Maximum temperature is outside the radiative'
635            print*,'transfer kmatrix bounds, exiting.'
636            print*,'t(',k,')=',tlevrad(k),' > ',tgasmax
637            ! call abort
638          endif
639         enddo
640
641!=======================================================================
642!     Calling the main radiative transfer subroutines
643
644!-----------------------------------------------------------------------
645!     Shortwave
646
647         IF(fract(ig) .GE. 1.0e-4) THEN ! only during daylight  IPM?! flux UV...
648
649            fluxtoplanet=0.
650            DO nw=1,L_NSPECTV
651               stel_fract(nw)= stel(nw) * fract(ig)
652               fluxtoplanet=fluxtoplanet + stel_fract(nw)
653            END DO
654
655            !print*, 'starting optcv'
656            call optcv_pluto(dtauv,tauv,taucumv,plevrad,  &
657                qxvaer,qsvaer,gvaer,wbarv,cosbv,tauray,tauaero, &
658                tmid,pmid,taugsurf,qvar)
659
660            call sfluxv_pluto(dtauv,tauv,taucumv,albv,dwnv,wbarv,cosbv,   &
661                acosz,stel_fract,nfluxtopv,nfluxtopv_nu,    &
662                fmnetv,fmnetv_nu,fluxupv,fluxdnv,fzerov,taugsurf)
663
664         ELSE                          ! during the night, fluxes = 0
665            nfluxtopv=0.0
666            fmnetv_nu(:,:)=0.0
667            nfluxtopv_nu(:)=0.0
668            DO l=1,L_NLAYRAD
669               fmnetv(l)=0.0
670               fluxupv(l)=0.0
671               fluxdnv(l)=0.0
672            END DO
673         END IF
674
675!-----------------------------------------------------------------------
676!     Longwave
677
678         call optci_pluto(plevrad,tlevrad,dtaui,taucumi,  &
679             qxiaer,qsiaer,giaer,cosbi,wbari,tauaero,tmid,pmid, &
680             taugsurfi,qvar)
681         call sfluxi_pluto(plevrad,tlevrad,dtaui,taucumi,ubari,albi,  &
682          wnoi,dwni,cosbi,wbari,nfluxtopi,nfluxtopi_nu, &
683          fmneti,fmneti_nu,fluxupi,fluxdni,fluxupi_nu,fzeroi,taugsurfi)
684
685
686!-----------------------------------------------------------------------
687!     Transformation of the correlated-k code outputs
688!     (into dtlw, dtsw, fluxsurf_lw, fluxsurf_sw, fluxtop_lw, fluxtop_sw)
689
690         fluxtop_lw(ig)  = fluxupi(1)
691         fluxsurf_lw(ig) = fluxdni(L_NLAYRAD)
692         fluxtop_sw(ig)  = fluxupv(1)
693         fluxsurf_sw(ig) = fluxdnv(L_NLAYRAD)
694
695!     Flux incident at the top of the atmosphere
696         fluxtop_dn(ig)=fluxdnv(1)
697
698!        IR spectral output from top of the atmosphere
699         if(specOLR)then
700            do nw=1,L_NSPECTI
701               OLR_nu(ig,nw)=nfluxtopi_nu(nw)
702            end do
703         endif
704
705! **********************************************************
706!     Finally, the heating rates
707!     g/cp*DF/DP
708! **********************************************************
709
710         DO l=2,L_NLAYRAD
711               dpp = g/(cpp*scalep*(plevrad(2*l+1)-plevrad(2*l-1)))
712
713               ! DTSW :
714               !dtsw(ig,L_NLAYRAD+1-l)=(fmnetv(l)-fmnetv(l-1))*dpp   !averaged dtlw on each wavelength
715               do nw=1,L_NSPECTV
716                 dtsw_nu(L_NLAYRAD+1-l,nw)= &
717                   (fmnetv_nu(l,nw)-fmnetv_nu(l-1,nw))*dpp
718               end do
719
720               ! DTLW : detailed with wavelength so that we can apply NLTE
721               !dtlw(ig,L_NLAYRAD+1-l)=(fmneti(l)-fmneti(l-1))*dpp   !averaged dtlw on each wavelength
722               do nw=1,L_NSPECTI
723                 dtlw_nu(L_NLAYRAD+1-l,nw)= &
724                   (fmneti_nu(l,nw)-fmneti_nu(l-1,nw))*dpp
725               end do
726         END DO
727         ! values at top of atmosphere
728         dpp = g/(cpp*scalep*(plevrad(3)-plevrad(1)))
729
730         ! SW
731         !dtsw(ig,L_NLAYRAD)=(fmnetv(1)-nfluxtopv)*dpp
732         do nw=1,L_NSPECTV
733            dtsw_nu(L_NLAYRAD,nw)=  &
734             (fmnetv_nu(1,nw)-nfluxtopv_nu(nw))*dpp
735         end do
736
737         ! LW
738!        dtlw(ig,L_NLAYRAD)=(fmneti(1)-nfluxtopi) *dpp
739         do nw=1,L_NSPECTI
740             dtlw_nu(L_NLAYRAD,nw)= &
741             (fmneti_nu(1,nw)-nfluxtopi_nu(nw))*dpp
742         end do
743
744! **********************************************************
745!     NON NLTE correction in Pluto atmosphere
746!     And conversion of LW spectral heating rates to total rates
747! **********************************************************
748
749         if (.not.nlte) then
750            eps_nlte_sw23(ig,:) =1. ! IF no NLTE
751            eps_nlte_sw33(ig,:) =1. ! IF no NLTE
752            eps_nlte_lw(ig,:) =1. ! IF no NLTE
753         endif
754
755         do l=1,nlayer
756
757            !LW
758            dtlw(ig,l) =0.
759!           dtlw_co(ig,l) =0.  ! only for diagnostic
760            do nw=1,L_NSPECTI
761              ! wewei : wavelength in micrometer
762              if ((wavei(nw).gt.6.).and.(wavei(nw).lt.9)) then
763                dtlw_nu(l,nw)=dtlw_nu(l,nw)*eps_nlte_lw(ig,l)
764              else
765                !dtlw_nu(l,nw)=1.*dtlw_nu(l,nw) ! no CO correction (Strobbel 1996)
766                dtlw_nu(l,nw)=0.33*dtlw_nu(l,nw) ! CO correction (Strobbel 1996)
767!               dtlw_co(ig,l)=dtlw_co(ig,l)+ dtlw_nu(l,nw) ! diagnostic
768              end if
769              dtlw(ig,l)=dtlw(ig,l)+ dtlw_nu(l,nw) !average now on each wavelength
770            end do
771            ! adding c2h2 if cooling active
772            dtlw(ig,l)=dtlw(ig,l)+dtlw_hcn_c2h2(ig,l)
773
774            !SW
775            dtsw(ig,l) =0.
776
777            if (strobel) then
778
779             do nw=1,L_NSPECTV
780              if ((wavev(nw).gt.2).and.(wavev(nw).lt.2.6)) then
781                dtsw_nu(l,nw)=dtsw_nu(l,nw)*eps_nlte_sw23(ig,l)
782              elseif ((wavev(nw).gt.3).and.(wavev(nw).lt.3.6)) then
783                dtsw_nu(l,nw)=dtsw_nu(l,nw)*eps_nlte_sw33(ig,l)
784              else
785                dtsw_nu(l,nw)=dtsw_nu(l,nw)
786              end if
787              dtsw(ig,l)=dtsw(ig,l)+ dtsw_nu(l,nw)
788             end do
789
790            else ! total heating rates multiplied by nlte coef
791
792             do nw=1,L_NSPECTV
793                dtsw_nu(l,nw)=dtsw_nu(l,nw)*eps_nlte_sw23(ig,l)
794                dtsw(ig,l)=dtsw(ig,l)+ dtsw_nu(l,nw)
795             enddo
796
797            endif
798
799
800         end do
801! **********************************************************
802
803!     Diagnotics for last call for each grid point
804         !if (lastcall) then
805
806          !print*,'albedi vis=',albv
807          !print*,'albedo ir=',albi
808          !print*,'fluxup ir (:)=',fluxupi
809          !print*,'flux ir net (:)=',fluxdni-fluxupi
810          !print*,'cumulative flux net ir (:)=',fmneti
811          !print*,'cumulative flux net vis (:)=',fmnetv
812          !print*,'fluxdn vis (:)=',fluxdnv
813          !print*,'fluxtop vis=',fluxtop_sw
814          !print*,'fluxsurf vis=',fluxsurf_sw
815          !print*,'fluxtop ir=',fluxtop_lw
816          !print*,'fluxsurf ir=',fluxsurf_lw
817
818!         write(*,*) 'pressure, eps_nlte_sw, eps_nlte_lw'
819!         do l=1,nlayer
820!           write(*,*)pplay(1,l),eps_nlte_sw(1,l),eps_nlte_lw(1,l)
821!         end do
822
823         !endif
824
825!     ---------------------------------------------------------------
826      end do                    ! end of big loop over every GCM column (ig = 1:ngrid)
827
828!-----------------------------------------------------------------------
829!     Additional diagnostics
830
831!     IR spectral output, for exoplanet observational comparison
832    !   if(specOLR)then
833    !     if(ngrid.ne.1)then
834    !       call writediagspec(ngrid,"OLR3D", &
835    !          "OLR(lon,lat,band)","W m^-2",3,OLR_nu)
836    !     endif
837    !   endif
838
839      if(lastcall)then
840
841        ! 1D Output
842        if(ngrid.eq.1)then
843
844!     surface diagnotics
845           diagrad_surf=.true.
846           if(diagrad_surf)then
847               open(116,file='surf_vals.out')
848               write(116,*) tsurf(1),pplev(1,1),    &
849                  fluxtop_dn(1) - fluxtop_sw(1),fluxtop_lw(1)
850               do nw=1,L_NSPECTV
851                 write(116,*) wavev(nw),fmnetv_nu(L_NLAYRAD,nw)
852               enddo
853               do nw=1,L_NSPECTI
854                 write(116,*) wavei(nw),fmneti_nu(L_NLAYRAD,nw)
855               enddo
856               close(116)
857           endif
858
859!     OLR by band
860           diagrad_OLR=.true.
861           if(diagrad_OLR)then
862               open(117,file='OLRnu.out')
863               write(117,*) 'IR wavel - band width - OLR'
864               do nw=1,L_NSPECTI
865                   write(117,*) wavei(nw),  &
866             abs(1.e4/bwnv(nw)-1.e4/bwnv(nw+1)),OLR_nu(1,nw)
867               enddo
868               close(117)
869           endif
870
871!     OLR vs altitude: in a .txt file
872           diagrad_OLRz=.true.
873           if(diagrad_OLRz)then
874               open(118,file='OLRz_plevs.out')
875               open(119,file='OLRz.out')
876               do l=1,L_NLAYRAD
877                 write(118,*) plevrad(2*l)
878                 do nw=1,L_NSPECTI
879                     write(119,*) fluxupi_nu(l,nw)
880                 enddo
881               enddo
882               close(118)
883               close(119)
884           endif
885
886!     Heating rates vs altitude in a .txt file
887           diagrad_rates=.true.
888           if(diagrad_rates)then
889             open(120,file='heating_rates.out')
890             write(120,*) "Pressure - Alt - HR tot - Rates (wavel SW)"
891             do l=1,nlayer
892               write(120,*) pplay(1,l),zzlay(1,l),dtsw(1,l),dtsw_nu(l,:)
893             enddo
894             close(120)
895
896             open(121,file='cooling_rates.out')
897             write(121,*) "Pressure - Alt - CR tot - Rates (wavel LW)"
898             do l=1,nlayer
899               write(121,*) pplay(1,l),zzlay(1,l),dtlw(1,l),dtlw_nu(l,:)
900             enddo
901             close(121)
902
903             open(122,file='bands.out')
904             write(122,*) "wavel - bands boundaries (microns)"
905             do nw=1,L_NSPECTV
906               write(122,*) wavev(nw),1.e4/bwnv(nw+1),1.e4/bwnv(nw)
907             enddo
908             do nw=1,L_NSPECTI
909               write(122,*) wavei(nw),1.e4/bwni(nw+1),1.e4/bwni(nw)
910             enddo
911             close(122)
912
913             open(123,file='c2h2_rates.out')
914             write(123,*) "Pressure - Alt - CR c2h2"
915             do l=1,nlayer
916               write(123,*) pplay(1,l),zzlay(1,l),dtlw_hcn_c2h2(1,l)
917             enddo
918             close(123)
919
920           endif
921
922        endif ! ngrid.eq.1
923
924      endif ! lastcall
925
926    end subroutine callcorrk_pluto
927
928END MODULE callcorrk_pluto_mod
Note: See TracBrowser for help on using the repository browser.