source: trunk/LMDZ.VENUS/libf/phyvenus/aeroptproperties.F90 @ 3556

Last change on this file since 3556 was 2560, checked in by slebonnois, 3 years ago

SL: Implementation of SW computation based on generic model. Switch between this new SW module or old module that reads R. Haus tables implemented with a key (solarchoice)

File size: 24.2 KB
RevLine 
[2560]1      SUBROUTINE aeroptproperties(ngrid,nlayer,reffrad,nueffrad,   &
2                                  QVISsQREF3d,omegaVIS3d,gVIS3d,QREFvis3d)
3
4      use radinc_h,    only: L_NSPECTV,nsizemax,naerkind
5      use radcommon_h, only: QVISsQREF,omegavis,gvis
6      use radcommon_h, only: qrefvis
7      use radcommon_h, only: radiustab,nsize
8
9      implicit none
10
11!     =============================================================
12!     Aerosol Optical Properties
13!
14!     Description:
15!       Compute the scattering parameters in each grid
16!       box, depending on aerosol grain sizes. Log-normal size
17!       distribution and Gauss-Legendre integration are used.
18
19!     Parameters:
20!       Don't forget to set the value of varyingnueff below; If
21!       the effective variance of the distribution for the given
22!       aerosol is considered homogeneous in the atmosphere, please
23!       set varyingnueff(iaer) to .false. Resulting computational
24!       time will be much better.
25
26!     Authors: J.-B. Madeleine, F. Forget, F. Montmessin
27!     Slightly modified and converted to F90 by R. Wordsworth (2009)
28!     Varying nueff section removed by R. Wordsworth for simplicity
29!     ==============================================================
30
31!     Local variables
32!     ---------------
33
34
35
36!     =============================================================
37      LOGICAL, PARAMETER :: varyingnueff(naerkind) = .false.
38!     =============================================================
39
40!     Min. and max radius of the interpolation grid (in METERS)
41      REAL, PARAMETER :: refftabmin = 2e-8 !2e-8
42!      REAL, PARAMETER :: refftabmax = 35e-6
43      REAL, PARAMETER :: refftabmax = 1e-3
44!     Log of the min and max variance of the interpolation grid
45      REAL, PARAMETER :: nuefftabmin = -4.6
46      REAL, PARAMETER :: nuefftabmax = 0.
47!     Number of effective radius of the interpolation grid
48      INTEGER, PARAMETER :: refftabsize = 200
49!     Number of effective variances of the interpolation grid
50!      INTEGER, PARAMETER :: nuefftabsize = 100
51      INTEGER, PARAMETER :: nuefftabsize = 1
52!     Interpolation grid indices (reff,nueff)
53      INTEGER :: grid_i,grid_j
54!     Intermediate variable
55      REAL :: var_tmp,var3d_tmp(ngrid,nlayer)
56!     Bilinear interpolation factors
57      REAL :: kx,ky,k1,k2,k3,k4
58!     Size distribution parameters
59      REAL :: sizedistk1,sizedistk2
60!     Pi!
61      REAL,SAVE :: pi
62!$OMP THREADPRIVATE(pi)
63!     Variables used by the Gauss-Legendre integration:
64      INTEGER radius_id,gausind
65      REAL kint
66      REAL drad
67      INTEGER, PARAMETER :: ngau = 10
68      REAL weightgaus(ngau),radgaus(ngau)
69      SAVE weightgaus,radgaus
70!      DATA weightgaus/.2955242247,.2692667193,.2190863625,.1494513491,.0666713443/
71!      DATA radgaus/.1488743389,.4333953941,.6794095682,.8650633666,.9739065285/
72      DATA    radgaus/0.07652652113350,0.22778585114165, &
73                      0.37370608871528,0.51086700195146, &
74                      0.63605368072468,0.74633190646476, &
75                      0.83911697181213,0.91223442826796, &
76                      0.96397192726078,0.99312859919241/
77
78      DATA weightgaus/0.15275338723120,0.14917298659407, &
79                      0.14209610937519,0.13168863843930, &
80                      0.11819453196154,0.10193011980823, &
81                      0.08327674160932,0.06267204829828, &
82                      0.04060142982019,0.01761400714091/
83!$OMP THREADPRIVATE(radgaus,weightgaus)
84!     Indices
85      INTEGER :: i,j,k,l,m,iaer,idomain
86      INTEGER :: ig,lg,chg
87
88!     Local saved variables
89!     ---------------------
90
91!     Radius axis of the interpolation grid
92      REAL,SAVE :: refftab(refftabsize)
93!     Variance axis of the interpolation grid
94      REAL,SAVE :: nuefftab(nuefftabsize)
95!     Volume ratio of the grid
96      REAL,SAVE :: logvratgrid,vratgrid
97!     Grid used to remember which calculation is done
98      LOGICAL,SAVE :: checkgrid(refftabsize,nuefftabsize,naerkind,2) = .false.
99!$OMP THREADPRIVATE(refftab,nuefftab,logvratgrid,vratgrid,checkgrid)
100!     Optical properties of the grid (VISIBLE)
101      REAL,SAVE :: qsqrefVISgrid(refftabsize,nuefftabsize,L_NSPECTV,naerkind)
102      REAL,SAVE :: qextVISgrid(refftabsize,nuefftabsize,L_NSPECTV,naerkind)
103      REAL,SAVE :: qscatVISgrid(refftabsize,nuefftabsize,L_NSPECTV,naerkind)
104      REAL,SAVE :: omegVISgrid(refftabsize,nuefftabsize,L_NSPECTV,naerkind)
105      REAL,SAVE :: gVISgrid(refftabsize,nuefftabsize,L_NSPECTV,naerkind)
106!$OMP THREADPRIVATE(qsqrefVISgrid,qextVISgrid,qscatVISgrid,omegVISgrid,gVISgrid)
107!     Optical properties of the grid (REFERENCE WAVELENGTHS)
108      REAL,SAVE :: qrefVISgrid(refftabsize,nuefftabsize,naerkind)
109      REAL,SAVE :: qscatrefVISgrid(refftabsize,nuefftabsize,naerkind)
110      REAL,SAVE :: omegrefVISgrid(refftabsize,nuefftabsize,naerkind)
111!$OMP THREADPRIVATE(qrefVISgrid,qscatrefVISgrid,omegrefVISgrid)
112!     Firstcall
113      LOGICAL,SAVE :: firstcall = .true.
114!$OMP THREADPRIVATE(firstcall)
115!     Variables used by the Gauss-Legendre integration:
116      REAL,SAVE :: normd(refftabsize,nuefftabsize,naerkind,2)
117      REAL,SAVE :: dista(refftabsize,nuefftabsize,naerkind,2,ngau)
118      REAL,SAVE :: distb(refftabsize,nuefftabsize,naerkind,2,ngau)
119!$OMP THREADPRIVATE(normd,dista,distb)
120
121      REAL,SAVE :: radGAUSa(ngau,naerkind,2)
122      REAL,SAVE :: radGAUSb(ngau,naerkind,2)
123!$OMP THREADPRIVATE(radGAUSa,radGAUSb)
124
125      REAL,SAVE :: qsqrefVISa(L_NSPECTV,ngau,naerkind)
126      REAL,SAVE :: qrefVISa(ngau,naerkind)
127      REAL,SAVE :: qsqrefVISb(L_NSPECTV,ngau,naerkind)
128      REAL,SAVE :: qrefVISb(ngau,naerkind)
129      REAL,SAVE :: omegVISa(L_NSPECTV,ngau,naerkind)
130      REAL,SAVE :: omegrefVISa(ngau,naerkind)
131      REAL,SAVE :: omegVISb(L_NSPECTV,ngau,naerkind)
132      REAL,SAVE :: omegrefVISb(ngau,naerkind)
133      REAL,SAVE :: gVISa(L_NSPECTV,ngau,naerkind)
134      REAL,SAVE :: gVISb(L_NSPECTV,ngau,naerkind)
135!$OMP THREADPRIVATE(qsqrefVISa,qrefVISa,qsqrefVISb,qrefVISb,omegVISa, &
136        !$OMP omegrefVISa,omegVISb,omegrefVISb,gVISa,gVISb)
137
138      REAL :: radiusm
139      REAL :: radiusr
140
141!     Inputs
142!     ------
143
144      INTEGER :: ngrid,nlayer
145!     Aerosol effective radius used for radiative transfer (meter)
146      REAL,INTENT(IN) :: reffrad(ngrid,nlayer,naerkind)
147!     Aerosol effective variance used for radiative transfer (n.u.)
148      REAL,INTENT(IN) :: nueffrad(ngrid,nlayer,naerkind)
149
150!     Outputs
151!     -------
152
153      REAL,INTENT(OUT) :: QVISsQREF3d(ngrid,nlayer,L_NSPECTV,naerkind)
154      REAL,INTENT(OUT) :: omegaVIS3d(ngrid,nlayer,L_NSPECTV,naerkind)
155      REAL,INTENT(OUT) :: gVIS3d(ngrid,nlayer,L_NSPECTV,naerkind)
156
157      REAL,INTENT(OUT) :: QREFvis3d(ngrid,nlayer,naerkind)
158
159      DO iaer = 1, naerkind ! Loop on aerosol kind
160!        IF ( (nsize(iaer,1).EQ.1).AND.(nsize(iaer,2).EQ.1) ) THEN
161        IF (nsize(iaer,1).EQ.1) THEN
162!==================================================================
163!       If there is one single particle size, optical
164!         properties of the considered aerosol are homogeneous
165          DO lg = 1, nlayer
166            DO ig = 1, ngrid
167              DO chg = 1, L_NSPECTV
168                QVISsQREF3d(ig,lg,chg,iaer)=QVISsQREF(chg,iaer,1)
169                omegaVIS3d(ig,lg,chg,iaer)=omegaVIS(chg,iaer,1)
170                gVIS3d(ig,lg,chg,iaer)=gVIS(chg,iaer,1)
171              ENDDO
172              QREFvis3d(ig,lg,iaer)=QREFvis(iaer,1)
173            ENDDO
174          ENDDO
175
176
177          if (firstcall) then
178             print*,'Optical prop. of the aerosol are homogenous for:'
179             print*,'iaer = ',iaer
180          endif
181
182        ELSE ! Varying effective radius and variance
183!      DO idomain = 1, 2 ! Loop on visible or infrared channel
184       idomain=1
185!==================================================================
186!     1. Creating the effective radius and variance grid
187!     --------------------------------------------------
188      IF (firstcall) THEN
189
190!       1.1 Pi!
191        pi = 2. * asin(1.e0)
192
193!       1.2 Effective radius
194        refftab(1)    = refftabmin
195        refftab(refftabsize) = refftabmax
196
197        logvratgrid = log(refftabmax/refftabmin) / float(refftabsize-1)*3.
198        vratgrid = exp(logvratgrid)
199
200        do i = 2, refftabsize-1
201          refftab(i) = refftab(i-1)*vratgrid**(1./3.)
202        enddo
203
204!       1.3 Effective variance
205        if(nuefftabsize.eq.1)then ! addded by RDW
206           print*,'Warning: no variance range in aeroptproperties'
207           nuefftab(1)=0.2
208        else
209           do i = 0, nuefftabsize-1
210              nuefftab(i+1) = exp( nuefftabmin + i*(nuefftabmax-nuefftabmin)/(nuefftabsize-1) )
211           enddo
212        endif
213
214        firstcall = .false.
215      ENDIF
216
217!       1.4 Radius middle point and range for Gauss integration
218        radiusm=0.5*(radiustab(iaer,idomain,nsize(iaer,idomain)) + radiustab(iaer,idomain,1))
219        radiusr=0.5*(radiustab(iaer,idomain,nsize(iaer,idomain)) - radiustab(iaer,idomain,1))
220
221!       1.5 Interpolating data at the Gauss quadrature points:
222        DO gausind=1,ngau
223          drad=radiusr*radgaus(gausind)
224          radGAUSa(gausind,iaer,idomain)=radiusm-drad
225
226          radius_id=minloc(abs(radiustab(iaer,idomain,:) - (radiusm-drad)),DIM=1)
227          IF ((radiustab(iaer,idomain,radius_id) - (radiusm-drad)).GT.0) THEN
228            radius_id=radius_id-1
229          ENDIF
230          IF (radius_id.GE.nsize(iaer,idomain)) THEN
231            radius_id=nsize(iaer,idomain)-1
232            kint = 1.
233          ELSEIF (radius_id.LT.1) THEN
234            radius_id=1
235            kint = 0.
236          ELSE
237          kint = ( (radiusm-drad) -                             &
238                   radiustab(iaer,idomain,radius_id) ) /        &
239                 ( radiustab(iaer,idomain,radius_id+1) -        &
240                   radiustab(iaer,idomain,radius_id) )
241          ENDIF
242!          IF (idomain.EQ.1) THEN ! VISIBLE DOMAIN -----------------
243            DO m=1,L_NSPECTV
244               qsqrefVISa(m,gausind,iaer)=                      &
245                    (1-kint)*QVISsQREF(m,iaer,radius_id) +      &
246                    kint*QVISsQREF(m,iaer,radius_id+1)
247            omegVISa(m,gausind,iaer)=                           &
248                    (1-kint)*omegaVIS(m,iaer,radius_id) +       &
249                    kint*omegaVIS(m,iaer,radius_id+1)
250            gVISa(m,gausind,iaer)=                              &
251                    (1-kint)*gVIS(m,iaer,radius_id) +           &
252                    kint*gVIS(m,iaer,radius_id+1)
253            ENDDO
254            qrefVISa(gausind,iaer)=                             &
255                    (1-kint)*QREFvis(iaer,radius_id) +          &
256                    kint*QREFvis(iaer,radius_id+1)
257            omegrefVISa(gausind,iaer)= 0
258!            omegrefVISa(gausind,iaer)=                          &
259!                    (1-kint)*omegaREFvis(iaer,radius_id) +      &
260!                    kint*omegaREFvis(iaer,radius_id+1)
261!          ENDIF
262        ENDDO
263
264        DO gausind=1,ngau
265          drad=radiusr*radgaus(gausind)
266          radGAUSb(gausind,iaer,idomain)=radiusm+drad
267
268          radius_id=minloc(abs(radiustab(iaer,idomain,:) -      &
269                               (radiusm+drad)),DIM=1)
270          IF ((radiustab(iaer,idomain,radius_id) -              &
271               (radiusm+drad)).GT.0) THEN
272            radius_id=radius_id-1
273          ENDIF
274          IF (radius_id.GE.nsize(iaer,idomain)) THEN
275            radius_id=nsize(iaer,idomain)-1
276            kint = 1.
277          ELSEIF (radius_id.LT.1) THEN
278            radius_id=1
279            kint = 0.
280          ELSE
281            kint = ( (radiusm+drad) -                           &
282                     radiustab(iaer,idomain,radius_id) ) /      &
283                   ( radiustab(iaer,idomain,radius_id+1) -      &
284                     radiustab(iaer,idomain,radius_id) )
285          ENDIF
286!          IF (idomain.EQ.1) THEN ! VISIBLE DOMAIN -----------------
287            DO m=1,L_NSPECTV
288            qsqrefVISb(m,gausind,iaer)=                         &
289                    (1-kint)*QVISsQREF(m,iaer,radius_id) +      &
290                    kint*QVISsQREF(m,iaer,radius_id+1)   
291            omegVISb(m,gausind,iaer)=                           &
292                    (1-kint)*omegaVIS(m,iaer,radius_id) +       &
293                    kint*omegaVIS(m,iaer,radius_id+1)
294            gVISb(m,gausind,iaer)=                              &
295                    (1-kint)*gVIS(m,iaer,radius_id) +           &
296                    kint*gVIS(m,iaer,radius_id+1)
297            ENDDO
298            qrefVISb(gausind,iaer)=                             &
299                    (1-kint)*QREFvis(iaer,radius_id) +          &
300                    kint*QREFvis(iaer,radius_id+1)
301            omegrefVISb(gausind,iaer)= 0
302!          ENDIF
303        ENDDO
304
305!==================================================================
306! CONSTANT NUEFF FROM HERE
307!==================================================================
308
309!     2. Compute the scattering parameters using linear
310!       interpolation over grain sizes and constant nueff
311!     ---------------------------------------------------
312
313      DO lg = 1,nlayer
314        DO ig = 1, ngrid
315!         2.1 Effective radius index and kx calculation
316          var_tmp=reffrad(ig,lg,iaer)/refftabmin
317          var_tmp=log(var_tmp)*3.
318          var_tmp=var_tmp/logvratgrid+1.
319          grid_i=floor(var_tmp)
320          IF (grid_i.GE.refftabsize) THEN
321!           WRITE(*,*) 'Warning: particle size in grid box #'
322!           WRITE(*,*) ig,' is too large to be used by the '
323!           WRITE(*,*) 'radiative transfer; please extend the '
324!           WRITE(*,*) 'interpolation grid to larger grain sizes.'
325            grid_i=refftabsize-1
326            kx = 1.
327          ELSEIF (grid_i.LT.1) THEN
328!           WRITE(*,*) 'Warning: particle size in grid box #'
329!           WRITE(*,*) ig,' is too small to be used by the '
330!           WRITE(*,*) 'radiative transfer; please extend the '
331!           WRITE(*,*) 'interpolation grid to smaller grain sizes.'
332            grid_i=1
333            kx = 0.
334          ELSE
335            kx = ( reffrad(ig,lg,iaer)-refftab(grid_i) ) /            &
336                 ( refftab(grid_i+1)-refftab(grid_i) )
337          ENDIF
338!         2.3 Integration
339          DO j=grid_i,grid_i+1
340!             2.3.1 Check if the calculation has been done
341              IF (.NOT.checkgrid(j,1,iaer,idomain)) THEN
342!               2.3.2 Log-normal dist., r_g and sigma_g are defined
343!                 in [hansen_1974], "Light scattering in planetary
344!                 atmospheres", Space Science Reviews 16 527-610.
345!                 Here, sizedistk1=r_g and sizedistk2=sigma_g^2
346                sizedistk2 = log(1.+nueffrad(1,1,iaer))
347                sizedistk1 = exp(2.5*sizedistk2)
348                sizedistk1 = refftab(j) / sizedistk1
349
350                normd(j,1,iaer,idomain) = 1e-30
351                DO gausind=1,ngau
352                  drad=radiusr*radgaus(gausind)
353                  dista(j,1,iaer,idomain,gausind) = LOG((radiusm-drad)/sizedistk1)
354                  dista(j,1,iaer,idomain,gausind) =                   &
355                    EXP(-dista(j,1,iaer,idomain,gausind) *            &
356                    dista(j,1,iaer,idomain,gausind) *                 &
357                    0.5e0/sizedistk2)/(radiusm-drad)                 
358                  dista(j,1,iaer,idomain,gausind) =                   &
359                    dista(j,1,iaer,idomain,gausind) /                 &
360                    (sqrt(2e0*pi*sizedistk2))
361
362                  distb(j,1,iaer,idomain,gausind) = LOG((radiusm+drad)/sizedistk1)
363                  distb(j,1,iaer,idomain,gausind) =                   &
364                    EXP(-distb(j,1,iaer,idomain,gausind) *            &
365                    distb(j,1,iaer,idomain,gausind) *                 &
366                    0.5e0/sizedistk2)/(radiusm+drad)
367                  distb(j,1,iaer,idomain,gausind) =                   &
368                    distb(j,1,iaer,idomain,gausind) /                 &
369                    (sqrt(2e0*pi*sizedistk2))
370
371                  normd(j,1,iaer,idomain)=normd(j,1,iaer,idomain) +   &
372                    weightgaus(gausind) *                             &
373                    (                                                 &
374                    distb(j,1,iaer,idomain,gausind) * pi *            &
375                    radGAUSb(gausind,iaer,idomain) *                  &
376                    radGAUSb(gausind,iaer,idomain) +                  &
377                    dista(j,1,iaer,idomain,gausind) * pi *            &
378                    radGAUSa(gausind,iaer,idomain) *                  &
379                    radGAUSa(gausind,iaer,idomain)                    &
380                    )
381                ENDDO
382!                IF (idomain.EQ.1) THEN ! VISIBLE DOMAIN -----------
383!                 2.3.3.vis Initialization
384                  qsqrefVISgrid(j,1,:,iaer)=0.
385                  qextVISgrid(j,1,:,iaer)=0.
386                  qscatVISgrid(j,1,:,iaer)=0.
387                  omegVISgrid(j,1,:,iaer)=0.
388                  gVISgrid(j,1,:,iaer)=0.
389                  qrefVISgrid(j,1,iaer)=0.
390                  qscatrefVISgrid(j,1,iaer)=0.
391                  omegrefVISgrid(j,1,iaer)=0.
392
393                  DO gausind=1,ngau
394                    DO m=1,L_NSPECTV
395!                     Convolution:
396                      qextVISgrid(j,1,m,iaer) =              &
397                        qextVISgrid(j,1,m,iaer) +            &
398                        weightgaus(gausind) *                &
399                        (                                    &
400                        qsqrefVISb(m,gausind,iaer) *         &
401                        qrefVISb(gausind,iaer) *             &
402                        pi*radGAUSb(gausind,iaer,idomain) *  &
403                        radGAUSb(gausind,iaer,idomain) *     &
404                        distb(j,1,iaer,idomain,gausind) +    &
405                        qsqrefVISa(m,gausind,iaer) *         &
406                        qrefVISa(gausind,iaer) *             &
407                        pi*radGAUSa(gausind,iaer,idomain) *  &
408                        radGAUSa(gausind,iaer,idomain) *     &
409                        dista(j,1,iaer,idomain,gausind)      &
410                        )
411                      qscatVISgrid(j,1,m,iaer) =             &
412                        qscatVISgrid(j,1,m,iaer) +           &
413                        weightgaus(gausind) *                &
414                        (                                    &
415                        omegVISb(m,gausind,iaer) *           &
416                        qsqrefVISb(m,gausind,iaer) *         &
417                        qrefVISb(gausind,iaer) *             &
418                        pi*radGAUSb(gausind,iaer,idomain) *  &
419                        radGAUSb(gausind,iaer,idomain) *     &
420                        distb(j,1,iaer,idomain,gausind) +    &
421                        omegVISa(m,gausind,iaer) *           &
422                        qsqrefVISa(m,gausind,iaer) *         &
423                        qrefVISa(gausind,iaer) *             &
424                        pi*radGAUSa(gausind,iaer,idomain) *  &
425                        radGAUSa(gausind,iaer,idomain) *     &
426                        dista(j,1,iaer,idomain,gausind)      &
427                        )
428                      gVISgrid(j,1,m,iaer) =                 &
429                        gVISgrid(j,1,m,iaer) +               &
430                        weightgaus(gausind) *                &
431                        (                                    &
432                        omegVISb(m,gausind,iaer) *           &
433                        qsqrefVISb(m,gausind,iaer) *         &
434                        qrefVISb(gausind,iaer) *             &
435                        gVISb(m,gausind,iaer) *              &
436                        pi*radGAUSb(gausind,iaer,idomain) *  &
437                        radGAUSb(gausind,iaer,idomain) *     &
438                        distb(j,1,iaer,idomain,gausind) +    &
439                        omegVISa(m,gausind,iaer) *           &
440                        qsqrefVISa(m,gausind,iaer) *         &
441                        qrefVISa(gausind,iaer) *             &
442                        gVISa(m,gausind,iaer) *              &
443                        pi*radGAUSa(gausind,iaer,idomain) *  &
444                        radGAUSa(gausind,iaer,idomain) *     &
445                        dista(j,1,iaer,idomain,gausind)      &
446                        )
447                    ENDDO
448                    qrefVISgrid(j,1,iaer) =                  &
449                      qrefVISgrid(j,1,iaer) +                &
450                      weightgaus(gausind) *                  &
451                      (                                      &
452                      qrefVISb(gausind,iaer) *               &
453                      pi*radGAUSb(gausind,iaer,idomain) *    &
454                      radGAUSb(gausind,iaer,idomain) *       &
455                      distb(j,1,iaer,idomain,gausind) +      &
456                      qrefVISa(gausind,iaer) *               &
457                      pi*radGAUSa(gausind,iaer,idomain) *    &
458                      radGAUSa(gausind,iaer,idomain) *       &
459                      dista(j,1,iaer,idomain,gausind)        &
460                      )
461                    qscatrefVISgrid(j,1,iaer) =              &
462                      qscatrefVISgrid(j,1,iaer) +            &
463                      weightgaus(gausind) *                  &
464                      (                                      &
465                      omegrefVISb(gausind,iaer) *            &
466                      qrefVISb(gausind,iaer) *               &
467                      pi*radGAUSb(gausind,iaer,idomain) *    &
468                      radGAUSb(gausind,iaer,idomain) *       &
469                      distb(j,1,iaer,idomain,gausind) +      &
470                      omegrefVISa(gausind,iaer) *            &
471                      qrefVISa(gausind,iaer) *               &
472                      pi*radGAUSa(gausind,iaer,idomain) *    &
473                      radGAUSa(gausind,iaer,idomain) *       &
474                      dista(j,1,iaer,idomain,gausind)        &
475                      )
476                  ENDDO
477
478                  qrefVISgrid(j,1,iaer)=qrefVISgrid(j,1,iaer) /          &
479                                normd(j,1,iaer,idomain)       
480                  qscatrefVISgrid(j,1,iaer)=qscatrefVISgrid(j,1,iaer) /  &
481                                normd(j,1,iaer,idomain)
482                  omegrefVISgrid(j,1,iaer)=qscatrefVISgrid(j,1,iaer) /   &
483                               qrefVISgrid(j,1,iaer)
484                  DO m=1,L_NSPECTV
485                    qextVISgrid(j,1,m,iaer)=qextVISgrid(j,1,m,iaer) /    &
486                                normd(j,1,iaer,idomain)
487                    qscatVISgrid(j,1,m,iaer)=qscatVISgrid(j,1,m,iaer) /  &
488                                normd(j,1,iaer,idomain)
489                    gVISgrid(j,1,m,iaer)=gVISgrid(j,1,m,iaer) /          &
490                                qscatVISgrid(j,1,m,iaer) /               &
491                                normd(j,1,iaer,idomain)
492
493                    qsqrefVISgrid(j,1,m,iaer)=qextVISgrid(j,1,m,iaer) /  &
494                                qrefVISgrid(j,1,iaer)
495                    omegVISgrid(j,1,m,iaer)=qscatVISgrid(j,1,m,iaer) /   &
496                                qextVISgrid(j,1,m,iaer)
497                  ENDDO
498!                ENDIF                  ! --------------------------
499                checkgrid(j,1,iaer,idomain) = .true.
500              ENDIF !checkgrid
501          ENDDO !grid_i
502!         2.4 Linear interpolation
503          k1 = (1-kx)
504          k2 = kx
505!          IF (idomain.EQ.1) THEN ! VISIBLE ------------------------
506          DO m=1,L_NSPECTV
507             QVISsQREF3d(ig,lg,m,iaer) =                           &
508                        k1*qsqrefVISgrid(grid_i,1,m,iaer) +        &
509                        k2*qsqrefVISgrid(grid_i+1,1,m,iaer)
510            omegaVIS3d(ig,lg,m,iaer) =                             &
511                        k1*omegVISgrid(grid_i,1,m,iaer) +          &
512                        k2*omegVISgrid(grid_i+1,1,m,iaer)
513            gVIS3d(ig,lg,m,iaer) =                                 &
514                        k1*gVISgrid(grid_i,1,m,iaer) +             &
515                        k2*gVISgrid(grid_i+1,1,m,iaer)
516          ENDDO !L_NSPECTV
517          QREFvis3d(ig,lg,iaer) =                                  &
518                        k1*qrefVISgrid(grid_i,1,iaer) +            &
519                        k2*qrefVISgrid(grid_i+1,1,iaer)
520!          ENDIF                  ! --------------------------------
521        ENDDO !nlayer
522      ENDDO !ngrid
523
524!==================================================================
525
526
527
528!      ENDDO ! idomain
529
530      ENDIF ! nsize = 1
531
532      ENDDO ! iaer (loop on aerosol kind)
533
534      RETURN
535    END SUBROUTINE aeroptproperties
536
537
538
539     
Note: See TracBrowser for help on using the repository browser.