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

Last change on this file since 3501 was 3501, checked in by afalco, 13 days ago

Pluto PCM: removed useless code.
Reindenting.
AF

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