source: trunk/LMDZ.GENERIC/libf/phystd/suaer_corrk.F90 @ 2871

Last change on this file since 2871 was 2831, checked in by emillour, 3 years ago

Generic PCM:
Add the possibility to include Venus-like aerosols (triggered by option
aerovenus=.true. in callphys.def); baseline is to use 5 distinct scatterers
but each may be turned on/off (via aerovenus1, aerovenus2, aerovenus2p,
aerovenus3, aerovenusUV flags which may be specified in callphys.def).
GG

File size: 21.6 KB
Line 
1      subroutine suaer_corrk
2
3      ! inputs
4      use radinc_h,    only: L_NSPECTI,L_NSPECTV,nsizemax,naerkind
5      use radcommon_h, only: blamv,blami,lamrefir,lamrefvis
6      use datafile_mod, only: datadir, aerdir
7
8      ! outputs
9      use radcommon_h, only: QVISsQREF,omegavis,gvis,QIRsQREF,omegair,gir
10      use radcommon_h, only: radiustab,nsize,tstellar
11      use radcommon_h, only: qrefvis,qrefir,omegarefir !,omegarefvis
12      use aerosol_mod
13      use callkeys_mod, only: tplanet, optprop_back2lay_vis, optprop_back2lay_ir, &
14                              optprop_aeronlay_vis, optprop_aeronlay_ir,          &
15                              aeronlay_lamref, nlayaero,aerogeneric
16      use tracer_h, only: noms
17      implicit none
18
19!==================================================================
20!     Purpose
21!     -------
22!     Initialize all radiative aerosol properties
23!
24!     Notes
25!     -----
26!     Reads the optical properties -> Mean  -> Variable assignment
27!     (ASCII files)                  (see radcommon_h.F90)
28!     wvl(nwvl)                      longsun
29!     ep(nwvl)                       epav     QVISsQREF(L_NSPECTV)
30!     omeg(nwvl)                     omegav   omegavis(L_NSPECTV)
31!     gfactor(nwvl)                  gav      gvis(L_NSPECTV)
32!     
33!     Authors
34!     -------
35!     Richard Fournier (1996) Francois Forget (1996)
36!     Frederic Hourdin
37!     Jean-jacques morcrette *ECMWF*
38!     MODIF Francois Forget (2000)
39!     MODIF Franck Montmessin (add water ice)
40!     MODIF J.-B. Madeleine 2008W27
41!     - Optical properties read in ASCII files
42!     - Add varying radius of the particules
43!     - Add water ice clouds
44!     MODIF R. Wordsworth (2009)
45!     - generalisation to arbitrary spectral bands
46!
47!==================================================================
48
49!     Optical properties (read in external ASCII files)
50      INTEGER,SAVE      :: nwvl  ! Number of wavelengths in
51                                ! the domain (VIS or IR), read by master
52
53!      REAL             :: solsir ! visible to infrared ratio
54!                                ! (tau_0.67um/tau_9um)
55
56      REAL, DIMENSION(:),&
57      ALLOCATABLE, SAVE :: wvl  ! Wavelength axis, read by master
58      REAL, DIMENSION(:),&
59      ALLOCATABLE, SAVE :: radiusdyn ! Particle size axis, read by master
60
61      REAL, DIMENSION(:,:),&
62      ALLOCATABLE, SAVE :: ep,& ! Extinction coefficient Qext, read by master
63      omeg,&                    ! Single Scattering Albedo, read by master
64      gfactor                   ! Assymetry Factor, read by master
65
66!     Local variables:
67
68      INTEGER :: iaer           ! Aerosol index
69      INTEGER :: idomain        ! Domain index (1=VIS,2=IR)
70      INTEGER :: ilw            ! longwave index
71      INTEGER :: isw            ! shortwave index
72      INTEGER :: isize          ! Particle size index
73      INTEGER :: jfile          ! ASCII file scan index
74      INTEGER :: file_unit = 60
75      LOGICAL :: file_ok, endwhile
76      CHARACTER(LEN=132) :: scanline ! ASCII file scanning line
77
78!     I/O  of "aerave" (subroutine that spectrally averages
79!     the single scattering parameters)
80
81      REAL lamref                      ! reference wavelengths
82      REAL epref                       ! reference extinction ep
83
84      REAL epavVI(L_NSPECTV)            ! Average ep (= <Qext>/Qext(lamrefvis) if epref=1)
85      REAL omegavVI(L_NSPECTV)          ! Average single scattering albedo
86      REAL gavVI(L_NSPECTV)             ! Average assymetry parameter
87
88      REAL epavIR(L_NSPECTI)            ! Average ep (= <Qext>/Qext(lamrefir) if epref=1)
89      REAL omegavIR(L_NSPECTI)          ! Average single scattering albedo
90      REAL gavIR(L_NSPECTI)             ! Average assymetry parameter
91     
92      logical forgetit                  ! use francois' data?
93      integer iwvl, ia
94
95!     Local saved variables:
96
97      CHARACTER(LEN=50), DIMENSION(naerkind,2), SAVE :: file_id
98!$OMP THREADPRIVATE(file_id)
99!---- Please indicate the names of the optical property files below
100!     Please also choose the reference wavelengths of each aerosol     
101
102!--------- README TO UNDERSTAND WHAT FOLLOWS (JVO 20) -------
103!     The lamref variables comes from the Martian model
104!     where the visible one is the one used for computing
105!     and the IR one is used to output scaled opacity to
106!     match instrumental data ... This is done (at least
107!     for now) in the generic, so lamrefir is dummy*!
108
109!     The important one is the VISIBLE one as it will be used
110!     to rescale the values in callcork.F90 assuming 'aerosol' is
111!     at this visible reference wavelenght.
112
113!     *Actually if you change lamrefir here there is a
114!     slight sensitvity in the outputs because of some
115!     machine precision differences that amplifys and lead
116!     up to 10-6 differences in the radiative balance...
117!     This could be good to clean the code but would require
118!     a lot of modifs and to take care !
119
120!--------------------------------------------------------------
121      if (noaero) then
122        print*, 'naerkind= 0'
123      endif
124      do iaer=1,naerkind
125       if (iaer.eq.iaero_co2) then
126        forgetit=.true.
127          if (.not.noaero) then
128              print*, 'naerkind= co2', iaer
129          end if
130!     visible
131        if(forgetit)then
132           file_id(iaer,1) = 'optprop_co2_vis_ff.dat' ! Francois' values
133        else
134           file_id(iaer,1) = 'optprop_co2ice_vis_n50.dat'
135        endif
136        lamrefvis(iaer)=1.5E-6   ! 1.5um OMEGA/MEx ???
137
138!     infrared
139        if(forgetit)then
140           file_id(iaer,2) = 'optprop_co2_ir_ff.dat' ! Francois' values
141        else
142           file_id(iaer,2) = 'optprop_co2ice_ir_n50.dat'
143        endif
144        lamrefir(iaer)=12.1E-6   ! Dummy in generic phys. (JVO 20)
145       endif ! CO2 aerosols 
146! NOTE: these lamref's are currently for dust, not CO2 ice.
147! JB thinks this shouldn't matter too much, but it needs to be
148! fixed / decided for the final version.
149
150       if (iaer.eq.iaero_h2o) then
151        print*, 'naerkind= h2o', iaer
152
153!     visible
154         file_id(iaer,1) = 'optprop_icevis_n50.dat'
155         lamrefvis(iaer)=1.5E-6   ! 1.5um OMEGA/MEx
156!     infrared
157         file_id(iaer,2) = 'optprop_iceir_n50.dat'
158         lamrefir(iaer)=12.1E-6   ! Dummy in generic phys. (JVO 20)
159       endif
160
161       if (iaer.eq.iaero_dust) then
162        print*, 'naerkind= dust', iaer
163
164!     visible
165         file_id(iaer,1) = 'optprop_dustvis_n50.dat'
166         !lamrefvis(iaer)=1.5E-6   ! 1.5um OMEGA/MEx
167         lamrefvis(iaer)=0.67e-6
168!     infrared
169         file_id(iaer,2) = 'optprop_dustir_n50.dat'
170         lamrefir(iaer)=9.3E-6     ! Dummy in generic phys. (JVO 20)
171       endif
172
173       if (iaer.eq.iaero_h2so4) then
174         print*, 'naerkind= h2so4', iaer
175
176!     visible
177         file_id(iaer,1) = 'optprop_h2so4vis_n20.dat'
178         lamrefvis(iaer)=1.5E-6   ! no idea, must find
179!     infrared
180         file_id(iaer,2) = 'optprop_h2so4ir_n20.dat'
181         lamrefir(iaer)=9.3E-6 ! ! Dummy in generic phys. (JVO 20)
182! added by LK
183       endif
184
185       if (iaer.eq.iaero_back2lay) then
186         print*, 'naerkind= back2lay', iaer
187         
188!     visible
189         file_id(iaer,1) = TRIM(optprop_back2lay_vis)
190         lamrefvis(iaer)=0.8E-6  ! This is the important one.
191!     infrared
192         file_id(iaer,2) = TRIM(optprop_back2lay_ir)
193         lamrefir(iaer)=6.E-6    ! This is dummy.
194! added by SG
195       endif
196     
197       if (iaer.eq.iaero_nh3) then
198         print*, 'naerkind= nh3', iaer
199
200!     visible
201         file_id(iaer,1) = 'optprop_nh3ice_vis.dat'
202         lamrefvis(iaer)=0.756E-6  !
203!     infrared
204         file_id(iaer,2) = 'optprop_nh3ice_ir.dat'
205         lamrefir(iaer)=6.E-6  ! dummy
206! added by SG
207       endif
208
209       do ia=1,nlayaero
210         if (iaer.eq.iaero_nlay(ia)) then
211           print*, 'naerkind= nlay', iaer
212           
213!       visible
214           file_id(iaer,1) = TRIM(optprop_aeronlay_vis(ia))
215           lamrefvis(iaer)=aeronlay_lamref(ia)
216!       infrared
217           file_id(iaer,2) = TRIM(optprop_aeronlay_ir(ia))
218           lamrefir(iaer)=6.E-6 ! Dummy value
219         endif
220       enddo
221! added by JVO
222     
223       if (iaer.eq.iaero_aurora) then
224         print*, 'naerkind= aurora', iaer
225
226!     visible
227         file_id(iaer,1) = 'optprop_aurora_vis.dat'
228         lamrefvis(iaer)=0.3E-6  !
229!     infrared
230         file_id(iaer,2) = 'optprop_aurora_ir.dat'
231         lamrefir(iaer)=6.E-6  ! dummy
232! added by SG
233       endif
234
235! VENUS CLOUDS
236
237       if (iaer.eq.iaero_venus1) then
238         print*, 'naerkind= venus1', iaer
239
240!     visible
241         file_id(iaer,1) = 'optprop_h2so4vis_n50.dat'
242         lamrefvis(iaer)=1.5E-6   ! no idea, must find
243!     infrared
244         file_id(iaer,2) = 'optprop_h2so4ir_n50.dat'
245         lamrefir(iaer)=9.3E-6 ! no idea, must find
246! added by SL
247       endif
248
249       if (iaer.eq.iaero_venus2) then
250         print*, 'naerkind= venus2', iaer
251
252!     visible
253         file_id(iaer,1) = 'optprop_h2so4vis_n50.dat'
254         lamrefvis(iaer)=1.5E-6   ! no idea, must find
255!     infrared
256         file_id(iaer,2) = 'optprop_h2so4ir_n50.dat'
257         lamrefir(iaer)=9.3E-6 ! no idea, must find
258! added by SL
259       endif
260
261       if (iaer.eq.iaero_venus2p) then
262         print*, 'naerkind= venus2p', iaer
263
264!     visible
265         file_id(iaer,1) = 'optprop_h2so4vis_n50.dat'
266         lamrefvis(iaer)=1.5E-6   ! no idea, must find
267!     infrared
268         file_id(iaer,2) = 'optprop_h2so4ir_n50.dat'
269         lamrefir(iaer)=9.3E-6 ! no idea, must find
270! added by SL
271       endif
272
273       if (iaer.eq.iaero_venus3) then
274         print*, 'naerkind= venus3', iaer
275
276!     visible
277         file_id(iaer,1) = 'optprop_h2so4vis_n50.dat'
278         lamrefvis(iaer)=1.5E-6   ! no idea, must find
279!     infrared
280         file_id(iaer,2) = 'optprop_h2so4ir_n50.dat'
281         lamrefir(iaer)=9.3E-6 ! no idea, must find
282! added by SL
283       endif
284
285       if (iaer.eq.iaero_venusUV) then
286         print*, 'naerkind= venusUV', iaer
287
288!     visible
289         file_id(iaer,1) = 'optprop_venusUVvis.dat'
290         lamrefvis(iaer)=3.5E-7   ! Haus et al. 2015
291!     infrared
292         file_id(iaer,2) = 'optprop_venusUVir.dat'
293         lamrefir(iaer)=9.3E-6 ! not used anyway
294! added by SL
295       endif
296
297! END VENUS CLOUDS
298       
299! the following was added by LT
300       do ia=1,aerogeneric ! Read Radiative Generic Condensable Species data
301         if (iaer .eq. iaero_generic(ia)) then
302            if (index(noms(i_rgcs_ice(ia)),'Fe') .ne. 0) then
303               print*,"Reading Fe file"
304               file_id(iaer,1)='Fe.ocst.txt'
305               file_id(iaer,2)='Fe.ocst.txt'
306               lamrefvis(iaer) = 1.0E-6 !random pick
307               lamrefir(iaer) = 1.0E-6 !dummy but random pick
308            else if (index(noms(i_rgcs_ice(ia)),'Mn') .ne. 0) then
309               print*,"Reading MnS file"
310               file_id(iaer,1)='MnS.ocst.txt'
311               file_id(iaer,2)='MnS.ocst.txt'
312               lamrefvis(iaer) = 1.0E-6 !random pick
313               lamrefir(iaer) = 1.0E-6 !dummy but random pick   
314            else if (index(noms(i_rgcs_ice(ia)),'Mg') .ne. 0) then 
315               print*,"Reading Mg2SiO4 file"
316               file_id(iaer,1)='Mg2SiO4.ocst.txt'
317               file_id(iaer,2)='Mg2SiO4.ocst.txt'
318               lamrefvis(iaer) = 1.0E-6 !random pick
319               lamrefir(iaer) = 1.0E-6 !dummy but random pick 
320            else if (index(noms(i_rgcs_ice(ia)),'Cr') .ne. 0) then
321               print*,"Reading Cr file"
322               file_id(iaer,1)='Cr.ocst.txt'
323               file_id(iaer,2)='Cr.ocst.txt'
324               lamrefvis(iaer) = 1.0E-6 !random pick
325               lamrefir(iaer) = 1.0E-6 !dummy but random pick
326            else
327! If you want to add another specie, copy,paste & adapt the elseif block up here to your new specie (LT 2022)
328               call abort_physic("suaer_corrk", "Unknown specie in radiative generic condensable species",1)
329            endif
330         endif
331       enddo ! ia=1,aerogeneric
332      enddo ! of do iaer=1,naerkind
333     
334!------------------------------------------------------------------
335
336!     Initializations:
337
338      radiustab(:,:,:) = 0.0
339      gvis(:,:,:)      = 0.0
340      omegavis(:,:,:)  = 0.0
341      QVISsQREF(:,:,:) = 0.0
342      gir(:,:,:)       = 0.0
343      omegair(:,:,:)   = 0.0
344      QIRsQREF(:,:,:)  = 0.0
345
346      DO iaer = 1, naerkind     ! Loop on aerosol kind
347         DO idomain = 1, 2      ! Loop on radiation domain (VIS or IR)
348!==================================================================
349!     1. READ OPTICAL PROPERTIES
350!==================================================================
351
352!     1.1 Open the ASCII file
353
354!$OMP MASTER
355
356            INQUIRE(FILE=TRIM(datadir)//'/'//TRIM(aerdir)//&
357                    '/'//TRIM(file_id(iaer,idomain)),&
358                    EXIST=file_ok)
359            IF (file_ok) THEN
360              OPEN(UNIT=file_unit,&
361                   FILE=TRIM(datadir)//'/'//TRIM(aerdir)//&
362                        '/'//TRIM(file_id(iaer,idomain)),&
363                   FORM='formatted')
364            ELSE
365             ! In ye old days these files were stored in datadir;
366             ! let's be retro-compatible
367              INQUIRE(FILE=TRIM(datadir)//&
368                      '/'//TRIM(file_id(iaer,idomain)),&
369                      EXIST=file_ok)
370              IF (file_ok) THEN
371                OPEN(UNIT=file_unit,&
372                     FILE=TRIM(datadir)//&
373                          '/'//TRIM(file_id(iaer,idomain)),&
374                     FORM='formatted')
375              ENDIF             
376            ENDIF
377            IF(.NOT.file_ok) THEN
378               write(*,*)'suaer_corrk: Problem opening ',&
379               TRIM(file_id(iaer,idomain))
380               write(*,*)'It should be in: ',TRIM(datadir)//'/'//TRIM(aerdir)
381               write(*,*)'1) You can set this directory address ',&
382               'in callphys.def with:'
383               write(*,*)' datadir = /absolute/path/to/datagcm'
384               write(*,*)'2) If ',&
385              TRIM(file_id(iaer,idomain)),&
386               ' is a LMD reference datafile, it'
387               write(*,*)' can be obtained online at:'
388               write(*,*)' http://www.lmd.jussieu.fr/',&
389               '~lmdz/planets/LMDZ.GENERIC/datagcm/'
390               CALL ABORT
391            ENDIF
392
393!     1.2 Allocate the optical property table
394
395            jfile = 1
396            endwhile = .false.
397            DO WHILE (.NOT.endwhile)
398               READ(file_unit,*) scanline
399               IF ((scanline(1:1) .ne. '#').and.&
400               (scanline(1:1) .ne. ' ')) THEN
401               BACKSPACE(file_unit)
402               reading1_seq: SELECT CASE (jfile) ! ====================
403            CASE(1) reading1_seq ! nwvl ----------------------------
404               read(file_unit,*) nwvl
405               jfile = jfile+1
406            CASE(2) reading1_seq ! nsize ---------------------------
407               read(file_unit,*) nsize(iaer,idomain)
408               endwhile = .true.
409            CASE DEFAULT reading1_seq ! ----------------------------
410               WRITE(*,*) 'readoptprop: ',&
411               'Error while loading optical properties.'
412               CALL ABORT
413            END SELECT reading1_seq ! ==============================
414         ENDIF
415      ENDDO
416
417      ALLOCATE(wvl(nwvl))       ! wvl
418      ALLOCATE(radiusdyn(nsize(iaer,idomain))) ! radiusdyn
419      ALLOCATE(ep(nwvl,nsize(iaer,idomain))) ! ep
420      ALLOCATE(omeg(nwvl,nsize(iaer,idomain))) ! omeg
421      ALLOCATE(gfactor(nwvl,nsize(iaer,idomain))) ! g
422
423
424!     1.3 Read the data
425
426      jfile = 1
427      endwhile = .false.
428      DO WHILE (.NOT.endwhile)
429         READ(file_unit,*) scanline
430         IF ((scanline(1:1) .ne. '#').and.&
431         (scanline(1:1) .ne. ' ')) THEN
432         BACKSPACE(file_unit)
433         reading2_seq: SELECT CASE (jfile) ! ====================
434      CASE(1) reading2_seq      ! wvl -----------------------------
435         read(file_unit,*) wvl
436         jfile = jfile+1
437      CASE(2) reading2_seq      ! radiusdyn -----------------------
438         read(file_unit,*) radiusdyn
439         jfile = jfile+1
440      CASE(3) reading2_seq      ! ep ------------------------------
441         isize = 1
442         DO WHILE (isize .le. nsize(iaer,idomain))
443            READ(file_unit,*) scanline
444            IF ((scanline(1:1) .ne. '#').and.&
445            (scanline(1:1) .ne. ' ')) THEN
446            BACKSPACE(file_unit)
447            read(file_unit,*) ep(:,isize)
448            isize = isize + 1
449         ENDIF
450      ENDDO
451
452      jfile = jfile+1
453      CASE(4) reading2_seq      ! omeg ----------------------------
454         isize = 1
455         DO WHILE (isize .le. nsize(iaer,idomain))
456            READ(file_unit,*) scanline
457            IF ((scanline(1:1) .ne. '#').and.&
458            (scanline(1:1) .ne. ' ')) THEN
459            BACKSPACE(file_unit)
460            read(file_unit,*) omeg(:,isize)
461            isize = isize + 1
462         ENDIF
463      ENDDO
464
465      jfile = jfile+1
466      CASE(5) reading2_seq      ! gfactor -------------------------
467         isize = 1
468         DO WHILE (isize .le. nsize(iaer,idomain))
469            READ(file_unit,*) scanline
470            IF ((scanline(1:1) .ne. '#').and.&
471            (scanline(1:1) .ne. ' ')) THEN
472            BACKSPACE(file_unit)
473            read(file_unit,*) gfactor(:,isize)
474            isize = isize + 1
475         ENDIF
476      ENDDO
477
478      jfile = jfile+1
479      IF ((idomain.NE.iaero_co2).OR.(iaer.NE.iaero_co2)) THEN
480         endwhile = .true.
481      ENDIF
482      CASE(6) reading2_seq
483         endwhile = .true.
484      CASE DEFAULT reading2_seq ! ----------------------------
485         WRITE(*,*) 'readoptprop: ',&
486         'Error while loading optical properties.'
487         CALL ABORT
488      END SELECT reading2_seq   ! ==============================
489      ENDIF
490      ENDDO
491
492!     1.4 Close the file
493
494      CLOSE(file_unit)
495
496!     1.5 If Francois' data, convert wvl to metres
497       if(iaer.eq.iaero_co2)then
498         if(forgetit)then
499         !   print*,'please sort out forgetit for naerkind>1'
500            do iwvl=1,nwvl
501               wvl(iwvl)=wvl(iwvl)*1.e-6
502            enddo
503         endif
504      endif
505
506!$OMP END MASTER
507!$OMP BARRIER
508
509
510
511
512!==================================================================
513!     2. AVERAGED PROPERTIES AND VARIABLE ASSIGNMENTS
514!==================================================================
515      domain: SELECT CASE (idomain)
516!==================================================================
517      CASE(1) domain            !       VISIBLE DOMAIN (idomain=1)
518!==================================================================
519
520         lamref=lamrefvis(iaer)
521         epref=1.E+0
522
523         DO isize=1,nsize(iaer,idomain)
524
525!     Save the particle sizes
526            radiustab(iaer,idomain,isize)=radiusdyn(isize)
527
528!     Averaged optical properties (GCM channels)
529
530!            CALL aerave ( nwvl,&
531!            wvl(:),ep(:,isize),omeg(:,isize),gfactor(:,isize),&
532!            lamref,epref,tstellar,&
533!            L_NSPECTV,blamv,epavVI,&
534!            omegavVI,gavVI,QREFvis(iaer,isize))!,omegaREFir(iaer,isize))
535            CALL aerave_new ( nwvl,&
536            wvl(:),ep(:,isize),omeg(:,isize),gfactor(:,isize),&
537            lamref,epref,tstellar,&
538            L_NSPECTV,blamv,epavVI,&
539            omegavVI,gavVI,QREFvis(iaer,isize),omegaREFir(iaer,isize))
540
541!     Variable assignments (declared in radcommon)
542            DO isw=1,L_NSPECTV
543               QVISsQREF(isw,iaer,isize)=epavVI(isw)
544               gvis(isw,iaer,isize)=gavVI(isw)
545               omegavis(isw,iaer,isize)=omegavVI(isw)
546            END DO
547
548         ENDDO
549!==================================================================
550      CASE(2) domain            !      INFRARED DOMAIN (idomain=2)
551!==================================================================
552
553         DO isize=1,nsize(iaer,idomain) ! ----------------------------------
554
555            lamref=lamrefir(iaer)
556            epref=1.E+0
557
558!     Save the particle sizes
559            radiustab(iaer,idomain,isize)=radiusdyn(isize)
560
561!     Averaged optical properties (GCM channels)
562
563!     epav is <QIR>/Qext(lamrefir) since epref=1
564!     Note: aerave also computes the extinction coefficient at
565!     the reference wavelength. This is called QREFvis or QREFir
566!     (not epref, which is a different parameter).
567!     Reference wavelengths SHOULD BE defined for each aerosol in
568!     radcommon_h.F90
569
570!            CALL aerave ( nwvl,&
571!            wvl(:),ep(:,isize),omeg(:,isize),gfactor(:,isize),&
572!            lamref,epref,tplanet,&
573!            L_NSPECTI,blami,epavIR,&
574!            omegavIR,gavIR,QREFir(iaer,isize))!,omegaREFir(iaer,isize))
575            CALL aerave_new ( nwvl,&
576            wvl(:),ep(:,isize),omeg(:,isize),gfactor(:,isize),&
577            lamref,epref,tplanet,&
578            L_NSPECTI,blami,epavIR,&
579            omegavIR,gavIR,QREFir(iaer,isize),omegaREFir(iaer,isize))
580
581
582!     Variable assignments (declared in radcommon)
583            DO ilw=1,L_NSPECTI
584               QIRsQREF(ilw,iaer,isize)=epavIR(ilw)
585               gir(ilw,iaer,isize)=gavIR(ilw)
586               omegair(ilw,iaer,isize)=omegavIR(ilw)
587            END DO
588
589
590      ENDDO                     ! isize (particle size) -------------------------------------
591
592      END SELECT domain
593
594!========================================================================
595!     3. Deallocate temporary variables that were read in the ASCII files
596!========================================================================
597
598!$OMP BARRIER
599!$OMP MASTER
600      IF (ALLOCATED(wvl)) DEALLOCATE(wvl)                 ! wvl
601      IF (ALLOCATED(radiusdyn)) DEALLOCATE(radiusdyn)     ! radiusdyn
602      IF (ALLOCATED(ep)) DEALLOCATE(ep)                   ! ep
603      IF (ALLOCATED(omeg)) DEALLOCATE(omeg)               ! omeg
604      IF (ALLOCATED(gfactor)) DEALLOCATE(gfactor)         ! g
605!$OMP END MASTER
606!$OMP BARRIER
607
608      END DO                    ! Loop on iaer
609      END DO                    ! Loop on idomain
610!========================================================================
611      RETURN
612
613
614
615    END subroutine suaer_corrk
616     
Note: See TracBrowser for help on using the repository browser.