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

Last change on this file since 1974 was 1725, checked in by jvatant, 7 years ago

Optimization of the optci/cv routines

  • The repeated calls to huge matrices gasi/v increased dramatically the execution time because of memory access
  • Added a tmpk variable
  • Save ~ 50% time on the RT, ~30% on the whole code on the tested simulations

JVO

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