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

Last change on this file since 3581 was 3238, checked in by mturbet, 11 months ago

minor change in mineral aerosol name format

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