source: trunk/LMDZ.PLUTO/libf/phypluto/optci.F90 @ 3879

Last change on this file since 3879 was 3695, checked in by debatzbr, 9 months ago

Calculate the visible opacity everywhere
BBT

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