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

Last change on this file since 3957 was 3889, checked in by debatzbr, 4 months ago

Pluto PCM: Pass haze production parameters as keys in callphys.def
BBT

File size: 17.4 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,callmufi,Fabs_aers_IR,Fabs_aerf_IR
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  ! Variables for aerosol absorption
97  real*8 Fabs_aer(NAERKIND)
98
99  integer igas, jgas
100
101  integer interm
102 
103  logical :: firstcall=.true.
104!$OMP THREADPRIVATE(firstcall)
105
106  !--- Kasting's CIA ----------------------------------------
107  !real*8, parameter :: Ci(L_NSPECTI)=[                         &
108  !     3.8E-5, 1.2E-5, 2.8E-6, 7.6E-7, 4.5E-7, 2.3E-7,    &
109  !     5.4E-7, 1.6E-6, 0.0,                               &
110  !     0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0,            &
111  !     0.0, 4.0E-7, 4.0E-6, 1.4E-5,    &
112  !     1.0E-5, 1.2E-6, 2.0E-7, 5.0E-8, 3.0E-8, 0.0 ]
113  !real*8, parameter :: Ti(L_NSPECTI)=[ -2.2, -1.9,             &
114  !     -1.7, -1.7, -1.7, -1.7, -1.7, -1.7,                &
115  !     0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, &
116  !     -1.7,-1.7,-1.7,-1.7,-1.7,-1.7,-1.7, -1.7,0.0 ]
117  !----------------------------------------------------------
118
119  if (firstcall) then
120    ! allocate local arrays of size "naerkind" (which are also
121    ! "saved" so that this is done only once in for all even if
122    ! we don't need to store the value from a time step to the next)
123    allocate(TAUAEROLK(L_LEVELS,L_NSPECTI,NAERKIND))
124    allocate(TAEROS(L_LEVELS,L_NSPECTI,NAERKIND))
125    firstcall=.false.
126  endif ! of if (firstcall)
127
128  !! AS: to save time in computing continuum (see bilinearbig)
129  IF (.not.ALLOCATED(indi)) THEN
130      ALLOCATE(indi(L_NSPECTI,ngasmx,ngasmx))
131      indi = -9999  ! this initial value means "to be calculated"
132  ENDIF
133
134  !=======================================================================
135  !     Determine the total gas opacity throughout the column, for each
136  !     spectral interval, NW, and each Gauss point, NG.
137
138  dtaui(:,:,:) = 0.0
139  taucumi(:,:,:) = 0.0
140
141  taugsurf(:,:) = 0.0
142  dpr(:)        = 0.0
143  lkcoef(:,:)   = 0.0
144
145  if(callmufi) then
146   Fabs_aer(1) = Fabs_aers_IR
147   Fabs_aer(2) = Fabs_aerf_IR
148  else
149   Fabs_aer(:) = 1.0
150  endif
151
152  do K=2,L_LEVELS
153     DPR(k) = PLEV(K)-PLEV(K-1)
154
155     !--- Kasting's CIA ----------------------------------------
156     !dz(k)=dpr(k)*189.02*TMID(K)/(0.03720*PMID(K))
157     ! this is N2 path length (in cm) as written by Francois
158     ! delta_z = delta_p * R_specific * T / (g * P)
159     ! But Kasting states that W is in units of _atmosphere_ cm
160     ! So we do
161     !dz(k)=dz(k)*(PMID(K)/1013.25)
162     !dz(k)=dz(k)/100.0 ! in m for SI calc
163     !----------------------------------------------------------
164
165     ! if we have continuum opacities, we need dz
166     if(kastprof)then
167        dz(k) = dpr(k)*(1000.0d0*8.3145d0/muvar(k))*TMID(K)/(g*PMID(K))
168        U(k)  = Cmk*DPR(k)*mugaz/muvar(k)
169     else
170        dz(k) = dpr(k)*R*TMID(K)/(glat_ig*PMID(K))*mugaz/muvar(k)
171        U(k)  = Cmk*DPR(k)*mugaz/muvar(k)     ! only Cmk line in optci.F 
172            !JL13 the mugaz/muvar factor takes into account water meanmolecular weight if water is present
173     endif
174
175     call tpindex(PMID(K),TMID(K),QVAR(K),pfgasref,tgasref,WREFVAR, &
176          LCOEF,MT(K),MP(K),NVAR(K),WRATIO(K))
177
178     do LK=1,4
179        LKCOEF(K,LK) = LCOEF(LK)
180     end do
181  end do                    ! levels
182
183  ! Spectral dependance of aerosol absorption
184  do iaer=1,naerkind
185     DO NW=1,L_NSPECTI
186        do K=2,L_LEVELS
187           TAEROS(K,NW,IAER) = TAUAERO(K,IAER) * QXIAER(K,NW,IAER)
188           ! Aerosol absorption correction --> impact on opacity.
189           if (callmufi .and. (TAEROS(K,NW,IAER).gt.0.d0)) then
190            TAEROS(K,NW,IAER) = TAEROS(K,NW,IAER) * &
191                                ((QSIAER(K,NW,IAER)/QXIAER(K,NW,IAER)) + Fabs_aer(IAER)*(1.-(QSIAER(K,NW,IAER)/QXIAER(K,NW,IAER))))
192           endif
193        end do ! L_LEVELS
194     END DO ! L_NSPECTI
195  end do ! naerkind
196
197  do NW=1,L_NSPECTI
198
199     do K=2,L_LEVELS
200     
201        DAERO=SUM(TAEROS(K,NW,1:naerkind)) ! aerosol absorption
202
203        DCONT = 0.0d0 ! continuum absorption
204
205        if(continuum.and.(.not.graybody))then
206           ! include continua if necessary
207           wn_cont = dble(wnoi(nw))
208           T_cont  = dble(TMID(k))
209           do igas=1,ngasmx
210
211              if(gfrac(igas).eq.-1)then ! variable
212                 p_cont  = dble(PMID(k)*scalep*QVAR(k)) ! qvar = mol/mol
213              else
214                 p_cont  = dble(PMID(k)*scalep*gfrac(igas)*(1.-QVAR(k)))
215              endif
216
217              dtemp=0.0d0
218              if(igas.eq.igas_N2)then
219
220                 interm = indi(nw,igas,igas)
221                 call interpolateN2N2(wn_cont,T_cont,p_cont,dtemp,.false.,interm)
222                 indi(nw,igas,igas) = interm
223
224            !   elseif(igas.eq.igas_H2)then !AF24: removed
225
226              elseif(igas.eq.igas_CH4)then
227
228                 ! first do self-induced absorption
229                 interm = indi(nw,igas,igas)
230                 call interpolateCH4CH4(wn_cont,T_cont,p_cont,dtemp,.false.,interm)
231                 indi(nw,igas,igas) = interm
232
233                 ! then cross-interactions with other gases !AF24: removed
234
235            !   elseif(igas.eq.igas_H2O.and.T_cont.gt.100.0)then !AF24: removed
236
237              endif
238              DCONT = DCONT + dtemp
239           enddo
240
241           ! Oobleck test
242           !rho = PMID(k)*scalep / (TMID(k)*286.99)
243           !if(WNOI(nw).gt.300.0 .and. WNOI(nw).lt.500.0)then
244           !   DCONT = rho * 0.125 * 4.6e-4
245           !elseif(WNOI(nw).gt.500.0 .and. WNOI(nw).lt.700.0)then
246           !   DCONT = 1000*dpr(k) * 1.0 * 4.6e-4 / g
247           !   DCONT = rho * 1.0 * 4.6e-4
248           !elseif(WNOI(nw).gt.700.0 .and. WNOI(nw).lt.900.0)then
249           !   DCONT = rho * 0.125 * 4.6e-4
250           !endif
251
252           DCONT = DCONT*dz(k)
253
254        endif
255
256        do ng=1,L_NGAUSS-1
257
258           ! Now compute TAUGAS
259
260           ! Interpolate between water mixing ratios
261           ! WRATIO = 0.0 if the requested water amount is equal to, or outside the
262           ! the water data range
263
264           if(L_REFVAR.eq.1)then ! added by RW for special no variable case
265           
266              ! JVO 2017 : added tmpk because the repeated calls to gasi/v increased dramatically
267              ! the execution time of optci/v -> ~ factor 2 on the whole radiative
268              ! transfer on the tested simulations !
269
270              IF (corrk_recombin) THEN ! added by JVO
271                tmpk = GASI_RECOMB(MT(K):MT(K)+1,MP(K):MP(K)+1,1,NW,NG) ! contains the mix of recombined species
272              ELSE
273                tmpk = GASI(MT(K):MT(K)+1,MP(K):MP(K)+1,1,NW,NG)
274              ENDIF
275
276              KCOEF(1) = tmpk(1,1) ! KCOEF(1) = GASI(MT(K),MP(K),1,NW,NG)
277              KCOEF(2) = tmpk(1,2) ! KCOEF(2) = GASI(MT(K),MP(K)+1,1,NW,NG)
278              KCOEF(3) = tmpk(2,2) ! KCOEF(3) = GASI(MT(K)+1,MP(K)+1,1,NW,NG)
279              KCOEF(4) = tmpk(2,1) ! KCOEF(4) = GASI(MT(K)+1,MP(K),1,NW,NG)
280
281           else
282
283              IF (corrk_recombin) THEN ! added by JVO
284                tmpkvar = GASI_RECOMB(MT(K):MT(K)+1,MP(K):MP(K)+1,NVAR(K):NVAR(K)+1,NW,NG)
285              ELSE
286                tmpkvar = GASI(MT(K):MT(K)+1,MP(K):MP(K)+1,NVAR(K):NVAR(K)+1,NW,NG)
287              ENDIF
288
289              KCOEF(1) = tmpkvar(1,1,1) + WRATIO(K) *  &
290                        ( tmpkvar(1,1,2)-tmpkvar(1,1,1) )
291
292              KCOEF(2) = tmpkvar(1,2,1) + WRATIO(K) *  &
293                        ( tmpkvar(1,2,2)-tmpkvar(1,2,1) )
294
295              KCOEF(3) = tmpkvar(2,2,1) + WRATIO(K) *  &
296                        ( tmpkvar(2,2,2)-tmpkvar(2,2,1) )
297             
298              KCOEF(4) = tmpkvar(2,1,1) + WRATIO(K) *  &
299                        ( tmpkvar(2,1,2)-tmpkvar(2,1,1) )
300
301           endif
302
303           ! Interpolate the gaseous k-coefficients to the requested T,P values
304
305           ANS = LKCOEF(K,1)*KCOEF(1) + LKCOEF(K,2)*KCOEF(2) +            &
306                LKCOEF(K,3)*KCOEF(3) + LKCOEF(K,4)*KCOEF(4)
307
308           TAUGAS  = U(k)*ANS
309
310           TAUGSURF(NW,NG) = TAUGSURF(NW,NG) + TAUGAS + DCONT
311           DTAUKI(K,nw,ng) = TAUGAS    &
312                             + DCONT   & ! For parameterized continuum absorption
313                             + DAERO     ! For aerosol absorption
314
315        end do
316
317        ! Now fill in the "clear" part of the spectrum (NG = L_NGAUSS),
318        ! which holds continuum opacity only
319
320        NG              = L_NGAUSS
321        DTAUKI(K,nw,ng) = 0.d0      &
322                          + DCONT   & ! For parameterized continuum absorption
323                          + DAERO     ! For aerosol absorption
324
325     end do
326  end do
327
328  !=======================================================================
329  !     Now the full treatment for the layers, where besides the opacity
330  !     we need to calculate the scattering albedo and asymmetry factors
331
332  do iaer=1,naerkind
333    DO NW=1,L_NSPECTI
334     DO K=2,L_LEVELS
335           TAUAEROLK(K,NW,IAER) = TAUAERO(K,IAER)*QSIAER(K,NW,IAER) ! effect of scattering albedo
336     ENDDO
337    ENDDO
338  end do
339
340  DO NW = 1, L_NSPECTI
341     DO L = 1, L_NLAYRAD-1
342        K = 2*L + 1
343        ! Aerosol absorption correction --> impact on single scattering albedo.
344        if (callmufi) then
345           if ((TAEROS(K,NW,1).gt.0.d0) .and. TAEROS(K+1,NW,1).gt.0.d0) then
346              btemp(L,NW) = TAUAEROLK(K,NW,1) / &
347                          (QSIAER(K,NW,1)/QXIAER(K,NW,1) + &
348                             Fabs_aer(1)*(1.-QSIAER(K,NW,1)/QXIAER(K,NW,1))) + &
349                          TAUAEROLK(K+1,NW,1) / &
350                          (QSIAER(K+1,NW,1)/QXIAER(K+1,NW,1) + &
351                             Fabs_aer(1)*(1.-QSIAER(K+1,NW,1)/QXIAER(K+1,NW,1)))
352           else
353              btemp(L,NW) = TAUAEROLK(K,NW,1) + TAUAEROLK(K+1,NW,1)
354           endif
355           if ((TAEROS(K,NW,2).gt.0.d0) .and. TAEROS(K+1,NW,2).gt.0.d0) then
356              btemp(L,NW) = btemp(L,NW) + &
357                          (TAUAEROLK(K,NW,2) / &
358                          (QSIAER(K,NW,2)/QXIAER(K,NW,2) + &
359                             Fabs_aer(2)*(1.-QSIAER(K,NW,2)/QXIAER(K,NW,2)))) + &
360                          (TAUAEROLK(K+1,NW,2) / &
361                          (QSIAER(K+1,NW,2)/QXIAER(K+1,NW,2) + &
362                             Fabs_aer(2)*(1.-QSIAER(K+1,NW,2)/QXIAER(K+1,NW,2))))
363           else
364              btemp(L,NW) = btemp(L,NW) + TAUAEROLK(K,NW,2) + TAUAEROLK(K+1,NW,2)
365           endif
366        else
367           btemp(L,NW) = SUM(TAUAEROLK(K,NW,1:naerkind)) + SUM(TAUAEROLK(K+1,NW,1:naerkind))
368        endif ! callmufi
369     END DO ! L vertical loop
370     
371     ! Last level
372     !-----------
373     L = L_NLAYRAD
374     K = 2*L+1
375     ! Aerosol absorption correction --> impact on single scattering albedo.
376     if (callmufi) then
377      if (TAEROS(K,NW,1).gt.0.d0) then
378         btemp(L,NW) = TAUAEROLK(K,NW,1) / &
379                       (QSIAER(K,NW,1)/QXIAER(K,NW,1) + &
380                        Fabs_aer(1)*(1.-QSIAER(K,NW,1)/QXIAER(K,NW,1)))
381      else
382         btemp(L,NW) = TAUAEROLK(K,NW,1)
383      endif
384      if (TAEROS(K,NW,2).gt.0.d0) then
385         btemp(L,NW) = btemp(L,NW) + &
386                       (TAUAEROLK(K,NW,2) / &
387                       (QSIAER(K,NW,2)/QXIAER(K,NW,2) + &
388                        Fabs_aer(2)*(1.-QSIAER(K,NW,2)/QXIAER(K,NW,2))))
389      else
390         btemp(L,NW) = btemp(L,NW) + TAUAEROLK(K,NW,2)
391      endif
392     else
393      btemp(L,NW) = SUM(TAUAEROLK(K,NW,1:naerkind))
394     endif ! callmufi
395  END DO ! NW spectral loop
396 
397
398  DO NW=1,L_NSPECTI
399     NG = L_NGAUSS
400     DO L=1,L_NLAYRAD-1
401
402        K              = 2*L+1
403        DTAUI(L,nw,ng) = DTAUKI(K,NW,NG) + DTAUKI(K+1,NW,NG)! + 1.e-50
404
405        atemp = 0.
406        if(DTAUI(L,NW,NG) .GT. 1.0D-9) then
407           do iaer=1,naerkind
408              atemp = atemp +                                     &
409                   GIAER(K,NW,IAER)   * TAUAEROLK(K,NW,IAER) +    &
410                   GIAER(K+1,NW,IAER) * TAUAEROLK(K+1,NW,IAER)
411           end do
412           WBARI(L,nw,ng) = btemp(L,nw)  / DTAUI(L,NW,NG)
413        else
414           WBARI(L,nw,ng) = 0.0D0
415           DTAUI(L,NW,NG) = 1.0D-9
416        endif
417
418        if(btemp(L,nw) .GT. 0.0d0) then
419           cosbi(L,NW,NG) = atemp/btemp(L,nw)
420        else
421           cosbi(L,NW,NG) = 0.0D0
422        end if
423
424     END DO ! L vertical loop
425     
426     ! Last level
427     !-----------
428     L              = L_NLAYRAD
429     K              = 2*L+1
430     DTAUI(L,nw,ng) = DTAUKI(K,NW,NG) ! + 1.e-50
431
432     atemp = 0.
433     if(DTAUI(L,NW,NG) .GT. 1.0D-9) then
434        do iaer=1,naerkind
435           atemp = atemp + GIAER(K,NW,IAER)   * TAUAEROLK(K,NW,IAER)
436        end do
437        WBARI(L,nw,ng) = btemp(L,nw)  / DTAUI(L,NW,NG)
438     else
439        WBARI(L,nw,ng) = 0.0D0
440        DTAUI(L,NW,NG) = 1.0D-9
441     endif
442
443     if(btemp(L,nw) .GT. 0.0d0) then
444        cosbi(L,NW,NG) = atemp/btemp(L,nw)
445     else
446        cosbi(L,NW,NG) = 0.0D0
447     end if
448     
449
450     ! Now the other Gauss points, if needed.
451
452     DO NG=1,L_NGAUSS-1
453        IF(TAUGSURF(NW,NG) .gt. TLIMIT) THEN
454
455           DO L=1,L_NLAYRAD-1
456              K              = 2*L+1
457              DTAUI(L,nw,ng) = DTAUKI(K,NW,NG)+DTAUKI(K+1,NW,NG)! + 1.e-50
458
459              if(DTAUI(L,NW,NG) .GT. 1.0D-9) then
460
461                 WBARI(L,nw,ng) = btemp(L,nw)  / DTAUI(L,NW,NG)
462
463              else
464                 WBARI(L,nw,ng) = 0.0D0
465                 DTAUI(L,NW,NG) = 1.0D-9
466              endif
467
468              cosbi(L,NW,NG) = cosbi(L,NW,L_NGAUSS)
469           END DO ! L vertical loop
470           
471           ! Last level
472           L              = L_NLAYRAD
473           K              = 2*L+1
474           DTAUI(L,nw,ng) = DTAUKI(K,NW,NG)! + 1.e-50
475
476           if(DTAUI(L,NW,NG) .GT. 1.0D-9) then
477
478              WBARI(L,nw,ng) = btemp(L,nw)  / DTAUI(L,NW,NG)
479
480           else
481              WBARI(L,nw,ng) = 0.0D0
482              DTAUI(L,NW,NG) = 1.0D-9
483           endif
484
485           cosbi(L,NW,NG) = cosbi(L,NW,L_NGAUSS)
486           
487        END IF
488
489     END DO                 ! NG Gauss loop
490  END DO                    ! NW spectral loop
491
492  ! Total extinction optical depths
493
494  DO NG=1,L_NGAUSS       ! full gauss loop
495     DO NW=1,L_NSPECTI       
496        TAUCUMI(1,NW,NG)=0.0D0
497        DO K=2,L_LEVELS
498           TAUCUMI(K,NW,NG)=TAUCUMI(K-1,NW,NG)+DTAUKI(K,NW,NG)
499        END DO
500     END DO                 ! end full gauss loop
501  END DO
502
503  ! be aware when comparing with textbook results
504  ! (e.g. Pierrehumbert p. 218) that
505  ! taucumi does not take the <cos theta>=0.5 factor into
506  ! account. It is the optical depth for a vertically
507  ! ascending ray with angle theta = 0.
508
509  !open(127,file='taucum.out')
510  !do nw=1,L_NSPECTI
511  !   write(127,*) taucumi(L_LEVELS,nw,L_NGAUSS)
512  !enddo
513  !close(127)
514 
515!  print*,'WBARI'
516!  print*,WBARI
517!  print*,'DTAUI'
518!  print*,DTAUI
519!  call abort
520
521end subroutine optci
522
523END MODULE optci_mod
524
Note: See TracBrowser for help on using the repository browser.