source: trunk/LMDZ.GENERIC/libf/phystd/optci.F90 @ 2582

Last change on this file since 2582 was 2582, checked in by emillour, 3 years ago

Generic GCM:
Fixed bug in tpindex (for low temperatures, between first and second
reference temperatures, temperature was wrongly set to tref(1)).
The input temperature was also allowed to be modified by the routine,
which is probably not a good idea and no longer the case.
Took this opportunity to turn tpindex into a module.
GM+EM

File size: 15.0 KB
Line 
1MODULE optci_mod
2
3IMPLICIT NONE
4
5CONTAINS
6
7subroutine optci(PLEV,TLEV,DTAUI,TAUCUMI,      &
8     QXIAER,QSIAER,GIAER,COSBI,WBARI,TAUAERO,  &
9     TMID,PMID,TAUGSURF,QVAR,MUVAR)
10
11  use radinc_h, only: L_LEVELS, L_NLAYRAD, L_NSPECTI, L_NGAUSS, &
12                      L_NLEVRAD, L_REFVAR, naerkind
13  use radcommon_h, only: gasi,tlimit,wrefVAR,Cmk,tgasref,pfgasref,wnoi,scalep,indi,glat_ig
14  use gases_h, only: gfrac, ngasmx, igas_N2, igas_He, igas_H2O, igas_H2
15  use comcstfi_mod, only: g, r, mugaz
16  use callkeys_mod, only: kastprof,continuum,graybody
17  use recombin_corrk_mod, only: corrk_recombin, gasi_recomb
18  use tpindex_mod, only: tpindex
19
20  implicit none
21
22  !==================================================================
23  !     
24  !     Purpose
25  !     -------
26  !     Calculates longwave optical constants at each level. For each
27  !     layer and spectral interval in the IR it calculates WBAR, DTAU
28  !     and COSBAR. For each level it calculates TAU.
29  !     
30  !     TAUCUMI(L,LW) is the cumulative optical depth at level L (or alternatively
31  !     at the *bottom* of layer L), LW is the spectral wavelength interval.
32  !     
33  !     TLEV(L) - Temperature at the layer boundary (i.e., level)
34  !     PLEV(L) - Pressure at the layer boundary (i.e., level)
35  !
36  !     Authors
37  !     -------
38  !     Adapted from the NASA Ames code by R. Wordsworth (2009)
39  !     
40  !==================================================================
41
42
43  real*8 DTAUI(L_NLAYRAD,L_NSPECTI,L_NGAUSS)
44  real*8 DTAUKI(L_LEVELS,L_NSPECTI,L_NGAUSS)
45  real*8 TAUI(L_NLEVRAD,L_NSPECTI,L_NGAUSS)
46  real*8 TAUCUMI(L_LEVELS,L_NSPECTI,L_NGAUSS)
47  real*8 PLEV(L_LEVELS)
48  real*8 TLEV(L_LEVELS)
49  real*8 TMID(L_LEVELS), PMID(L_LEVELS)
50  real*8 COSBI(L_NLAYRAD,L_NSPECTI,L_NGAUSS)
51  real*8 WBARI(L_NLAYRAD,L_NSPECTI,L_NGAUSS)
52
53  ! for aerosols
54  real*8  QXIAER(L_LEVELS,L_NSPECTI,NAERKIND)
55  real*8  QSIAER(L_LEVELS,L_NSPECTI,NAERKIND)
56  real*8  GIAER(L_LEVELS,L_NSPECTI,NAERKIND)
57  real*8  TAUAERO(L_LEVELS,NAERKIND)
58  real*8  TAUAEROLK(L_LEVELS,L_NSPECTI,NAERKIND)
59  real*8  TAEROS(L_LEVELS,L_NSPECTI,NAERKIND)
60
61  integer L, NW, NG, K, LK, IAER
62  integer MT(L_LEVELS), MP(L_LEVELS), NP(L_LEVELS)
63  real*8  ANS, TAUGAS
64  real*8  DPR(L_LEVELS), U(L_LEVELS)
65  real*8  LCOEF(4), LKCOEF(L_LEVELS,4)
66
67  real*8 taugsurf(L_NSPECTI,L_NGAUSS-1)
68  real*8 DCONT,DAERO
69  double precision wn_cont, p_cont, p_air, T_cont, dtemp, dtempc
70  double precision p_cross
71
72  ! variable species mixing ratio variables
73  real*8  QVAR(L_LEVELS), WRATIO(L_LEVELS), MUVAR(L_LEVELS)
74  real*8  KCOEF(4)
75  integer NVAR(L_LEVELS)
76 
77  ! temporary variables to reduce memory access time to gasi
78  real*8 tmpk(2,2)
79  real*8 tmpkvar(2,2,2)
80
81  ! temporary variables for multiple aerosol calculation
82  real*8 atemp
83  real*8 btemp(L_NLAYRAD,L_NSPECTI)
84
85  ! variables for k in units m^-1
86  real*8 dz(L_LEVELS)
87  !real*8 rho !! see test below
88
89  integer igas, jgas
90
91  integer interm
92
93  !--- Kasting's CIA ----------------------------------------
94  !real*8, parameter :: Ci(L_NSPECTI)=[                         &
95  !     3.8E-5, 1.2E-5, 2.8E-6, 7.6E-7, 4.5E-7, 2.3E-7,    &
96  !     5.4E-7, 1.6E-6, 0.0,                               &
97  !     0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0,            &
98  !     0.0, 4.0E-7, 4.0E-6, 1.4E-5,    &
99  !     1.0E-5, 1.2E-6, 2.0E-7, 5.0E-8, 3.0E-8, 0.0 ]
100  !real*8, parameter :: Ti(L_NSPECTI)=[ -2.2, -1.9,             &
101  !     -1.7, -1.7, -1.7, -1.7, -1.7, -1.7,                &
102  !     0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, &
103  !     -1.7,-1.7,-1.7,-1.7,-1.7,-1.7,-1.7, -1.7,0.0 ]
104  !----------------------------------------------------------
105
106  !! AS: to save time in computing continuum (see bilinearbig)
107  IF (.not.ALLOCATED(indi)) THEN
108      ALLOCATE(indi(L_NSPECTI,ngasmx,ngasmx))
109      indi = -9999  ! this initial value means "to be calculated"
110  ENDIF
111
112  !=======================================================================
113  !     Determine the total gas opacity throughout the column, for each
114  !     spectral interval, NW, and each Gauss point, NG.
115
116  taugsurf(:,:) = 0.0
117  dpr(:)        = 0.0
118  lkcoef(:,:)   = 0.0
119
120  do K=2,L_LEVELS
121     DPR(k) = PLEV(K)-PLEV(K-1)
122
123     !--- Kasting's CIA ----------------------------------------
124     !dz(k)=dpr(k)*189.02*TMID(K)/(0.03720*PMID(K))
125     ! this is CO2 path length (in cm) as written by Francois
126     ! delta_z = delta_p * R_specific * T / (g * P)
127     ! But Kasting states that W is in units of _atmosphere_ cm
128     ! So we do
129     !dz(k)=dz(k)*(PMID(K)/1013.25)
130     !dz(k)=dz(k)/100.0 ! in m for SI calc
131     !----------------------------------------------------------
132
133     ! if we have continuum opacities, we need dz
134     if(kastprof)then
135        dz(k) = dpr(k)*(1000.0d0*8.3145d0/muvar(k))*TMID(K)/(g*PMID(K))
136        U(k)  = Cmk*DPR(k)*mugaz/muvar(k)
137     else
138        dz(k) = dpr(k)*R*TMID(K)/(glat_ig*PMID(K))*mugaz/muvar(k)
139        U(k)  = Cmk*DPR(k)*mugaz/muvar(k)     ! only Cmk line in optci.F 
140            !JL13 the mugaz/muvar factor takes into account water meanmolecular weight if water is present
141     endif
142
143     call tpindex(PMID(K),TMID(K),QVAR(K),pfgasref,tgasref,WREFVAR, &
144          LCOEF,MT(K),MP(K),NVAR(K),WRATIO(K))
145
146     do LK=1,4
147        LKCOEF(K,LK) = LCOEF(LK)
148     end do
149  end do                    ! levels
150
151  ! Spectral dependance of aerosol absorption
152  do iaer=1,naerkind
153     DO NW=1,L_NSPECTI
154        do K=2,L_LEVELS
155           TAEROS(K,NW,IAER) = TAUAERO(K,IAER) * QXIAER(K,NW,IAER)
156        end do                    ! levels
157     END DO
158  end do
159
160  do NW=1,L_NSPECTI
161
162     do K=2,L_LEVELS
163     
164        DAERO=SUM(TAEROS(K,NW,1:naerkind)) ! aerosol absorption
165
166        DCONT = 0.0d0 ! continuum absorption
167
168        if(continuum.and.(.not.graybody))then
169           ! include continua if necessary
170           wn_cont = dble(wnoi(nw))
171           T_cont  = dble(TMID(k))
172           do igas=1,ngasmx
173
174              if(gfrac(igas).eq.-1)then ! variable
175                 p_cont  = dble(PMID(k)*scalep*QVAR(k)) ! qvar = mol/mol
176              else
177                 p_cont  = dble(PMID(k)*scalep*gfrac(igas)*(1.-QVAR(k)))
178              endif
179
180              dtemp=0.0d0
181              if(igas.eq.igas_N2)then
182
183                 interm = indi(nw,igas,igas)
184                 call interpolateN2N2(wn_cont,T_cont,p_cont,dtemp,.false.,interm)
185                 indi(nw,igas,igas) = interm
186
187              elseif(igas.eq.igas_H2)then
188
189                 ! first do self-induced absorption
190                 interm = indi(nw,igas,igas)
191                 call interpolateH2H2(wn_cont,T_cont,p_cont,dtemp,.false.,interm)
192                 indi(nw,igas,igas) = interm
193
194                 ! then cross-interactions with other gases
195                 do jgas=1,ngasmx
196                    p_cross = dble(PMID(k)*scalep*gfrac(jgas)*(1.-QVAR(k)))
197                    dtempc  = 0.0d0
198                    if(jgas.eq.igas_N2)then
199                       interm = indi(nw,igas,jgas)
200                       call interpolateN2H2(wn_cont,T_cont,p_cross,p_cont,dtempc,.false.,interm)
201                       indi(nw,igas,jgas) = interm
202                    elseif(jgas.eq.igas_He)then
203                       interm = indi(nw,igas,jgas)
204                       call interpolateH2He(wn_cont,T_cont,p_cross,p_cont,dtempc,.false.,interm)
205                       indi(nw,igas,jgas) = interm
206                    endif
207                    dtemp = dtemp + dtempc
208                 enddo
209
210              elseif(igas.eq.igas_H2O.and.T_cont.gt.100.0)then
211                 ! Compute self and foreign (with air) continuum of H2O
212                 p_air = dble(PMID(k)*scalep) - p_cont ! note assumes background is air!
213                 interm = indi(nw,igas,igas)
214                 call interpolateH2O_self_foreign(wn_cont,T_cont,p_cont,p_air,dtemp,.false.,interm) ! MTCKD v3.3
215                 indi(nw,igas,igas) = interm
216
217              endif
218
219              DCONT = DCONT + dtemp
220
221           enddo
222
223           ! Oobleck test
224           !rho = PMID(k)*scalep / (TMID(k)*286.99)
225           !if(WNOI(nw).gt.300.0 .and. WNOI(nw).lt.500.0)then
226           !   DCONT = rho * 0.125 * 4.6e-4
227           !elseif(WNOI(nw).gt.500.0 .and. WNOI(nw).lt.700.0)then
228           !   DCONT = 1000*dpr(k) * 1.0 * 4.6e-4 / g
229           !   DCONT = rho * 1.0 * 4.6e-4
230           !elseif(WNOI(nw).gt.700.0 .and. WNOI(nw).lt.900.0)then
231           !   DCONT = rho * 0.125 * 4.6e-4
232           !endif
233
234           DCONT = DCONT*dz(k)
235
236        endif
237
238        do ng=1,L_NGAUSS-1
239
240           ! Now compute TAUGAS
241
242           ! Interpolate between water mixing ratios
243           ! WRATIO = 0.0 if the requested water amount is equal to, or outside the
244           ! the water data range
245
246           if(L_REFVAR.eq.1)then ! added by RW for special no variable case
247           
248              ! JVO 2017 : added tmpk because the repeated calls to gasi/v increased dramatically
249              ! the execution time of optci/v -> ~ factor 2 on the whole radiative
250              ! transfer on the tested simulations !
251
252              IF (corrk_recombin) THEN ! added by JVO
253                tmpk = GASI_RECOMB(MT(K):MT(K)+1,MP(K):MP(K)+1,1,NW,NG) ! contains the mix of recombined species
254              ELSE
255                tmpk = GASI(MT(K):MT(K)+1,MP(K):MP(K)+1,1,NW,NG)
256              ENDIF
257
258              KCOEF(1) = tmpk(1,1) ! KCOEF(1) = GASI(MT(K),MP(K),1,NW,NG)
259              KCOEF(2) = tmpk(1,2) ! KCOEF(2) = GASI(MT(K),MP(K)+1,1,NW,NG)
260              KCOEF(3) = tmpk(2,2) ! KCOEF(3) = GASI(MT(K)+1,MP(K)+1,1,NW,NG)
261              KCOEF(4) = tmpk(2,1) ! KCOEF(4) = GASI(MT(K)+1,MP(K),1,NW,NG)
262
263           else
264
265              IF (corrk_recombin) THEN ! added by JVO
266                tmpkvar = GASI_RECOMB(MT(K):MT(K)+1,MP(K):MP(K)+1,NVAR(K):NVAR(K)+1,NW,NG)
267              ELSE
268                tmpkvar = GASI(MT(K):MT(K)+1,MP(K):MP(K)+1,NVAR(K):NVAR(K)+1,NW,NG)
269              ENDIF
270
271              KCOEF(1) = tmpkvar(1,1,1) + WRATIO(K) *  &
272                        ( tmpkvar(1,1,2)-tmpkvar(1,1,1) )
273
274              KCOEF(2) = tmpkvar(1,2,1) + WRATIO(K) *  &
275                        ( tmpkvar(1,2,2)-tmpkvar(1,2,1) )
276
277              KCOEF(3) = tmpkvar(2,2,1) + WRATIO(K) *  &
278                        ( tmpkvar(2,2,2)-tmpkvar(2,2,1) )
279             
280              KCOEF(4) = tmpkvar(2,1,1) + WRATIO(K) *  &
281                        ( tmpkvar(2,1,2)-tmpkvar(2,1,1) )
282
283           endif
284
285           ! Interpolate the gaseous k-coefficients to the requested T,P values
286
287           ANS = LKCOEF(K,1)*KCOEF(1) + LKCOEF(K,2)*KCOEF(2) +            &
288                LKCOEF(K,3)*KCOEF(3) + LKCOEF(K,4)*KCOEF(4)
289
290           TAUGAS  = U(k)*ANS
291
292           TAUGSURF(NW,NG) = TAUGSURF(NW,NG) + TAUGAS + DCONT
293           DTAUKI(K,nw,ng) = TAUGAS    &
294                             + DCONT   & ! For parameterized continuum absorption
295                             + DAERO     ! For aerosol absorption
296
297        end do
298
299        ! Now fill in the "clear" part of the spectrum (NG = L_NGAUSS),
300        ! which holds continuum opacity only
301
302        NG              = L_NGAUSS
303        DTAUKI(K,nw,ng) = 0.d0      &
304                          + DCONT   & ! For parameterized continuum absorption
305                          + DAERO     ! For aerosol absorption
306
307     end do
308  end do
309
310  !=======================================================================
311  !     Now the full treatment for the layers, where besides the opacity
312  !     we need to calculate the scattering albedo and asymmetry factors
313
314  do iaer=1,naerkind
315    DO NW=1,L_NSPECTI
316     DO K=2,L_LEVELS
317           TAUAEROLK(K,NW,IAER) = TAUAERO(K,IAER)*QSIAER(K,NW,IAER) ! effect of scattering albedo
318     ENDDO
319    ENDDO
320  end do
321 
322  DO NW=1,L_NSPECTI
323     DO L=1,L_NLAYRAD-1
324        K              = 2*L+1
325        btemp(L,NW) = SUM(TAUAEROLK(K,NW,1:naerkind)) + SUM(TAUAEROLK(K+1,NW,1:naerkind))
326     END DO ! L vertical loop
327     
328     ! Last level
329     L           = L_NLAYRAD
330     K           = 2*L+1   
331     btemp(L,NW) = SUM(TAUAEROLK(K,NW,1:naerkind))
332     
333  END DO                    ! NW spectral loop
334 
335
336  DO NW=1,L_NSPECTI
337     NG = L_NGAUSS
338     DO L=1,L_NLAYRAD-1
339
340        K              = 2*L+1
341        DTAUI(L,nw,ng) = DTAUKI(K,NW,NG) + DTAUKI(K+1,NW,NG)! + 1.e-50
342
343        atemp = 0.
344        if(DTAUI(L,NW,NG) .GT. 1.0D-9) then
345           do iaer=1,naerkind
346              atemp = atemp +                                     &
347                   GIAER(K,NW,IAER)   * TAUAEROLK(K,NW,IAER) +    &
348                   GIAER(K+1,NW,IAER) * TAUAEROLK(K+1,NW,IAER)
349           end do
350           WBARI(L,nw,ng) = btemp(L,nw)  / DTAUI(L,NW,NG)
351        else
352           WBARI(L,nw,ng) = 0.0D0
353           DTAUI(L,NW,NG) = 1.0D-9
354        endif
355
356        if(btemp(L,nw) .GT. 0.0d0) then
357           cosbi(L,NW,NG) = atemp/btemp(L,nw)
358        else
359           cosbi(L,NW,NG) = 0.0D0
360        end if
361
362     END DO ! L vertical loop
363     
364     ! Last level
365     
366     L              = L_NLAYRAD
367     K              = 2*L+1
368     DTAUI(L,nw,ng) = DTAUKI(K,NW,NG) ! + 1.e-50
369
370     atemp = 0.
371     if(DTAUI(L,NW,NG) .GT. 1.0D-9) then
372        do iaer=1,naerkind
373           atemp = atemp + GIAER(K,NW,IAER)   * TAUAEROLK(K,NW,IAER)
374        end do
375        WBARI(L,nw,ng) = btemp(L,nw)  / DTAUI(L,NW,NG)
376     else
377        WBARI(L,nw,ng) = 0.0D0
378        DTAUI(L,NW,NG) = 1.0D-9
379     endif
380
381     if(btemp(L,nw) .GT. 0.0d0) then
382        cosbi(L,NW,NG) = atemp/btemp(L,nw)
383     else
384        cosbi(L,NW,NG) = 0.0D0
385     end if
386     
387
388     ! Now the other Gauss points, if needed.
389
390     DO NG=1,L_NGAUSS-1
391        IF(TAUGSURF(NW,NG) .gt. TLIMIT) THEN
392
393           DO L=1,L_NLAYRAD-1
394              K              = 2*L+1
395              DTAUI(L,nw,ng) = DTAUKI(K,NW,NG)+DTAUKI(K+1,NW,NG)! + 1.e-50
396
397              if(DTAUI(L,NW,NG) .GT. 1.0D-9) then
398
399                 WBARI(L,nw,ng) = btemp(L,nw)  / DTAUI(L,NW,NG)
400
401              else
402                 WBARI(L,nw,ng) = 0.0D0
403                 DTAUI(L,NW,NG) = 1.0D-9
404              endif
405
406              cosbi(L,NW,NG) = cosbi(L,NW,L_NGAUSS)
407           END DO ! L vertical loop
408           
409           ! Last level
410           L              = L_NLAYRAD
411           K              = 2*L+1
412           DTAUI(L,nw,ng) = DTAUKI(K,NW,NG)! + 1.e-50
413
414           if(DTAUI(L,NW,NG) .GT. 1.0D-9) then
415
416              WBARI(L,nw,ng) = btemp(L,nw)  / DTAUI(L,NW,NG)
417
418           else
419              WBARI(L,nw,ng) = 0.0D0
420              DTAUI(L,NW,NG) = 1.0D-9
421           endif
422
423           cosbi(L,NW,NG) = cosbi(L,NW,L_NGAUSS)
424           
425        END IF
426
427     END DO                 ! NG Gauss loop
428  END DO                    ! NW spectral loop
429
430  ! Total extinction optical depths
431
432  DO NG=1,L_NGAUSS       ! full gauss loop
433     DO NW=1,L_NSPECTI       
434        TAUCUMI(1,NW,NG)=0.0D0
435        DO K=2,L_LEVELS
436           TAUCUMI(K,NW,NG)=TAUCUMI(K-1,NW,NG)+DTAUKI(K,NW,NG)
437        END DO
438     END DO                 ! end full gauss loop
439  END DO
440
441  ! be aware when comparing with textbook results
442  ! (e.g. Pierrehumbert p. 218) that
443  ! taucumi does not take the <cos theta>=0.5 factor into
444  ! account. It is the optical depth for a vertically
445  ! ascending ray with angle theta = 0.
446
447  !open(127,file='taucum.out')
448  !do nw=1,L_NSPECTI
449  !   write(127,*) taucumi(L_LEVELS,nw,L_NGAUSS)
450  !enddo
451  !close(127)
452 
453!  print*,'WBARI'
454!  print*,WBARI
455!  print*,'DTAUI'
456!  print*,DTAUI
457!  call abort
458
459end subroutine optci
460
461END MODULE optci_mod
462
Note: See TracBrowser for help on using the repository browser.