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

Last change on this file since 537 was 374, checked in by emillour, 13 years ago

Generic GCM: Upgrade: The location of the 'datagcm' directory can now be given in the callphys.def file ( datadir = /absolute/path/to/datagcm ). Changed "datafile.h" into a F90 module "datafile_mod.F90" and spread this change to all routines that used to use "datafile.h".
EM

File size: 14.2 KB
Line 
1      subroutine suaer_corrk
2
3      ! inputs
4      use radinc_h,    only: L_NSPECTI,L_NSPECTV,naerkind,nsizemax,iim,jjm
5      use radcommon_h, only: blamv,blami,lamrefir,lamrefvis
6      use datafile_mod, only: datadir
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,omegarefvis,omegarefir
12
13      implicit none
14
15!==================================================================
16!     Purpose
17!     -------
18!     Initialize all radiative aerosol properties
19!
20!     Notes
21!     -----
22!     Reads the optical properties -> Mean  -> Variable assignment
23!     (ASCII files)                  (see radcommon_h.F90)
24!     wvl(nwvl)                      longsun
25!     ep(nwvl)                       epav     QVISsQREF(L_NSPECTV)
26!     omeg(nwvl)                     omegav   omegavis(L_NSPECTV)
27!     gfactor(nwvl)                  gav      gvis(L_NSPECTV)
28!     
29!     Authors
30!     -------
31!     Richard Fournier (1996) Francois Forget (1996)
32!     Frederic Hourdin
33!     Jean-jacques morcrette *ECMWF*
34!     MODIF Francois Forget (2000)
35!     MODIF Franck Montmessin (add water ice)
36!     MODIF J.-B. Madeleine 2008W27
37!     - Optical properties read in ASCII files
38!     - Add varying radius of the particules
39!     - Add water ice clouds
40!     MODIF R. Wordsworth (2009)
41!     - generalisation to arbitrary spectral bands
42!
43!==================================================================
44
45#include "callkeys.h"
46
47!     Optical properties (read in external ASCII files)
48      INTEGER          :: nwvl  ! Number of wavelengths in
49                                ! the domain (VIS or IR)
50
51!      REAL             :: solsir ! visible to infrared ratio
52!                                ! (tau_0.67um/tau_9um)
53
54      REAL, DIMENSION(:),&
55      ALLOCATABLE, SAVE :: wvl  ! Wavelength axis
56      REAL, DIMENSION(:),&
57      ALLOCATABLE, SAVE :: radiusdyn ! Particle size axis
58
59      REAL, DIMENSION(:,:),&
60      ALLOCATABLE, SAVE :: ep,& ! Extinction coefficient Qext
61      omeg,&                    ! Single Scattering Albedo
62      gfactor                   ! Assymetry Factor
63
64!     Local variables:
65
66      INTEGER :: iaer           ! Aerosol index
67      INTEGER :: idomain        ! Domain index (1=VIS,2=IR)
68      INTEGER :: ilw            ! longwave index
69      INTEGER :: isw           ! shortwave index
70      INTEGER :: isize          ! Particle size index
71      INTEGER :: jfile          ! ASCII file scan index
72      INTEGER :: file_unit = 60
73      LOGICAL :: file_ok, endwhile
74      CHARACTER(LEN=132) :: scanline ! ASCII file scanning line
75
76!     I/O  of "aerave" (subroutine that spectrally averages
77!     the single scattering parameters)
78
79      REAL lamref                      ! reference wavelengths
80      REAL epref                       ! reference extinction ep
81
82!      REAL epav(L_NSPECTI)            ! Average ep (= <Qext>/Qext(lamref) if epref=1)
83!      REAL omegav(L_NSPECTI)          ! Average single scattering albedo
84!      REAL gav(L_NSPECTI)             ! Average assymetry parameter
85
86      REAL epavVI(L_NSPECTV)            ! Average ep (= <Qext>/Qext(lamref) if epref=1)
87      REAL omegavVI(L_NSPECTV)          ! Average single scattering albedo
88      REAL gavVI(L_NSPECTV)             ! Average assymetry parameter
89
90      REAL epavIR(L_NSPECTI)            ! Average ep (= <Qext>/Qext(lamref) if epref=1)
91      REAL omegavIR(L_NSPECTI)          ! Average single scattering albedo
92      REAL gavIR(L_NSPECTI)             ! Average assymetry parameter
93     
94      logical forgetit                  ! use francois' data?
95      integer iwvl
96
97!     Local saved variables:
98
99      CHARACTER(LEN=30), DIMENSION(naerkind,2), SAVE :: file_id
100!---- Please indicate the names of the optical property files below
101!     Please also choose the reference wavelengths of each aerosol
102
103      forgetit=.true.
104
105!     naerkind=1
106!     visible
107      if(forgetit)then
108         file_id(1,1) = 'optprop_co2_vis_ff.dat' ! Francois' values
109      else
110         file_id(1,1) = 'optprop_co2ice_vis_n50.dat'
111      endif
112      lamrefvis(1)=1.5E-6   ! 1.5um OMEGA/MEx ???
113
114!     infrared
115      if(forgetit)then
116         file_id(1,2) = 'optprop_co2_ir_ff.dat' ! Francois' values
117      else
118         file_id(1,2) = 'optprop_co2ice_ir_n50.dat'
119      endif
120      lamrefir(1)=12.1E-6   ! 825cm-1 TES/MGS ???
121         
122! NOTE: these lamref's are currently for dust, not CO2 ice.
123! JB thinks this shouldn't matter too much, but it needs to be
124! fixed / decided for the final version.
125
126      IF (naerkind .gt. 1) THEN
127!     naerkind=2
128!     visible
129         file_id(2,1) = 'optprop_icevis_n50.dat'
130         lamrefvis(2)=1.5E-6   ! 1.5um OMEGA/MEx
131!     infrared
132         file_id(2,2) = 'optprop_iceir_n50.dat'
133         lamrefir(2)=12.1E-6   ! 825cm-1 TES/MGS
134      ENDIF
135
136      IF (naerkind .eq. 3) THEN
137!     naerkind=3
138!     visible
139         file_id(naerkind,1) = 'optprop_dustvis_n50.dat'
140         !lamrefvis(3)=1.5E-6   ! 1.5um OMEGA/MEx
141         lamrefvis(naerkind)=0.67e-6
142!     infrared
143         file_id(naerkind,2) = 'optprop_dustir_n50.dat'
144         !lamrefir(3)=12.1E-6   ! 825cm-1 TES/MGS
145         lamrefir(naerkind)=9.3E-6
146
147      ENDIF
148
149      IF (naerkind .gt. 3) THEN
150         print*,'naerkind = ',naerkind
151         print*,'but we only have data for 3 types, exiting.'
152         call abort
153      ENDIF
154
155!------------------------------------------------------------------
156
157!     Initializations:
158
159      radiustab(:,:,:) = 0.0
160      gvis(:,:,:)      = 0.0
161      omegavis(:,:,:)  = 0.0
162      QVISsQREF(:,:,:) = 0.0
163      gir(:,:,:)       = 0.0
164      omegair(:,:,:)   = 0.0
165      QIRsQREF(:,:,:)  = 0.0
166
167      DO iaer = 1, naerkind     ! Loop on aerosol kind
168         DO idomain = 1, 2      ! Loop on radiation domain (VIS or IR)
169!==================================================================
170!     1. READ OPTICAL PROPERTIES
171!==================================================================
172
173!     1.1 Open the ASCII file
174
175            INQUIRE(FILE=TRIM(datadir)//&
176         '/'//TRIM(file_id(iaer,idomain)),&
177            EXIST=file_ok)
178            IF(.NOT.file_ok) THEN
179               write(*,*)'suaer_corrk: Problem opening ',&
180               TRIM(file_id(iaer,idomain))
181               write(*,*)'It should be in: ',TRIM(datadir)
182               write(*,*)'1) You can set this directory address ',&
183               'in callphys.def with:'
184               write(*,*)' datadir = /absolute/path/to/datagcm'
185               write(*,*)'2) If ',&
186              TRIM(file_id(iaer,idomain)),&
187               ' is a LMD reference datafile, it'
188               write(*,*)' can be obtained online at:'
189               write(*,*)' http://www.lmd.jussieu.fr/',&
190               '~forget/datagcm/datafile'
191               CALL ABORT
192            ENDIF
193            OPEN(UNIT=file_unit,&
194            FILE=TRIM(datadir)//&
195         '/'//TRIM(file_id(iaer,idomain)),&
196            FORM='formatted')
197
198!     1.2 Allocate the optical property table
199
200            jfile = 1
201            endwhile = .false.
202            DO WHILE (.NOT.endwhile)
203               READ(file_unit,*) scanline
204               IF ((scanline(1:1) .ne. '#').and.&
205               (scanline(1:1) .ne. ' ')) THEN
206               BACKSPACE(file_unit)
207               reading1_seq: SELECT CASE (jfile) ! ====================
208            CASE(1) reading1_seq ! nwvl ----------------------------
209               read(file_unit,*) nwvl
210               jfile = jfile+1
211            CASE(2) reading1_seq ! nsize ---------------------------
212               read(file_unit,*) nsize(iaer,idomain)
213               endwhile = .true.
214            CASE DEFAULT reading1_seq ! ----------------------------
215               WRITE(*,*) 'readoptprop: ',&
216               'Error while loading optical properties.'
217               CALL ABORT
218            END SELECT reading1_seq ! ==============================
219         ENDIF
220      ENDDO
221
222      ALLOCATE(wvl(nwvl))       ! wvl
223      ALLOCATE(radiusdyn(nsize(iaer,idomain))) ! radiusdyn
224      ALLOCATE(ep(nwvl,nsize(iaer,idomain))) ! ep
225      ALLOCATE(omeg(nwvl,nsize(iaer,idomain))) ! omeg
226      ALLOCATE(gfactor(nwvl,nsize(iaer,idomain))) ! g
227
228
229!     1.3 Read the data
230
231      jfile = 1
232      endwhile = .false.
233      DO WHILE (.NOT.endwhile)
234         READ(file_unit,*) scanline
235         IF ((scanline(1:1) .ne. '#').and.&
236         (scanline(1:1) .ne. ' ')) THEN
237         BACKSPACE(file_unit)
238         reading2_seq: SELECT CASE (jfile) ! ====================
239      CASE(1) reading2_seq      ! wvl -----------------------------
240         read(file_unit,*) wvl
241         jfile = jfile+1
242      CASE(2) reading2_seq      ! radiusdyn -----------------------
243         read(file_unit,*) radiusdyn
244         jfile = jfile+1
245      CASE(3) reading2_seq      ! ep ------------------------------
246         isize = 1
247         DO WHILE (isize .le. nsize(iaer,idomain))
248            READ(file_unit,*) scanline
249            IF ((scanline(1:1) .ne. '#').and.&
250            (scanline(1:1) .ne. ' ')) THEN
251            BACKSPACE(file_unit)
252            read(file_unit,*) ep(:,isize)
253            isize = isize + 1
254         ENDIF
255      ENDDO
256
257      jfile = jfile+1
258      CASE(4) reading2_seq      ! omeg ----------------------------
259         isize = 1
260         DO WHILE (isize .le. nsize(iaer,idomain))
261            READ(file_unit,*) scanline
262            IF ((scanline(1:1) .ne. '#').and.&
263            (scanline(1:1) .ne. ' ')) THEN
264            BACKSPACE(file_unit)
265            read(file_unit,*) omeg(:,isize)
266            isize = isize + 1
267         ENDIF
268      ENDDO
269
270      jfile = jfile+1
271      CASE(5) reading2_seq      ! gfactor -------------------------
272         isize = 1
273         DO WHILE (isize .le. nsize(iaer,idomain))
274            READ(file_unit,*) scanline
275            IF ((scanline(1:1) .ne. '#').and.&
276            (scanline(1:1) .ne. ' ')) THEN
277            BACKSPACE(file_unit)
278            read(file_unit,*) gfactor(:,isize)
279            isize = isize + 1
280         ENDIF
281      ENDDO
282
283      jfile = jfile+1
284      IF ((idomain.NE.1).OR.(iaer.NE.1)) THEN
285         endwhile = .true.
286      ENDIF
287      CASE(6) reading2_seq
288         endwhile = .true.
289      CASE DEFAULT reading2_seq ! ----------------------------
290         WRITE(*,*) 'readoptprop: ',&
291         'Error while loading optical properties.'
292         CALL ABORT
293      END SELECT reading2_seq   ! ==============================
294      ENDIF
295      ENDDO
296
297!     1.4 Close the file
298
299      CLOSE(file_unit)
300
301!     1.5 If Francois' data, convert wvl to metres
302      if(forgetit)then
303         !   print*,'please sort out forgetit for naerkind>1'
304         if(iaer.eq.1)then
305            do iwvl=1,nwvl
306               wvl(iwvl)=wvl(iwvl)*1.e-6
307            enddo
308         endif
309      endif
310
311
312
313
314
315
316!==================================================================
317!     2. AVERAGED PROPERTIES AND VARIABLE ASSIGNMENTS
318!==================================================================
319      domain: SELECT CASE (idomain)
320!==================================================================
321      CASE(1) domain            !       VISIBLE DOMAIN (idomain=1)
322!==================================================================
323
324         lamref=lamrefvis(iaer)
325         epref=1.E+0
326
327         DO isize=1,nsize(iaer,idomain)
328
329!     Save the particle sizes
330            radiustab(iaer,idomain,isize)=radiusdyn(isize)
331
332!     Averaged optical properties (GCM channels)
333
334!            CALL aerave ( nwvl,&
335!            wvl(:),ep(:,isize),omeg(:,isize),gfactor(:,isize),&
336!            lamref,epref,tstellar,&
337!            L_NSPECTV,blamv,epavVI,&
338!            omegavVI,gavVI,QREFvis(iaer,isize))!,omegaREFir(iaer,isize))
339            CALL aerave_new ( nwvl,&
340            wvl(:),ep(:,isize),omeg(:,isize),gfactor(:,isize),&
341            lamref,epref,tstellar,&
342            L_NSPECTV,blamv,epavVI,&
343            omegavVI,gavVI,QREFvis(iaer,isize),omegaREFir(iaer,isize))
344
345!     Variable assignments (declared in radcommon)
346            DO isw=1,L_NSPECTV
347               QVISsQREF(isw,iaer,isize)=epavVI(isw)
348               gvis(isw,iaer,isize)=gavVI(isw)
349               omegavis(isw,iaer,isize)=omegavVI(isw)
350            END DO
351
352         ENDDO
353
354!==================================================================
355      CASE(2) domain            !      INFRARED DOMAIN (idomain=2)
356!==================================================================
357
358
359         DO isize=1,nsize(iaer,idomain) ! ----------------------------------
360
361            lamref=lamrefir(iaer)
362            epref=1.E+0
363
364!     Save the particle sizes
365            radiustab(iaer,idomain,isize)=radiusdyn(isize)
366
367!     Averaged optical properties (GCM channels)
368
369!     epav is <QIR>/Qext(lamrefir) since epref=1
370!     Note: aerave also computes the extinction coefficient at
371!     the reference wavelength. This is called QREFvis or QREFir
372!     (not epref, which is a different parameter).
373!     Reference wavelengths SHOULD BE defined for each aerosol in
374!     radcommon_h.F90
375
376!            CALL aerave ( nwvl,&
377!            wvl(:),ep(:,isize),omeg(:,isize),gfactor(:,isize),&
378!            lamref,epref,tplanet,&
379!            L_NSPECTI,blami,epavIR,&
380!            omegavIR,gavIR,QREFir(iaer,isize))!,omegaREFir(iaer,isize))
381            CALL aerave_new ( nwvl,&
382            wvl(:),ep(:,isize),omeg(:,isize),gfactor(:,isize),&
383            lamref,epref,tplanet,&
384            L_NSPECTI,blami,epavIR,&
385            omegavIR,gavIR,QREFir(iaer,isize),omegaREFir(iaer,isize))
386
387
388!     Variable assignments (declared in radcommon)
389            DO ilw=1,L_NSPECTI
390               QIRsQREF(ilw,iaer,isize)=epavIR(ilw)
391               gir(ilw,iaer,isize)=gavIR(ilw)
392               omegair(ilw,iaer,isize)=omegavIR(ilw)
393            END DO
394
395
396      ENDDO                     ! isize (particle size) -------------------------------------
397
398      END SELECT domain
399
400!========================================================================
401!     3. Deallocate temporary variables that were read in the ASCII files
402!========================================================================
403
404      DEALLOCATE(wvl)           ! wvl
405      DEALLOCATE(radiusdyn)     ! radiusdyn
406      DEALLOCATE(ep)            ! ep
407      DEALLOCATE(omeg)          ! omeg
408      DEALLOCATE(gfactor)       ! g
409
410      END DO                    ! Loop on iaer
411      END DO                    ! Loop on idomain
412!========================================================================
413
414      RETURN
415
416
417
418    END subroutine suaer_corrk
419     
Note: See TracBrowser for help on using the repository browser.