source: trunk/LMDZ.PLUTO/libf/phypluto/aeroptproperties.F90 @ 3210

Last change on this file since 3210 was 3196, checked in by afalco, 22 months ago

Pluto PCM:

Included N2 condensation from PLUTO.old & read haze aerosols.

AF

File size: 40.0 KB
Line 
1module aeroptproperties_mod
2
3implicit none
4
5contains
6
7      SUBROUTINE aeroptproperties(ngrid,nlayer,reffrad,nueffrad,   &
8                                  QVISsQREF3d,omegaVIS3d,gVIS3d,   &
9                                  QIRsQREF3d,omegaIR3d,gIR3d,      &
10                                  QREFvis3d,QREFir3d)!,            &
11!                                  omegaREFvis3d,omegaREFir3d)
12
13      use radinc_h,    only: L_NSPECTI,L_NSPECTV,nsizemax,naerkind
14      use radcommon_h, only: QVISsQREF,omegavis,gvis,QIRsQREF,omegair,gir
15      use radcommon_h, only: qrefvis,qrefir,omegarefir !,omegarefvis
16      use radcommon_h, only: radiustab,nsize
17
18      implicit none
19
20!     =============================================================
21!     Aerosol Optical Properties
22!
23!     Description:
24!       Compute the scattering parameters in each grid
25!       box, depending on aerosol grain sizes. Log-normal size
26!       distribution and Gauss-Legendre integration are used.
27
28!     Parameters:
29!       Don't forget to set the value of varyingnueff below; If
30!       the effective variance of the distribution for the given
31!       aerosol is considered homogeneous in the atmosphere, please
32!       set varyingnueff(iaer) to .false. Resulting computational
33!       time will be much better.
34
35!     Authors: J.-B. Madeleine, F. Forget, F. Montmessin
36!     Slightly modified and converted to F90 by R. Wordsworth (2009)
37!     Varying nueff section removed by R. Wordsworth for simplicity
38!     ==============================================================
39
40!     Local variables
41!     ---------------
42
43
44
45!     =============================================================
46!      LOGICAL, PARAMETER :: varyingnueff(naerkind) = .false. ! not used!
47!     =============================================================
48
49!     Min. and max radius of the interpolation grid (in METERS)
50      REAL, PARAMETER :: refftabmin = 2e-8 !2e-8
51!      REAL, PARAMETER :: refftabmax = 35e-6
52      REAL, PARAMETER :: refftabmax = 1e-3
53!     Log of the min and max variance of the interpolation grid
54      REAL, PARAMETER :: nuefftabmin = -4.6
55      REAL, PARAMETER :: nuefftabmax = 0.
56!     Number of effective radius of the interpolation grid
57      INTEGER, PARAMETER :: refftabsize = 200
58!     Number of effective variances of the interpolation grid
59!      INTEGER, PARAMETER :: nuefftabsize = 100
60      INTEGER, PARAMETER :: nuefftabsize = 1
61!     Interpolation grid indices (reff,nueff)
62      INTEGER :: grid_i,grid_j
63!     Intermediate variable
64      REAL :: var_tmp,var3d_tmp(ngrid,nlayer)
65!     Bilinear interpolation factors
66      REAL :: kx,ky,k1,k2,k3,k4
67!     Size distribution parameters
68      REAL :: sizedistk1,sizedistk2
69!     Pi!
70      REAL,SAVE :: pi
71!$OMP THREADPRIVATE(pi)
72!     Variables used by the Gauss-Legendre integration:
73      INTEGER radius_id,gausind
74      REAL kint
75      REAL drad
76      INTEGER, PARAMETER :: ngau = 10
77      REAL weightgaus(ngau),radgaus(ngau)
78      SAVE weightgaus,radgaus
79!      DATA weightgaus/.2955242247,.2692667193,.2190863625,.1494513491,.0666713443/
80!      DATA radgaus/.1488743389,.4333953941,.6794095682,.8650633666,.9739065285/
81      DATA    radgaus/0.07652652113350,0.22778585114165, &
82                      0.37370608871528,0.51086700195146, &
83                      0.63605368072468,0.74633190646476, &
84                      0.83911697181213,0.91223442826796, &
85                      0.96397192726078,0.99312859919241/
86
87      DATA weightgaus/0.15275338723120,0.14917298659407, &
88                      0.14209610937519,0.13168863843930, &
89                      0.11819453196154,0.10193011980823, &
90                      0.08327674160932,0.06267204829828, &
91                      0.04060142982019,0.01761400714091/
92!$OMP THREADPRIVATE(radgaus,weightgaus)
93!     Indices
94      INTEGER :: i,j,k,l,m,iaer,idomain
95      INTEGER :: ig,lg,chg
96
97!     Local saved variables
98!     ---------------------
99
100!     Radius axis of the interpolation grid
101      REAL,SAVE :: refftab(refftabsize)
102!     Variance axis of the interpolation grid
103      REAL,SAVE :: nuefftab(nuefftabsize)
104!     Volume ratio of the grid
105      REAL,SAVE :: logvratgrid,vratgrid
106!     Grid used to remember which calculation is done
107      LOGICAL,SAVE,ALLOCATABLE :: checkgrid(:,:,:,:)
108!$OMP THREADPRIVATE(refftab,nuefftab,logvratgrid,vratgrid,checkgrid)
109!     Optical properties of the grid (VISIBLE)
110      REAL,SAVE,ALLOCATABLE :: qsqrefVISgrid(:,:,:,:)
111      REAL,SAVE,ALLOCATABLE :: qextVISgrid(:,:,:,:)
112      REAL,SAVE,ALLOCATABLE :: qscatVISgrid(:,:,:,:)
113      REAL,SAVE,ALLOCATABLE :: omegVISgrid(:,:,:,:)
114      REAL,SAVE,ALLOCATABLE :: gVISgrid(:,:,:,:)
115!$OMP THREADPRIVATE(qsqrefVISgrid,qextVISgrid,qscatVISgrid,omegVISgrid,gVISgrid)
116!     Optical properties of the grid (INFRARED)
117      REAL,SAVE,ALLOCATABLE :: qsqrefIRgrid(:,:,:,:)
118      REAL,SAVE,ALLOCATABLE :: qextIRgrid(:,:,:,:)
119      REAL,SAVE,ALLOCATABLE :: qscatIRgrid(:,:,:,:)
120      REAL,SAVE,ALLOCATABLE :: omegIRgrid(:,:,:,:)
121      REAL,SAVE,ALLOCATABLE :: gIRgrid(:,:,:,:)
122!$OMP THREADPRIVATE(qsqrefIRgrid,qextIRgrid,qscatIRgrid,omegIRgrid,gIRgrid)
123!     Optical properties of the grid (REFERENCE WAVELENGTHS)
124      REAL,SAVE,ALLOCATABLE :: qrefVISgrid(:,:,:)
125      REAL,SAVE,ALLOCATABLE :: qscatrefVISgrid(:,:,:)
126      REAL,SAVE,ALLOCATABLE :: qrefIRgrid(:,:,:)
127      REAL,SAVE,ALLOCATABLE :: qscatrefIRgrid(:,:,:)
128      REAL,SAVE,ALLOCATABLE :: omegrefVISgrid(:,:,:)
129      REAL,SAVE,ALLOCATABLE :: omegrefIRgrid(:,:,:)
130!$OMP THREADPRIVATE(qrefVISgrid,qscatrefVISgrid,qrefIRgrid,qscatrefIRgrid,omegrefVISgrid,&
131!$OMP omegrefIRgrid)
132!     Firstcall
133      LOGICAL,SAVE :: firstcall = .true.
134      LOGICAL,SAVE :: first_allocate=.true.
135!$OMP THREADPRIVATE(firstcall,first_allocate)
136!     Variables used by the Gauss-Legendre integration:
137      REAL,SAVE,ALLOCATABLE :: normd(:,:,:,:)
138      REAL,SAVE,ALLOCATABLE :: dista(:,:,:,:,:)
139      REAL,SAVE,ALLOCATABLE :: distb(:,:,:,:,:)
140!$OMP THREADPRIVATE(normd,dista,distb)
141
142      REAL,SAVE,ALLOCATABLE :: radGAUSa(:,:,:)
143      REAL,SAVE,ALLOCATABLE :: radGAUSb(:,:,:)
144!$OMP THREADPRIVATE(radGAUSa,radGAUSb)
145
146      REAL,SAVE,ALLOCATABLE :: qsqrefVISa(:,:,:)
147      REAL,SAVE,ALLOCATABLE :: qrefVISa(:,:)
148      REAL,SAVE,ALLOCATABLE :: qsqrefVISb(:,:,:)
149      REAL,SAVE,ALLOCATABLE :: qrefVISb(:,:)
150!$OMP THREADPRIVATE(qsqrefVISa,qrefVISa,qsqrefVISb,qrefVISb)
151      REAL,SAVE,ALLOCATABLE :: omegVISa(:,:,:)
152      REAL,SAVE,ALLOCATABLE :: omegrefVISa(:,:)
153      REAL,SAVE,ALLOCATABLE :: omegVISb(:,:,:)
154      REAL,SAVE,ALLOCATABLE :: omegrefVISb(:,:)
155      REAL,SAVE,ALLOCATABLE :: gVISa(:,:,:)
156      REAL,SAVE,ALLOCATABLE :: gVISb(:,:,:)
157!$OMP THREADPRIVATE(omegVISa,omegrefVISa,omegVISb,omegrefVISb,gVISa,gVISb)
158
159      REAL,SAVE,ALLOCATABLE :: qsqrefIRa(:,:,:)
160      REAL,SAVE,ALLOCATABLE :: qrefIRa(:,:)
161      REAL,SAVE,ALLOCATABLE :: qsqrefIRb(:,:,:)
162      REAL,SAVE,ALLOCATABLE :: qrefIRb(:,:)
163!$OMP THREADPRIVATE(qsqrefIRa,qrefIRa,qsqrefIRb,qrefIRb)
164      REAL,SAVE,ALLOCATABLE :: omegIRa(:,:,:)
165      REAL,SAVE,ALLOCATABLE :: omegrefIRa(:,:)
166      REAL,SAVE,ALLOCATABLE :: omegIRb(:,:,:)
167      REAL,SAVE,ALLOCATABLE :: omegrefIRb(:,:)
168      REAL,SAVE,ALLOCATABLE :: gIRa(:,:,:)
169      REAL,SAVE,ALLOCATABLE :: gIRb(:,:,:)
170!$OMP THREADPRIVATE(omegIRa,omegrefIRa,omegIRb,omegrefIRb,gIRa,gIRb)
171
172      REAL :: radiusm
173      REAL :: radiusr
174
175!     Inputs
176!     ------
177
178      INTEGER :: ngrid,nlayer
179!     Aerosol effective radius used for radiative transfer (meter)
180      REAL,INTENT(IN) :: reffrad(ngrid,nlayer,naerkind)
181!     Aerosol effective variance used for radiative transfer (n.u.)
182      REAL,INTENT(IN) :: nueffrad(ngrid,nlayer,naerkind)
183
184!     Outputs
185!     -------
186
187      REAL,INTENT(OUT) :: QVISsQREF3d(ngrid,nlayer,L_NSPECTV,naerkind)
188      REAL,INTENT(OUT) :: omegaVIS3d(ngrid,nlayer,L_NSPECTV,naerkind)
189      REAL,INTENT(OUT) :: gVIS3d(ngrid,nlayer,L_NSPECTV,naerkind)
190
191      REAL,INTENT(OUT) :: QIRsQREF3d(ngrid,nlayer,L_NSPECTI,naerkind)
192      REAL,INTENT(OUT) :: omegaIR3d(ngrid,nlayer,L_NSPECTI,naerkind)
193      REAL,INTENT(OUT) :: gIR3d(ngrid,nlayer,L_NSPECTI,naerkind)
194
195      REAL,INTENT(OUT) :: QREFvis3d(ngrid,nlayer,naerkind)
196      REAL,INTENT(OUT) :: QREFir3d(ngrid,nlayer,naerkind)
197
198!      REAL :: omegaREFvis3d(ngrid,nlayer,naerkind)
199!      REAL :: omegaREFir3d(ngrid,nlayer,naerkind)
200
201      REAL :: minrad    ! minimal radius in table .dat radiustab (outside 0)
202      REAL :: maxrad
203
204!     0. Allocate local saved arrays at firstcall
205!     --------------------------------------------------
206      IF (first_allocate) THEN
207        ! Grid used to remember computations already done at previous calls
208        ALLOCATE(checkgrid(refftabsize,nuefftabsize,naerkind,2))
209        checkgrid(:,:,:,:)=.false.
210        ! Optical properties of the grid (VISIBLE)
211        ALLOCATE(qsqrefVISgrid(refftabsize,nuefftabsize,L_NSPECTV,naerkind))
212        ALLOCATE(qextVISgrid(refftabsize,nuefftabsize,L_NSPECTV,naerkind))
213        ALLOCATE(qscatVISgrid(refftabsize,nuefftabsize,L_NSPECTV,naerkind))
214        ALLOCATE(omegVISgrid(refftabsize,nuefftabsize,L_NSPECTV,naerkind))
215        ALLOCATE(gVISgrid(refftabsize,nuefftabsize,L_NSPECTV,naerkind))
216        ! Optical properties of the grid (INFRARED)
217        ALLOCATE(qsqrefIRgrid(refftabsize,nuefftabsize,L_NSPECTI,naerkind))
218        ALLOCATE(qextIRgrid(refftabsize,nuefftabsize,L_NSPECTI,naerkind))
219        ALLOCATE(qscatIRgrid(refftabsize,nuefftabsize,L_NSPECTI,naerkind))
220        ALLOCATE(omegIRgrid(refftabsize,nuefftabsize,L_NSPECTI,naerkind))
221        ALLOCATE(gIRgrid(refftabsize,nuefftabsize,L_NSPECTI,naerkind))
222        ! Optical properties of the grid (REFERENCE WAVELENGTHS)
223        ALLOCATE(qrefVISgrid(refftabsize,nuefftabsize,naerkind))
224        ALLOCATE(qscatrefVISgrid(refftabsize,nuefftabsize,naerkind))
225        ALLOCATE(qrefIRgrid(refftabsize,nuefftabsize,naerkind))
226        ALLOCATE(qscatrefIRgrid(refftabsize,nuefftabsize,naerkind))
227        ALLOCATE(omegrefVISgrid(refftabsize,nuefftabsize,naerkind))
228        ALLOCATE(omegrefIRgrid(refftabsize,nuefftabsize,naerkind))
229        ! Variables used by the Gauss-Legendre integration:
230        ALLOCATE(normd(refftabsize,nuefftabsize,naerkind,2))
231        ALLOCATE(dista(refftabsize,nuefftabsize,naerkind,2,ngau))
232        ALLOCATE(distb(refftabsize,nuefftabsize,naerkind,2,ngau))
233        ALLOCATE(radGAUSa(ngau,naerkind,2))
234        ALLOCATE(radGAUSb(ngau,naerkind,2))
235        !
236        ALLOCATE(qsqrefVISa(L_NSPECTV,ngau,naerkind))
237        ALLOCATE(qrefVISa(ngau,naerkind))
238        ALLOCATE(qsqrefVISb(L_NSPECTV,ngau,naerkind))
239        ALLOCATE(qrefVISb(ngau,naerkind))
240        ALLOCATE(omegVISa(L_NSPECTV,ngau,naerkind))
241        ALLOCATE(omegrefVISa(ngau,naerkind))
242        ALLOCATE(omegVISb(L_NSPECTV,ngau,naerkind))
243        ALLOCATE(omegrefVISb(ngau,naerkind))
244        ALLOCATE(gVISa(L_NSPECTV,ngau,naerkind))
245        ALLOCATE(gVISb(L_NSPECTV,ngau,naerkind))
246        !
247        ALLOCATE(qsqrefIRa(L_NSPECTI,ngau,naerkind))
248        ALLOCATE(qrefIRa(ngau,naerkind))
249        ALLOCATE(qsqrefIRb(L_NSPECTI,ngau,naerkind))
250        ALLOCATE(qrefIRb(ngau,naerkind))
251
252        ALLOCATE(omegIRa(L_NSPECTI,ngau,naerkind))
253        ALLOCATE(omegrefIRa(ngau,naerkind))
254        ALLOCATE(omegIRb(L_NSPECTI,ngau,naerkind))
255        ALLOCATE(omegrefIRb(ngau,naerkind))
256        ALLOCATE(gIRa(L_NSPECTI,ngau,naerkind))
257        ALLOCATE(gIRb(L_NSPECTI,ngau,naerkind))
258
259        first_allocate=.false.
260      ENDIF ! of IF (first_allocate)
261
262      DO iaer = 1, naerkind ! Loop on aerosol kind
263        IF ( (nsize(iaer,1).EQ.1).AND.(nsize(iaer,2).EQ.1) ) THEN
264!==================================================================
265!       If there is one single particle size, optical
266!         properties of the considered aerosol are homogeneous
267          DO lg = 1, nlayer
268            DO ig = 1, ngrid
269              DO chg = 1, L_NSPECTV
270                QVISsQREF3d(ig,lg,chg,iaer)=QVISsQREF(chg,iaer,1)
271                omegaVIS3d(ig,lg,chg,iaer)=omegaVIS(chg,iaer,1)
272                gVIS3d(ig,lg,chg,iaer)=gVIS(chg,iaer,1)
273              ENDDO
274              DO chg = 1, L_NSPECTI
275                QIRsQREF3d(ig,lg,chg,iaer)=QIRsQREF(chg,iaer,1)
276                omegaIR3d(ig,lg,chg,iaer)=omegaIR(chg,iaer,1)
277                gIR3d(ig,lg,chg,iaer)=gIR(chg,iaer,1)
278              ENDDO
279              QREFvis3d(ig,lg,iaer)=QREFvis(iaer,1)
280              QREFir3d(ig,lg,iaer)=QREFir(iaer,1)
281!              omegaREFvis3d(ig,lg,iaer)=omegaREFvis(iaer,1)
282!              omegaREFir3d(ig,lg,iaer)=omegaREFir(iaer,1)
283            ENDDO
284          ENDDO
285
286
287          if (firstcall) then
288             print*,'Optical prop. of the aerosol are homogenous for:'
289             print*,'iaer = ',iaer
290          endif
291
292        ELSE ! Varying effective radius and variance
293      DO idomain = 1, 2 ! Loop on visible or infrared channel
294!==================================================================
295!     1. Creating the effective radius and variance grid
296!     --------------------------------------------------
297      IF (firstcall) THEN
298
299!       1.1 Pi!
300        pi = 2. * asin(1.e0)
301
302!       1.2 Effective radius
303        refftab(1)    = refftabmin
304        refftab(refftabsize) = refftabmax
305
306        logvratgrid = log(refftabmax/refftabmin) / float(refftabsize-1)*3.
307        vratgrid = exp(logvratgrid)
308
309        do i = 2, refftabsize-1
310          refftab(i) = refftab(i-1)*vratgrid**(1./3.)
311        enddo
312
313!       1.3 Effective variance
314        if(nuefftabsize.eq.1)then ! addded by RDW
315           print*,'Warning: no variance range in aeroptproperties'
316           nuefftab(1)=0.2
317        else
318           do i = 0, nuefftabsize-1
319              nuefftab(i+1) = exp( nuefftabmin + i*(nuefftabmax-nuefftabmin)/(nuefftabsize-1) )
320           enddo
321        endif
322
323        ! check that radiustab (table.dat) cover the reffrad used for haze : iaer=1
324        minrad=min(MINVAL(radiustab(1,1,1:nsize(1,1))),MINVAL(radiustab(1,2,1:nsize(1,2))))
325        maxrad=min(MAXVAL(radiustab(1,1,1:nsize(1,1))),MAXVAL(radiustab(1,2,1:nsize(1,2))))
326        IF ((MINVAL(reffrad).LE.minrad).OR.(MAXVAL(reffrad).GE.maxrad)) then
327           WRITE(*,*) 'Warning: particle size in grid box #'
328           WRITE(*,*) ig,' is too large to be used by the '
329           WRITE(*,*) 'radiative transfer; please extend the '
330           WRITE(*,*) 'interpolation grid to larger grain sizes.'
331           WRITE(*,*) 'radiustab=',minrad,'-',maxrad
332           WRITE(*,*) 'reffrad=',MINVAL(reffrad),'-',MAXVAL(reffrad)
333           stop
334        ENDIF
335
336        firstcall = .false.
337      ENDIF ! of IF (firstcall)
338
339!       1.4 Radius middle point and range for Gauss integration
340        radiusm=0.5*(radiustab(iaer,idomain,nsize(iaer,idomain)) + radiustab(iaer,idomain,1))
341        radiusr=0.5*(radiustab(iaer,idomain,nsize(iaer,idomain)) - radiustab(iaer,idomain,1))
342
343!       1.5 Interpolating data at the Gauss quadrature points:
344        DO gausind=1,ngau
345          drad=radiusr*radgaus(gausind)
346          radGAUSa(gausind,iaer,idomain)=radiusm-drad
347
348          radius_id=minloc(abs(radiustab(iaer,idomain,:) - (radiusm-drad)),DIM=1)
349          IF ((radiustab(iaer,idomain,radius_id) - (radiusm-drad)).GT.0) THEN
350            radius_id=radius_id-1
351          ENDIF
352          IF (radius_id.GE.nsize(iaer,idomain)) THEN
353            radius_id=nsize(iaer,idomain)-1
354            kint = 1.
355          ELSEIF (radius_id.LT.1) THEN
356            radius_id=1
357            kint = 0.
358          ELSE
359          kint = ( (radiusm-drad) -                             &
360                   radiustab(iaer,idomain,radius_id) ) /        &
361                 ( radiustab(iaer,idomain,radius_id+1) -        &
362                   radiustab(iaer,idomain,radius_id) )
363          ENDIF
364          IF (idomain.EQ.1) THEN ! VISIBLE DOMAIN -----------------
365            DO m=1,L_NSPECTV
366               qsqrefVISa(m,gausind,iaer)=                      &
367                    (1-kint)*QVISsQREF(m,iaer,radius_id) +      &
368                    kint*QVISsQREF(m,iaer,radius_id+1)
369            omegVISa(m,gausind,iaer)=                           &
370                    (1-kint)*omegaVIS(m,iaer,radius_id) +       &
371                    kint*omegaVIS(m,iaer,radius_id+1)
372            gVISa(m,gausind,iaer)=                              &
373                    (1-kint)*gVIS(m,iaer,radius_id) +           &
374                    kint*gVIS(m,iaer,radius_id+1)
375            ENDDO
376            qrefVISa(gausind,iaer)=                             &
377                    (1-kint)*QREFvis(iaer,radius_id) +          &
378                    kint*QREFvis(iaer,radius_id+1)
379            omegrefVISa(gausind,iaer)= 0
380!            omegrefVISa(gausind,iaer)=                          &
381!                    (1-kint)*omegaREFvis(iaer,radius_id) +      &
382!                    kint*omegaREFvis(iaer,radius_id+1)
383          ELSE ! INFRARED DOMAIN ----------------------------------
384            DO m=1,L_NSPECTI
385            qsqrefIRa(m,gausind,iaer)=                          &
386                    (1-kint)*QIRsQREF(m,iaer,radius_id) +       &
387                    kint*QIRsQREF(m,iaer,radius_id+1)
388            omegIRa(m,gausind,iaer)=                            &
389                    (1-kint)*omegaIR(m,iaer,radius_id) +        &
390                    kint*omegaIR(m,iaer,radius_id+1)
391            gIRa(m,gausind,iaer)=                               &
392                    (1-kint)*gIR(m,iaer,radius_id) +            &
393                    kint*gIR(m,iaer,radius_id+1)
394            ENDDO
395            qrefIRa(gausind,iaer)=                              &
396                    (1-kint)*QREFir(iaer,radius_id) +           &
397                    kint*QREFir(iaer,radius_id+1)
398            omegrefIRa(gausind,iaer)=                           &
399                    (1-kint)*omegaREFir(iaer,radius_id) +       &
400                    kint*omegaREFir(iaer,radius_id+1)
401          ENDIF
402        ENDDO
403
404        DO gausind=1,ngau
405          drad=radiusr*radgaus(gausind)
406          radGAUSb(gausind,iaer,idomain)=radiusm+drad
407
408          radius_id=minloc(abs(radiustab(iaer,idomain,:) -      &
409                               (radiusm+drad)),DIM=1)
410          IF ((radiustab(iaer,idomain,radius_id) -              &
411               (radiusm+drad)).GT.0) THEN
412            radius_id=radius_id-1
413          ENDIF
414          IF (radius_id.GE.nsize(iaer,idomain)) THEN
415            radius_id=nsize(iaer,idomain)-1
416            kint = 1.
417          ELSEIF (radius_id.LT.1) THEN
418            radius_id=1
419            kint = 0.
420          ELSE
421            kint = ( (radiusm+drad) -                           &
422                     radiustab(iaer,idomain,radius_id) ) /      &
423                   ( radiustab(iaer,idomain,radius_id+1) -      &
424                     radiustab(iaer,idomain,radius_id) )
425          ENDIF
426          IF (idomain.EQ.1) THEN ! VISIBLE DOMAIN -----------------
427            DO m=1,L_NSPECTV
428            qsqrefVISb(m,gausind,iaer)=                         &
429                    (1-kint)*QVISsQREF(m,iaer,radius_id) +      &
430                    kint*QVISsQREF(m,iaer,radius_id+1)
431            omegVISb(m,gausind,iaer)=                           &
432                    (1-kint)*omegaVIS(m,iaer,radius_id) +       &
433                    kint*omegaVIS(m,iaer,radius_id+1)
434            gVISb(m,gausind,iaer)=                              &
435                    (1-kint)*gVIS(m,iaer,radius_id) +           &
436                    kint*gVIS(m,iaer,radius_id+1)
437            ENDDO
438            qrefVISb(gausind,iaer)=                             &
439                    (1-kint)*QREFvis(iaer,radius_id) +          &
440                    kint*QREFvis(iaer,radius_id+1)
441            omegrefVISb(gausind,iaer)= 0
442!            omegrefVISb(gausind,iaer)=                          &
443!                    (1-kint)*omegaREFvis(iaer,radius_id) +      &
444!                    kint*omegaREFvis(iaer,radius_id+1)
445          ELSE ! INFRARED DOMAIN ----------------------------------
446            DO m=1,L_NSPECTI
447            qsqrefIRb(m,gausind,iaer)=                          &
448                    (1-kint)*QIRsQREF(m,iaer,radius_id) +       &
449                    kint*QIRsQREF(m,iaer,radius_id+1)
450            omegIRb(m,gausind,iaer)=                            &
451                    (1-kint)*omegaIR(m,iaer,radius_id) +        &
452                    kint*omegaIR(m,iaer,radius_id+1)
453            gIRb(m,gausind,iaer)=                               &
454                    (1-kint)*gIR(m,iaer,radius_id) +            &
455                    kint*gIR(m,iaer,radius_id+1)
456            ENDDO
457            qrefIRb(gausind,iaer)=                              &
458                    (1-kint)*QREFir(iaer,radius_id) +           &
459                    kint*QREFir(iaer,radius_id+1)
460            omegrefIRb(gausind,iaer)=                           &
461                    (1-kint)*omegaREFir(iaer,radius_id) +       &
462                    kint*omegaREFir(iaer,radius_id+1)
463          ENDIF
464        ENDDO
465
466!==================================================================
467! CONSTANT NUEFF FROM HERE
468!==================================================================
469
470!     2. Compute the scattering parameters using linear
471!       interpolation over grain sizes and constant nueff
472!     ---------------------------------------------------
473
474      DO lg = 1,nlayer
475        DO ig = 1, ngrid
476!         2.1 Effective radius index and kx calculation
477          var_tmp=reffrad(ig,lg,iaer)/refftabmin
478          var_tmp=log(var_tmp)*3.
479          var_tmp=var_tmp/logvratgrid+1.
480          grid_i=floor(var_tmp)
481          IF (grid_i.GE.refftabsize) THEN
482!           WRITE(*,*) 'Warning: particle size in grid box #'
483!           WRITE(*,*) ig,' is too large to be used by the '
484!           WRITE(*,*) 'radiative transfer; please extend the '
485!           WRITE(*,*) 'interpolation grid to larger grain sizes.'
486            grid_i=refftabsize-1
487            kx = 1.
488          ELSEIF (grid_i.LT.1) THEN
489!           WRITE(*,*) 'Warning: particle size in grid box #'
490!           WRITE(*,*) ig,' is too small to be used by the '
491!           WRITE(*,*) 'radiative transfer; please extend the '
492!           WRITE(*,*) 'interpolation grid to smaller grain sizes.'
493            grid_i=1
494            kx = 0.
495          ELSE
496            kx = ( reffrad(ig,lg,iaer)-refftab(grid_i) ) /            &
497                 ( refftab(grid_i+1)-refftab(grid_i) )
498          ENDIF
499!         2.3 Integration
500          DO j=grid_i,grid_i+1
501!             2.3.1 Check if the calculation has been done
502              IF (.NOT.checkgrid(j,1,iaer,idomain)) THEN
503!               2.3.2 Log-normal dist., r_g and sigma_g are defined
504!                 in [hansen_1974], "Light scattering in planetary
505!                 atmospheres", Space Science Reviews 16 527-610.
506!                 Here, sizedistk1=r_g and sizedistk2=sigma_g^2
507                sizedistk2 = log(1.+nueffrad(1,1,iaer))
508                sizedistk1 = exp(2.5*sizedistk2)
509                sizedistk1 = refftab(j) / sizedistk1
510
511                normd(j,1,iaer,idomain) = 1e-30
512                DO gausind=1,ngau
513                  drad=radiusr*radgaus(gausind)
514                  dista(j,1,iaer,idomain,gausind) = LOG((radiusm-drad)/sizedistk1)
515                  dista(j,1,iaer,idomain,gausind) =                   &
516                    EXP(-dista(j,1,iaer,idomain,gausind) *            &
517                    dista(j,1,iaer,idomain,gausind) *                 &
518                    0.5e0/sizedistk2)/(radiusm-drad)
519                  dista(j,1,iaer,idomain,gausind) =                   &
520                    dista(j,1,iaer,idomain,gausind) /                 &
521                    (sqrt(2e0*pi*sizedistk2))
522
523                  distb(j,1,iaer,idomain,gausind) = LOG((radiusm+drad)/sizedistk1)
524                  distb(j,1,iaer,idomain,gausind) =                   &
525                    EXP(-distb(j,1,iaer,idomain,gausind) *            &
526                    distb(j,1,iaer,idomain,gausind) *                 &
527                    0.5e0/sizedistk2)/(radiusm+drad)
528                  distb(j,1,iaer,idomain,gausind) =                   &
529                    distb(j,1,iaer,idomain,gausind) /                 &
530                    (sqrt(2e0*pi*sizedistk2))
531
532                  normd(j,1,iaer,idomain)=normd(j,1,iaer,idomain) +   &
533                    weightgaus(gausind) *                             &
534                    (                                                 &
535                    distb(j,1,iaer,idomain,gausind) * pi *            &
536                    radGAUSb(gausind,iaer,idomain) *                  &
537                    radGAUSb(gausind,iaer,idomain) +                  &
538                    dista(j,1,iaer,idomain,gausind) * pi *            &
539                    radGAUSa(gausind,iaer,idomain) *                  &
540                    radGAUSa(gausind,iaer,idomain)                    &
541                    )
542                ENDDO
543                IF (idomain.EQ.1) THEN ! VISIBLE DOMAIN -----------
544!                 2.3.3.vis Initialization
545                  qsqrefVISgrid(j,1,:,iaer)=0.
546                  qextVISgrid(j,1,:,iaer)=0.
547                  qscatVISgrid(j,1,:,iaer)=0.
548                  omegVISgrid(j,1,:,iaer)=0.
549                  gVISgrid(j,1,:,iaer)=0.
550                  qrefVISgrid(j,1,iaer)=0.
551                  qscatrefVISgrid(j,1,iaer)=0.
552                  omegrefVISgrid(j,1,iaer)=0.
553
554                  DO gausind=1,ngau
555                    DO m=1,L_NSPECTV
556!                     Convolution:
557                      qextVISgrid(j,1,m,iaer) =              &
558                        qextVISgrid(j,1,m,iaer) +            &
559                        weightgaus(gausind) *                &
560                        (                                    &
561                        qsqrefVISb(m,gausind,iaer) *         &
562                        qrefVISb(gausind,iaer) *             &
563                        pi*radGAUSb(gausind,iaer,idomain) *  &
564                        radGAUSb(gausind,iaer,idomain) *     &
565                        distb(j,1,iaer,idomain,gausind) +    &
566                        qsqrefVISa(m,gausind,iaer) *         &
567                        qrefVISa(gausind,iaer) *             &
568                        pi*radGAUSa(gausind,iaer,idomain) *  &
569                        radGAUSa(gausind,iaer,idomain) *     &
570                        dista(j,1,iaer,idomain,gausind)      &
571                        )
572                      qscatVISgrid(j,1,m,iaer) =             &
573                        qscatVISgrid(j,1,m,iaer) +           &
574                        weightgaus(gausind) *                &
575                        (                                    &
576                        omegVISb(m,gausind,iaer) *           &
577                        qsqrefVISb(m,gausind,iaer) *         &
578                        qrefVISb(gausind,iaer) *             &
579                        pi*radGAUSb(gausind,iaer,idomain) *  &
580                        radGAUSb(gausind,iaer,idomain) *     &
581                        distb(j,1,iaer,idomain,gausind) +    &
582                        omegVISa(m,gausind,iaer) *           &
583                        qsqrefVISa(m,gausind,iaer) *         &
584                        qrefVISa(gausind,iaer) *             &
585                        pi*radGAUSa(gausind,iaer,idomain) *  &
586                        radGAUSa(gausind,iaer,idomain) *     &
587                        dista(j,1,iaer,idomain,gausind)      &
588                        )
589                      gVISgrid(j,1,m,iaer) =                 &
590                        gVISgrid(j,1,m,iaer) +               &
591                        weightgaus(gausind) *                &
592                        (                                    &
593                        omegVISb(m,gausind,iaer) *           &
594                        qsqrefVISb(m,gausind,iaer) *         &
595                        qrefVISb(gausind,iaer) *             &
596                        gVISb(m,gausind,iaer) *              &
597                        pi*radGAUSb(gausind,iaer,idomain) *  &
598                        radGAUSb(gausind,iaer,idomain) *     &
599                        distb(j,1,iaer,idomain,gausind) +    &
600                        omegVISa(m,gausind,iaer) *           &
601                        qsqrefVISa(m,gausind,iaer) *         &
602                        qrefVISa(gausind,iaer) *             &
603                        gVISa(m,gausind,iaer) *              &
604                        pi*radGAUSa(gausind,iaer,idomain) *  &
605                        radGAUSa(gausind,iaer,idomain) *     &
606                        dista(j,1,iaer,idomain,gausind)      &
607                        )
608                    ENDDO
609                    qrefVISgrid(j,1,iaer) =                  &
610                      qrefVISgrid(j,1,iaer) +                &
611                      weightgaus(gausind) *                  &
612                      (                                      &
613                      qrefVISb(gausind,iaer) *               &
614                      pi*radGAUSb(gausind,iaer,idomain) *    &
615                      radGAUSb(gausind,iaer,idomain) *       &
616                      distb(j,1,iaer,idomain,gausind) +      &
617                      qrefVISa(gausind,iaer) *               &
618                      pi*radGAUSa(gausind,iaer,idomain) *    &
619                      radGAUSa(gausind,iaer,idomain) *       &
620                      dista(j,1,iaer,idomain,gausind)        &
621                      )
622                    qscatrefVISgrid(j,1,iaer) =              &
623                      qscatrefVISgrid(j,1,iaer) +            &
624                      weightgaus(gausind) *                  &
625                      (                                      &
626                      omegrefVISb(gausind,iaer) *            &
627                      qrefVISb(gausind,iaer) *               &
628                      pi*radGAUSb(gausind,iaer,idomain) *    &
629                      radGAUSb(gausind,iaer,idomain) *       &
630                      distb(j,1,iaer,idomain,gausind) +      &
631                      omegrefVISa(gausind,iaer) *            &
632                      qrefVISa(gausind,iaer) *               &
633                      pi*radGAUSa(gausind,iaer,idomain) *    &
634                      radGAUSa(gausind,iaer,idomain) *       &
635                      dista(j,1,iaer,idomain,gausind)        &
636                      )
637                  ENDDO
638
639                  qrefVISgrid(j,1,iaer)=qrefVISgrid(j,1,iaer) /          &
640                                normd(j,1,iaer,idomain)
641                  qscatrefVISgrid(j,1,iaer)=qscatrefVISgrid(j,1,iaer) /  &
642                                normd(j,1,iaer,idomain)
643                  omegrefVISgrid(j,1,iaer)=qscatrefVISgrid(j,1,iaer) /   &
644                               qrefVISgrid(j,1,iaer)
645                  DO m=1,L_NSPECTV
646                    qextVISgrid(j,1,m,iaer)=qextVISgrid(j,1,m,iaer) /    &
647                                normd(j,1,iaer,idomain)
648                    qscatVISgrid(j,1,m,iaer)=qscatVISgrid(j,1,m,iaer) /  &
649                                normd(j,1,iaer,idomain)
650                    gVISgrid(j,1,m,iaer)=gVISgrid(j,1,m,iaer) /          &
651                                qscatVISgrid(j,1,m,iaer) /               &
652                                normd(j,1,iaer,idomain)
653
654                    qsqrefVISgrid(j,1,m,iaer)=qextVISgrid(j,1,m,iaer) /  &
655                                qrefVISgrid(j,1,iaer)
656                    omegVISgrid(j,1,m,iaer)=qscatVISgrid(j,1,m,iaer) /   &
657                                qextVISgrid(j,1,m,iaer)
658                  ENDDO
659                ELSE                   ! INFRARED DOMAIN ----------
660!                 2.3.3.ir Initialization
661                  qsqrefIRgrid(j,1,:,iaer)=0.
662                  qextIRgrid(j,1,:,iaer)=0.
663                  qscatIRgrid(j,1,:,iaer)=0.
664                  omegIRgrid(j,1,:,iaer)=0.
665                  gIRgrid(j,1,:,iaer)=0.
666                  qrefIRgrid(j,1,iaer)=0.
667                  qscatrefIRgrid(j,1,iaer)=0.
668                  omegrefIRgrid(j,1,iaer)=0.
669
670                  DO gausind=1,ngau
671                    DO m=1,L_NSPECTI
672!                     Convolution:
673                      qextIRgrid(j,1,m,iaer) =                  &
674                        qextIRgrid(j,1,m,iaer) +                &
675                        weightgaus(gausind) *                   &
676                        (                                       &
677                        qsqrefIRb(m,gausind,iaer) *             &
678                        qrefVISb(gausind,iaer) *                &
679                        pi*radGAUSb(gausind,iaer,idomain) *     &
680                        radGAUSb(gausind,iaer,idomain) *        &
681                        distb(j,1,iaer,idomain,gausind) +       &
682                        qsqrefIRa(m,gausind,iaer) *             &
683                        qrefVISa(gausind,iaer) *                &
684                        pi*radGAUSa(gausind,iaer,idomain) *     &
685                        radGAUSa(gausind,iaer,idomain) *        &
686                        dista(j,1,iaer,idomain,gausind)         &
687                        )
688                      qscatIRgrid(j,1,m,iaer) =                 &
689                        qscatIRgrid(j,1,m,iaer) +               &
690                        weightgaus(gausind) *                   &
691                        (                                       &
692                        omegIRb(m,gausind,iaer) *               &
693                        qsqrefIRb(m,gausind,iaer) *             &
694                        qrefVISb(gausind,iaer) *                &
695                        pi*radGAUSb(gausind,iaer,idomain) *     &
696                        radGAUSb(gausind,iaer,idomain) *        &
697                        distb(j,1,iaer,idomain,gausind) +       &
698                        omegIRa(m,gausind,iaer) *               &
699                        qsqrefIRa(m,gausind,iaer) *             &
700                        qrefVISa(gausind,iaer) *                &
701                        pi*radGAUSa(gausind,iaer,idomain) *     &
702                        radGAUSa(gausind,iaer,idomain) *        &
703                        dista(j,1,iaer,idomain,gausind)         &
704                        )
705                      gIRgrid(j,1,m,iaer) =                     &
706                        gIRgrid(j,1,m,iaer) +                   &
707                        weightgaus(gausind) *                   &
708                        (                                       &
709                        omegIRb(m,gausind,iaer) *               &
710                        qsqrefIRb(m,gausind,iaer) *             &
711                        qrefVISb(gausind,iaer) *                &
712                        gIRb(m,gausind,iaer) *                  &
713                        pi*radGAUSb(gausind,iaer,idomain) *     &
714                        radGAUSb(gausind,iaer,idomain) *        &
715                        distb(j,1,iaer,idomain,gausind) +       &
716                        omegIRa(m,gausind,iaer) *               &
717                        qsqrefIRa(m,gausind,iaer) *             &
718                        qrefVISa(gausind,iaer) *                &
719                        gIRa(m,gausind,iaer) *                  &
720                        pi*radGAUSa(gausind,iaer,idomain) *     &
721                        radGAUSa(gausind,iaer,idomain) *        &
722                        dista(j,1,iaer,idomain,gausind)         &
723                        )
724                    ENDDO
725                    qrefIRgrid(j,1,iaer) =                      &
726                      qrefIRgrid(j,1,iaer) +                    &
727                      weightgaus(gausind) *                     &
728                      (                                         &
729                      qrefIRb(gausind,iaer) *                   &
730                      pi*radGAUSb(gausind,iaer,idomain) *       &
731                      radGAUSb(gausind,iaer,idomain) *          &
732                      distb(j,1,iaer,idomain,gausind) +         &
733                      qrefIRa(gausind,iaer) *                   &
734                      pi*radGAUSa(gausind,iaer,idomain) *       &
735                      radGAUSa(gausind,iaer,idomain) *          &
736                      dista(j,1,iaer,idomain,gausind)           &
737                      )
738                    qscatrefIRgrid(j,1,iaer) =                  &
739                      qscatrefIRgrid(j,1,iaer) +                &
740                      weightgaus(gausind) *                     &
741                      (                                         &
742                      omegrefIRb(gausind,iaer) *                &
743                      qrefIRb(gausind,iaer) *                   &
744                      pi*radGAUSb(gausind,iaer,idomain) *       &
745                      radGAUSb(gausind,iaer,idomain) *          &
746                      distb(j,1,iaer,idomain,gausind) +         &
747                      omegrefIRa(gausind,iaer) *                &
748                      qrefIRa(gausind,iaer) *                   &
749                      pi*radGAUSa(gausind,iaer,idomain) *       &
750                      radGAUSa(gausind,iaer,idomain) *          &
751                      dista(j,1,iaer,idomain,gausind)           &
752                      )
753                  ENDDO
754
755                  qrefIRgrid(j,1,iaer)=qrefIRgrid(j,1,iaer) /          &
756                                normd(j,1,iaer,idomain)
757                  qscatrefIRgrid(j,1,iaer)=qscatrefIRgrid(j,1,iaer) /  &
758                                normd(j,1,iaer,idomain)
759                  omegrefIRgrid(j,1,iaer)=qscatrefIRgrid(j,1,iaer) /   &
760                               qrefIRgrid(j,1,iaer)
761                  DO m=1,L_NSPECTI
762                    qextIRgrid(j,1,m,iaer)=qextIRgrid(j,1,m,iaer) /    &
763                                normd(j,1,iaer,idomain)
764                    qscatIRgrid(j,1,m,iaer)=qscatIRgrid(j,1,m,iaer) /  &
765                                normd(j,1,iaer,idomain)
766                    gIRgrid(j,1,m,iaer)=gIRgrid(j,1,m,iaer) /          &
767                                qscatIRgrid(j,1,m,iaer) /              &
768                                normd(j,1,iaer,idomain)
769
770                    qsqrefIRgrid(j,1,m,iaer)=qextIRgrid(j,1,m,iaer) /  &
771                                qrefVISgrid(j,1,iaer)
772                    omegIRgrid(j,1,m,iaer)=qscatIRgrid(j,1,m,iaer) /   &
773                                qextIRgrid(j,1,m,iaer)
774                  ENDDO
775                ENDIF                  ! --------------------------
776                checkgrid(j,1,iaer,idomain) = .true.
777              ENDIF !checkgrid
778          ENDDO !grid_i
779!         2.4 Linear interpolation
780          k1 = (1-kx)
781          k2 = kx
782          IF (idomain.EQ.1) THEN ! VISIBLE ------------------------
783          DO m=1,L_NSPECTV
784             QVISsQREF3d(ig,lg,m,iaer) =                           &
785                        k1*qsqrefVISgrid(grid_i,1,m,iaer) +        &
786                        k2*qsqrefVISgrid(grid_i+1,1,m,iaer)
787            omegaVIS3d(ig,lg,m,iaer) =                             &
788                        k1*omegVISgrid(grid_i,1,m,iaer) +          &
789                        k2*omegVISgrid(grid_i+1,1,m,iaer)
790            gVIS3d(ig,lg,m,iaer) =                                 &
791                        k1*gVISgrid(grid_i,1,m,iaer) +             &
792                        k2*gVISgrid(grid_i+1,1,m,iaer)
793          ENDDO !L_NSPECTV
794          QREFvis3d(ig,lg,iaer) =                                  &
795                        k1*qrefVISgrid(grid_i,1,iaer) +            &
796                        k2*qrefVISgrid(grid_i+1,1,iaer)
797!          omegaREFvis3d(ig,lg,iaer) =                              &
798!                        k1*omegrefVISgrid(grid_i,1,iaer) +         &
799!                        k2*omegrefVISgrid(grid_i+1,1,iaer)
800          ELSE                   ! INFRARED -----------------------
801          DO m=1,L_NSPECTI
802            QIRsQREF3d(ig,lg,m,iaer) =                             &
803                        k1*qsqrefIRgrid(grid_i,1,m,iaer) +         &
804                        k2*qsqrefIRgrid(grid_i+1,1,m,iaer)
805            omegaIR3d(ig,lg,m,iaer) =                              &
806                        k1*omegIRgrid(grid_i,1,m,iaer) +           &
807                        k2*omegIRgrid(grid_i+1,1,m,iaer)
808            gIR3d(ig,lg,m,iaer) =                                  &
809                        k1*gIRgrid(grid_i,1,m,iaer) +              &
810                        k2*gIRgrid(grid_i+1,1,m,iaer)
811          ENDDO !L_NSPECTI
812          QREFir3d(ig,lg,iaer) =                                   &
813                        k1*qrefIRgrid(grid_i,1,iaer) +             &
814                        k2*qrefIRgrid(grid_i+1,1,iaer)
815!          omegaREFir3d(ig,lg,iaer) =                               &
816!                        k1*omegrefIRgrid(grid_i,1,iaer) +          &
817!                        k2*omegrefIRgrid(grid_i+1,1,iaer)
818          ENDIF                  ! --------------------------------
819        ENDDO !nlayer
820      ENDDO !ngrid
821
822!==================================================================
823
824
825
826      ENDDO ! idomain
827
828      ENDIF ! nsize = 1
829
830      ENDDO ! iaer (loop on aerosol kind)
831
832    END SUBROUTINE aeroptproperties
833
834
835end module aeroptproperties_mod
Note: See TracBrowser for help on using the repository browser.