[3175] | 1 | subroutine suaer_corrk |
---|
| 2 | |
---|
| 3 | ! inputs |
---|
| 4 | use radinc_h, only: L_NSPECTI,L_NSPECTV,naerkind,nsizemax |
---|
| 5 | use radcommon_h, only: blamv,blami,lamrefir,lamrefvis |
---|
| 6 | use datafile_mod, only: datadir, aerdir, hazeprop |
---|
| 7 | use aerosol_mod |
---|
| 8 | |
---|
| 9 | ! outputs |
---|
| 10 | use radcommon_h, only: QVISsQREF,omegavis,gvis,QIRsQREF,omegair,gir |
---|
| 11 | use radcommon_h, only: radiustab,nsize,tstellar |
---|
| 12 | use radcommon_h, only: qrefvis,qrefir,omegarefvis,omegarefir |
---|
| 13 | |
---|
| 14 | implicit none |
---|
| 15 | |
---|
| 16 | !================================================================== |
---|
| 17 | ! Purpose |
---|
| 18 | ! ------- |
---|
| 19 | ! Initialize all radiative aerosol properties |
---|
| 20 | ! |
---|
| 21 | ! Notes |
---|
| 22 | ! ----- |
---|
| 23 | ! Reads the optical properties -> Mean -> Variable assignment |
---|
| 24 | ! (ASCII files) (see radcommon_h.F90) |
---|
| 25 | ! wvl(nwvl) longsun |
---|
| 26 | ! ep(nwvl) epav QVISsQREF(L_NSPECTV) |
---|
| 27 | ! omeg(nwvl) omegav omegavis(L_NSPECTV) |
---|
| 28 | ! gfactor(nwvl) gav gvis(L_NSPECTV) |
---|
| 29 | ! |
---|
| 30 | ! Authors |
---|
| 31 | ! ------- |
---|
| 32 | ! Richard Fournier (1996) Francois Forget (1996) |
---|
| 33 | ! Frederic Hourdin |
---|
| 34 | ! Jean-jacques morcrette *ECMWF* |
---|
| 35 | ! MODIF Francois Forget (2000) |
---|
| 36 | ! MODIF Franck Montmessin (add water ice) |
---|
| 37 | ! MODIF J.-B. Madeleine 2008W27 |
---|
| 38 | ! - Optical properties read in ASCII files |
---|
| 39 | ! - Add varying radius of the particules |
---|
| 40 | ! - Add water ice clouds |
---|
| 41 | ! MODIF R. Wordsworth (2009) |
---|
| 42 | ! - generalisation to arbitrary spectral bands |
---|
| 43 | ! MODIF Tanguy BERTRAND (2017) : adaptation Pluto |
---|
| 44 | !================================================================== |
---|
| 45 | |
---|
| 46 | #include "dimensions.h" |
---|
| 47 | #include "dimphys.h" |
---|
| 48 | #include "callkeys.h" |
---|
| 49 | #include "tracer.h" |
---|
| 50 | |
---|
| 51 | ! Optical properties (read in external ASCII files) |
---|
| 52 | INTEGER,SAVE :: nwvl ! Number of wavelengths in |
---|
| 53 | ! the domain (VIS or IR) |
---|
| 54 | |
---|
| 55 | REAL, DIMENSION(:),& |
---|
| 56 | ALLOCATABLE, SAVE :: wvl ! Wavelength axis |
---|
| 57 | REAL, DIMENSION(:),& |
---|
| 58 | ALLOCATABLE, SAVE :: radiusdyn ! Particle size axis |
---|
| 59 | |
---|
| 60 | REAL, DIMENSION(:,:),& |
---|
| 61 | ALLOCATABLE, SAVE :: ep,& ! Extinction coefficient Qext |
---|
| 62 | omeg,& ! Single Scattering Albedo |
---|
| 63 | gfactor ! Assymetry Factor |
---|
| 64 | |
---|
| 65 | ! Local variables: |
---|
| 66 | |
---|
| 67 | INTEGER :: iaer ! Aerosol index |
---|
| 68 | INTEGER :: idomain ! Domain index (1=VIS,2=IR) |
---|
| 69 | INTEGER :: ilw ! longwave index |
---|
| 70 | INTEGER :: isw ! shortwave index |
---|
| 71 | INTEGER :: isize ! Particle size index |
---|
| 72 | INTEGER :: jfile ! ASCII file scan index |
---|
| 73 | INTEGER :: file_unit = 60 |
---|
| 74 | LOGICAL :: file_ok, endwhile |
---|
| 75 | CHARACTER(LEN=132) :: scanline ! ASCII file scanning line |
---|
| 76 | |
---|
| 77 | ! I/O of "aerave" (subroutine that spectrally averages |
---|
| 78 | ! the single scattering parameters) |
---|
| 79 | |
---|
| 80 | REAL lamref ! reference wavelengths |
---|
| 81 | REAL epref ! reference extinction ep |
---|
| 82 | |
---|
| 83 | REAL epavVI(L_NSPECTV) ! Average ep (= <Qext>/Qext(lamref) if epref=1) |
---|
| 84 | REAL omegavVI(L_NSPECTV) ! Average single scattering albedo |
---|
| 85 | REAL gavVI(L_NSPECTV) ! Average assymetry parameter |
---|
| 86 | |
---|
| 87 | REAL epavIR(L_NSPECTI) ! Average ep (= <Qext>/Qext(lamref) if epref=1) |
---|
| 88 | REAL omegavIR(L_NSPECTI) ! Average single scattering albedo |
---|
| 89 | REAL gavIR(L_NSPECTI) ! Average assymetry parameter |
---|
| 90 | |
---|
| 91 | logical forgetit ! use francois' data? |
---|
| 92 | integer iwvl |
---|
| 93 | |
---|
| 94 | ! Local saved variables: |
---|
| 95 | |
---|
| 96 | CHARACTER(LEN=60), DIMENSION(naerkind,2), SAVE :: file_id |
---|
| 97 | ! --------------------------------------------------------------------------- |
---|
| 98 | |
---|
| 99 | !---- Please indicate the names of the optical property files below |
---|
| 100 | ! Please also choose the reference wavelengths of each aerosol |
---|
| 101 | |
---|
| 102 | ! NAERKIND dans radinc_h : = 1 |
---|
| 103 | ! Flag haze doit etre dans callcorrk... |
---|
| 104 | |
---|
| 105 | write(*,*) 'naerkind=',naerkind |
---|
| 106 | do iaer=1,naerkind |
---|
| 107 | ! naerkind=1: haze |
---|
| 108 | if (iaer.eq.1) then |
---|
| 109 | |
---|
| 110 | ! Only one table of optical properties : |
---|
| 111 | write(*,*)'Suaer haze optical properties, using: ', & |
---|
| 112 | TRIM(hazeprop) |
---|
| 113 | |
---|
| 114 | ! visible |
---|
| 115 | file_id(iaer,1)=TRIM(hazeprop) |
---|
| 116 | ! infrared |
---|
| 117 | file_id(iaer,2)=file_id(iaer,1) |
---|
| 118 | |
---|
| 119 | lamrefvis(iaer)=0.607E-6 ! reference wavelength for opacity vis pivot wavelength Cheng et al 2017 |
---|
| 120 | lamrefir(iaer)=2.E-6 ! reference wavelength for opacity IR (in the LEISA range) |
---|
| 121 | |
---|
| 122 | endif |
---|
| 123 | enddo |
---|
| 124 | |
---|
| 125 | IF (naerkind .gt. 1) THEN |
---|
| 126 | print*,'naerkind = ',naerkind |
---|
| 127 | print*,'but we only have data for 1 type, exiting.' |
---|
| 128 | call abort |
---|
| 129 | ENDIF |
---|
| 130 | |
---|
| 131 | !------------------------------------------------------------------ |
---|
| 132 | ! Initializations: |
---|
| 133 | |
---|
| 134 | radiustab(:,:,:) = 0.0 |
---|
| 135 | gvis(:,:,:) = 0.0 |
---|
| 136 | omegavis(:,:,:) = 0.0 |
---|
| 137 | QVISsQREF(:,:,:) = 0.0 |
---|
| 138 | gir(:,:,:) = 0.0 |
---|
| 139 | omegair(:,:,:) = 0.0 |
---|
| 140 | QIRsQREF(:,:,:) = 0.0 |
---|
| 141 | |
---|
| 142 | DO iaer = 1, naerkind ! Loop on aerosol kind |
---|
| 143 | |
---|
| 144 | DO idomain = 1, 2 ! Loop on radiation domain (VIS or IR) |
---|
| 145 | |
---|
| 146 | !================================================================== |
---|
| 147 | ! 1. READ OPTICAL PROPERTIES |
---|
| 148 | !================================================================== |
---|
| 149 | |
---|
| 150 | ! 1.1 Open the ASCII file |
---|
| 151 | INQUIRE(FILE=TRIM(datadir)//'/'//TRIM(aerdir)//& |
---|
| 152 | '/'//TRIM(file_id(iaer,idomain)),& |
---|
| 153 | EXIST=file_ok) |
---|
| 154 | write(*,*)' opening file=',TRIM(datadir)//'/'//TRIM(aerdir)//& |
---|
| 155 | '/'//TRIM(file_id(iaer,idomain)) |
---|
| 156 | IF (file_ok) THEN |
---|
| 157 | OPEN(UNIT=file_unit,& |
---|
| 158 | FILE=TRIM(datadir)//'/'//TRIM(aerdir)//& |
---|
| 159 | '/'//TRIM(file_id(iaer,idomain)),& |
---|
| 160 | FORM='formatted') |
---|
| 161 | ENDIF |
---|
| 162 | IF(.NOT.file_ok) THEN |
---|
| 163 | write(*,*)'suaer_corrk: Pb opening ',& |
---|
| 164 | TRIM(file_id(iaer,idomain)) |
---|
| 165 | write(*,*)'It should be in: ',& |
---|
| 166 | TRIM(datadir)//'/'//TRIM(aerdir) |
---|
| 167 | write(*,*)'1) You can set this directory address ',& |
---|
| 168 | 'in callphys.def with:' |
---|
| 169 | write(*,*)' datadir = /absolute/path/to/datagcm' |
---|
| 170 | CALL ABORT |
---|
| 171 | ENDIF |
---|
| 172 | |
---|
| 173 | ! 1.2 Allocate the optical property table |
---|
| 174 | |
---|
| 175 | jfile = 1 |
---|
| 176 | endwhile = .false. |
---|
| 177 | DO WHILE (.NOT.endwhile) |
---|
| 178 | READ(file_unit,*) scanline |
---|
| 179 | IF ((scanline(1:1) .ne. '#').and.& |
---|
| 180 | (scanline(1:1) .ne. ' ')) THEN |
---|
| 181 | BACKSPACE(file_unit) |
---|
| 182 | reading1_seq: SELECT CASE (jfile) ! ==================== |
---|
| 183 | CASE(1) reading1_seq ! nwvl ---------------------------- |
---|
| 184 | read(file_unit,*) nwvl |
---|
| 185 | jfile = jfile+1 |
---|
| 186 | CASE(2) reading1_seq ! nsize --------------------------- |
---|
| 187 | read(file_unit,*) nsize(iaer,idomain) |
---|
| 188 | endwhile = .true. |
---|
| 189 | CASE DEFAULT reading1_seq ! ---------------------------- |
---|
| 190 | WRITE(*,*) 'readoptprop: ',& |
---|
| 191 | 'Error while loading optical properties.' |
---|
| 192 | CALL ABORT |
---|
| 193 | END SELECT reading1_seq ! ============================== |
---|
| 194 | ENDIF |
---|
| 195 | ENDDO |
---|
| 196 | |
---|
| 197 | ALLOCATE(wvl(nwvl)) ! wvl |
---|
| 198 | ALLOCATE(radiusdyn(nsize(iaer,idomain))) ! radiusdyn |
---|
| 199 | ALLOCATE(ep(nwvl,nsize(iaer,idomain))) ! ep |
---|
| 200 | ALLOCATE(omeg(nwvl,nsize(iaer,idomain))) ! omeg |
---|
| 201 | ALLOCATE(gfactor(nwvl,nsize(iaer,idomain))) ! g |
---|
| 202 | |
---|
| 203 | ! 1.3 Read the data |
---|
| 204 | |
---|
| 205 | jfile = 1 |
---|
| 206 | endwhile = .false. |
---|
| 207 | DO WHILE (.NOT.endwhile) |
---|
| 208 | READ(file_unit,*) scanline |
---|
| 209 | IF ((scanline(1:1) .ne. '#').and.& |
---|
| 210 | (scanline(1:1) .ne. ' ')) THEN |
---|
| 211 | BACKSPACE(file_unit) |
---|
| 212 | reading2_seq: SELECT CASE (jfile) ! ==================== |
---|
| 213 | CASE(1) reading2_seq ! wvl ----------------------------- |
---|
| 214 | read(file_unit,*) wvl |
---|
| 215 | jfile = jfile+1 |
---|
| 216 | CASE(2) reading2_seq ! radiusdyn ----------------------- |
---|
| 217 | read(file_unit,*) radiusdyn |
---|
| 218 | jfile = jfile+1 |
---|
| 219 | CASE(3) reading2_seq ! ep ------------------------------ |
---|
| 220 | isize = 1 |
---|
| 221 | DO WHILE (isize .le. nsize(iaer,idomain)) |
---|
| 222 | READ(file_unit,*) scanline |
---|
| 223 | IF ((scanline(1:1) .ne. '#').and.& |
---|
| 224 | (scanline(1:1) .ne. ' ')) THEN |
---|
| 225 | BACKSPACE(file_unit) |
---|
| 226 | read(file_unit,*) ep(:,isize) |
---|
| 227 | isize = isize + 1 |
---|
| 228 | ENDIF |
---|
| 229 | ENDDO |
---|
| 230 | |
---|
| 231 | jfile = jfile+1 |
---|
| 232 | CASE(4) reading2_seq ! omeg ---------------------------- |
---|
| 233 | isize = 1 |
---|
| 234 | DO WHILE (isize .le. nsize(iaer,idomain)) |
---|
| 235 | READ(file_unit,*) scanline |
---|
| 236 | IF ((scanline(1:1) .ne. '#').and.& |
---|
| 237 | (scanline(1:1) .ne. ' ')) THEN |
---|
| 238 | BACKSPACE(file_unit) |
---|
| 239 | read(file_unit,*) omeg(:,isize) |
---|
| 240 | isize = isize + 1 |
---|
| 241 | ENDIF |
---|
| 242 | ENDDO |
---|
| 243 | |
---|
| 244 | jfile = jfile+1 |
---|
| 245 | CASE(5) reading2_seq ! gfactor ------------------------- |
---|
| 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,*) gfactor(:,isize) |
---|
| 253 | isize = isize + 1 |
---|
| 254 | ENDIF |
---|
| 255 | ENDDO |
---|
| 256 | |
---|
| 257 | jfile = jfile+1 |
---|
| 258 | IF ((idomain.NE.1).OR.(iaer.NE.1).OR.(jfile.EQ.6)) THEN |
---|
| 259 | endwhile = .true. |
---|
| 260 | ENDIF |
---|
| 261 | CASE(6) reading2_seq |
---|
| 262 | endwhile = .true. |
---|
| 263 | CASE DEFAULT reading2_seq ! ---------------------------- |
---|
| 264 | WRITE(*,*) 'readoptprop: ',& |
---|
| 265 | 'Error while loading optical properties.' |
---|
| 266 | CALL ABORT |
---|
| 267 | END SELECT reading2_seq ! ============================== |
---|
| 268 | ENDIF |
---|
| 269 | ENDDO |
---|
| 270 | |
---|
| 271 | ! 1.4 Close the file |
---|
| 272 | |
---|
| 273 | CLOSE(file_unit) |
---|
| 274 | |
---|
| 275 | ! 1.5 Convert wvl to metres (needed in aerave_new) |
---|
| 276 | ! not needed if already the case |
---|
| 277 | ! do iwvl=1,nwvl |
---|
| 278 | ! wvl(iwvl)=wvl(iwvl)*1.e-6 |
---|
| 279 | ! enddo |
---|
| 280 | |
---|
| 281 | !================================================================== |
---|
| 282 | ! 2. AVERAGED PROPERTIES AND VARIABLE ASSIGNMENTS |
---|
| 283 | !================================================================== |
---|
| 284 | domain: SELECT CASE (idomain) |
---|
| 285 | !================================================================== |
---|
| 286 | CASE(1) domain ! VISIBLE DOMAIN (idomain=1) |
---|
| 287 | !================================================================== |
---|
| 288 | |
---|
| 289 | lamref=lamrefvis(iaer) |
---|
| 290 | epref=1.E+0 |
---|
| 291 | |
---|
| 292 | DO isize=1,nsize(iaer,idomain) |
---|
| 293 | |
---|
| 294 | ! Save the particle sizes |
---|
| 295 | radiustab(iaer,idomain,isize)=radiusdyn(isize) |
---|
| 296 | ! Averaged optical properties (GCM channels) |
---|
| 297 | |
---|
| 298 | CALL aerave_new ( nwvl,& |
---|
| 299 | wvl(:),ep(:,isize),omeg(:,isize),gfactor(:,isize),& |
---|
| 300 | lamref,epref,tstellar,& |
---|
| 301 | L_NSPECTV,blamv,epavVI,& |
---|
| 302 | omegavVI,gavVI,QREFvis(iaer,isize),omegaREFir(iaer,isize)) |
---|
| 303 | |
---|
| 304 | ! Variable assignments (declared in radcommon) |
---|
| 305 | DO isw=1,L_NSPECTV |
---|
| 306 | QVISsQREF(isw,iaer,isize)=epavVI(isw) |
---|
| 307 | gvis(isw,iaer,isize)=gavVI(isw) |
---|
| 308 | omegavis(isw,iaer,isize)=omegavVI(isw) |
---|
| 309 | END DO |
---|
| 310 | ENDDO |
---|
| 311 | |
---|
| 312 | !================================================================== |
---|
| 313 | CASE(2) domain ! INFRARED DOMAIN (idomain=2) |
---|
| 314 | !================================================================== |
---|
| 315 | |
---|
| 316 | |
---|
| 317 | DO isize=1,nsize(iaer,idomain) ! ---------------------------------- |
---|
| 318 | |
---|
| 319 | lamref=lamrefir(iaer) |
---|
| 320 | epref=1.E+0 |
---|
| 321 | |
---|
| 322 | ! Save the particle sizes |
---|
| 323 | radiustab(iaer,idomain,isize)=radiusdyn(isize) |
---|
| 324 | ! radiustab is the table of radius from the table (iaer,idomain,isize) |
---|
| 325 | |
---|
| 326 | ! Averaged optical properties (GCM channels) |
---|
| 327 | |
---|
| 328 | ! epav is <QIR>/Qext(lamrefir) since epref=1 |
---|
| 329 | ! Note: aerave also computes the extinction coefficient at |
---|
| 330 | ! the reference wavelength. This is called QREFvis or QREFir |
---|
| 331 | ! (not epref, which is a different parameter). |
---|
| 332 | ! Reference wavelengths SHOULD BE defined for each aerosol in |
---|
| 333 | ! radcommon_h.F90 |
---|
| 334 | |
---|
| 335 | CALL aerave_new ( nwvl,& |
---|
| 336 | wvl(:),ep(:,isize),omeg(:,isize),gfactor(:,isize),& |
---|
| 337 | lamref,epref,tplanet,& |
---|
| 338 | L_NSPECTI,blami,epavIR,& |
---|
| 339 | omegavIR,gavIR,QREFir(iaer,isize),omegaREFir(iaer,isize)) |
---|
| 340 | |
---|
| 341 | ! Variable assignments (declared in radcommon) |
---|
| 342 | DO ilw=1,L_NSPECTI |
---|
| 343 | QIRsQREF(ilw,iaer,isize)=epavIR(ilw) |
---|
| 344 | gir(ilw,iaer,isize)=gavIR(ilw) |
---|
| 345 | omegair(ilw,iaer,isize)=omegavIR(ilw) |
---|
| 346 | END DO |
---|
| 347 | |
---|
| 348 | ENDDO ! isize (particle size) ------------------------------------- |
---|
| 349 | |
---|
| 350 | END SELECT domain |
---|
| 351 | |
---|
| 352 | |
---|
| 353 | !======================================================================== |
---|
| 354 | ! 3. Deallocate temporary variables that were read in the ASCII files |
---|
| 355 | !======================================================================== |
---|
| 356 | |
---|
| 357 | DEALLOCATE(wvl) ! wvl |
---|
| 358 | DEALLOCATE(radiusdyn) ! radiusdyn |
---|
| 359 | DEALLOCATE(ep) ! ep |
---|
| 360 | DEALLOCATE(omeg) ! omeg |
---|
| 361 | DEALLOCATE(gfactor) ! g |
---|
| 362 | |
---|
| 363 | END DO ! Loop on iaer |
---|
| 364 | END DO ! Loop on idomain |
---|
| 365 | !======================================================================== |
---|
| 366 | RETURN |
---|
| 367 | END subroutine suaer_corrk |
---|
| 368 | |
---|