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

Last change on this file since 2176 was 2133, checked in by jvatant, 6 years ago

Fic a bug in r2131 the writediagfi cannot be called
in RT routines because of the iradia it goes crazy with outputs freq ...
--JVO

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