[607] | 1 | subroutine read_dust_scenario(ngrid,nlayer,zday,pplev,tauref) |
---|
| 2 | |
---|
| 3 | ! Reading of the dust scenario file |
---|
| 4 | |
---|
| 5 | use netcdf |
---|
[1543] | 6 | use geometry_mod, only: latitude, longitude ! in radians |
---|
[1918] | 7 | use datafile_mod, only: datadir |
---|
[607] | 8 | implicit none |
---|
| 9 | |
---|
[1861] | 10 | include "callkeys.h" |
---|
[607] | 11 | |
---|
| 12 | integer, intent(in) :: ngrid,nlayer |
---|
| 13 | real, intent(in) :: zday ! date in martian days |
---|
| 14 | real, dimension(ngrid,nlayer+1), intent(in) :: pplev |
---|
| 15 | real, dimension(ngrid), intent(out) :: tauref |
---|
| 16 | |
---|
| 17 | ! Local variables |
---|
| 18 | |
---|
| 19 | real :: realday |
---|
| 20 | integer nid,nvarid,ierr |
---|
| 21 | integer ltloop,lsloop,iloop,jloop,varloop,ig |
---|
| 22 | real, dimension(2) :: taubuf |
---|
| 23 | real tau1(4),tau |
---|
| 24 | real alt(4) |
---|
| 25 | integer latp(4),lonp(4) |
---|
| 26 | integer yinf,ysup,xinf,xsup,tinf,tsup |
---|
| 27 | real latinf,latsup,loninf,lonsup |
---|
| 28 | real latintmp,lonintmp,latdeg,londeg |
---|
| 29 | real colat,dlat,dlon,colattmp |
---|
| 30 | logical, save :: firstcall=.true. |
---|
| 31 | logical :: timeflag |
---|
| 32 | real,save :: radeg,pi |
---|
| 33 | integer :: timedim,londim,latdim |
---|
| 34 | real, dimension(:), allocatable, save :: lat,lon,time |
---|
| 35 | real, dimension(:,:,:), allocatable, save :: tautes |
---|
| 36 | integer, save :: timelen,lonlen,latlen |
---|
| 37 | character(len=33),save :: filename |
---|
| 38 | |
---|
| 39 | realday=mod(zday,669.) |
---|
| 40 | |
---|
[1779] | 41 | ! AS: firstcall OK absolute |
---|
[607] | 42 | if (firstcall) then |
---|
| 43 | firstcall=.false. |
---|
| 44 | |
---|
| 45 | pi=acos(-1.) |
---|
| 46 | radeg=180/pi |
---|
| 47 | |
---|
| 48 | ! assimilated dust file: (NB: iaervar is a common in "callkeys.h") |
---|
| 49 | ! iaervar=4 means read dust_tes.nc file |
---|
[677] | 50 | ! iaervar=6 means read dust_cold.nc file |
---|
| 51 | ! iaervar=7 means read dust_warm.nc file |
---|
[1278] | 52 | ! iaervar=8 means read dust_clim.nc file |
---|
[607] | 53 | ! iaervar=24 means read dust_MY24.nc file |
---|
| 54 | ! iaervar=25 means read dust_MY25.nc file |
---|
| 55 | ! iaervar=26 means read dust_MY26.nc file, etc. |
---|
| 56 | if (iaervar.eq.4) then |
---|
| 57 | ! NB: 4: old TES assimilated MY24 dust scenarios (at 700Pa ref pressure!) |
---|
| 58 | filename="dust_tes.nc" |
---|
[677] | 59 | else if (iaervar.eq.6) then |
---|
| 60 | filename="dust_cold.nc" |
---|
| 61 | else if (iaervar.eq.7) then |
---|
| 62 | filename="dust_warm.nc" |
---|
[1278] | 63 | else if (iaervar.eq.8) then |
---|
| 64 | filename="dust_clim.nc" |
---|
[607] | 65 | else if (iaervar.eq.24) then |
---|
| 66 | filename="dust_MY24.nc" |
---|
| 67 | else if (iaervar.eq.25) then |
---|
| 68 | filename="dust_MY25.nc" |
---|
| 69 | else if (iaervar.eq.26) then |
---|
| 70 | filename="dust_MY26.nc" |
---|
| 71 | else if (iaervar.eq.27) then |
---|
| 72 | filename="dust_MY27.nc" |
---|
| 73 | else if (iaervar.eq.28) then |
---|
| 74 | filename="dust_MY28.nc" |
---|
| 75 | else if (iaervar.eq.29) then |
---|
| 76 | filename="dust_MY29.nc" |
---|
| 77 | else if (iaervar.eq.30) then |
---|
| 78 | filename="dust_MY30.nc" |
---|
[1278] | 79 | else if (iaervar.eq.31) then |
---|
| 80 | filename="dust_MY31.nc" |
---|
[1500] | 81 | else if (iaervar.eq.32) then |
---|
| 82 | filename="dust_MY32.nc" |
---|
[1861] | 83 | else if (iaervar.eq.33) then |
---|
| 84 | filename="dust_MY33.nc" |
---|
[2137] | 85 | else if (iaervar.eq.34) then |
---|
| 86 | filename="dust_MY34.nc" |
---|
[607] | 87 | ! 124,125,126: old TES assimilated dust scenarios (at 700Pa ref pressure!) |
---|
| 88 | else if (iaervar.eq.124) then |
---|
| 89 | filename="dust_tes_MY24.nc" |
---|
| 90 | else if (iaervar.eq.125) then |
---|
| 91 | filename="dust_tes_MY25.nc" |
---|
| 92 | else if (iaervar.eq.126) then |
---|
| 93 | filename="dust_tes_MY26.nc" |
---|
| 94 | endif |
---|
| 95 | |
---|
[1918] | 96 | ierr=nf90_open(trim(datadir)//"/"//trim(filename),NF90_NOWRITE,nid) |
---|
[607] | 97 | IF (ierr.NE.nf90_noerr) THEN |
---|
| 98 | write(*,*)'Problem opening ',trim(filename),' (in phymars/read_dust_scenario.F90)' |
---|
[1918] | 99 | write(*,*)'It should be in :',trim(datadir),'/' |
---|
[607] | 100 | write(*,*)'1) You can change this directory address in callfis.def with' |
---|
| 101 | write(*,*)' datadir=/path/to/datafiles' |
---|
| 102 | write(*,*)'2) If necessary, dust*.nc (and other datafiles)' |
---|
| 103 | write(*,*)' can be obtained online on:' |
---|
[1381] | 104 | write(*,*)' http://www.lmd.jussieu.fr/~lmdz/planets/mars/datadir' |
---|
[607] | 105 | CALL ABORT |
---|
| 106 | ENDIF |
---|
| 107 | |
---|
| 108 | ierr=nf90_inq_dimid(nid,"Time",timedim) |
---|
| 109 | if (ierr.ne.nf90_noerr) then |
---|
| 110 | ! 'Time' dimension not found, look for 'time' |
---|
| 111 | ierr=nf90_inq_dimid(nid,"time",timedim) |
---|
| 112 | if (ierr.ne.nf90_noerr) then |
---|
| 113 | write(*,*)"Error: read_dust_scenario <time> not found" |
---|
| 114 | endif |
---|
| 115 | endif |
---|
| 116 | ierr=nf90_inquire_dimension(nid,timedim,len=timelen) |
---|
| 117 | ierr=nf90_inq_dimid(nid,"latitude",latdim) |
---|
| 118 | ierr=nf90_inquire_dimension(nid,latdim,len=latlen) |
---|
| 119 | ierr=nf90_inq_dimid(nid,"longitude",londim) |
---|
| 120 | ierr=nf90_inquire_dimension(nid,londim,len=lonlen) |
---|
| 121 | |
---|
| 122 | |
---|
[1775] | 123 | if (.not.allocated(tautes)) allocate(tautes(lonlen,latlen,timelen)) |
---|
| 124 | if (.not.allocated(lat)) allocate(lat(latlen), lon(lonlen), time(timelen)) |
---|
[607] | 125 | |
---|
[1156] | 126 | ! "dustop" if loading visible extinction opacity |
---|
| 127 | ! "cdod" if loading IR absorption opacity |
---|
[607] | 128 | ierr=nf90_inq_varid(nid,"dustop",nvarid) |
---|
[1156] | 129 | if (ierr.eq.nf90_noerr) then |
---|
| 130 | ierr=nf90_get_var(nid,nvarid,tautes) |
---|
| 131 | IF (ierr .NE. nf90_noerr) THEN |
---|
| 132 | PRINT*, "Error: read_dust_scenario <dustop> not found" |
---|
| 133 | write(*,*)trim(nf90_strerror(ierr)) |
---|
| 134 | stop |
---|
| 135 | ENDIF |
---|
| 136 | else |
---|
| 137 | ! did not find "dustop" , look for "cdod" |
---|
| 138 | ierr=nf90_inq_varid(nid,"cdod",nvarid) |
---|
| 139 | ierr=nf90_get_var(nid,nvarid,tautes) |
---|
| 140 | IF (ierr .NE. nf90_noerr) THEN |
---|
| 141 | PRINT*, "Error: read_dust_scenario <cdod> not found" |
---|
| 142 | write(*,*)trim(nf90_strerror(ierr)) |
---|
| 143 | stop |
---|
| 144 | ENDIF |
---|
| 145 | ! and multiply by 2*1.3=2.6 to convert from IR absorption |
---|
| 146 | ! to visible extinction opacity |
---|
| 147 | tautes(:,:,:)=2.6*tautes(:,:,:) |
---|
| 148 | endif |
---|
[607] | 149 | |
---|
| 150 | ierr=nf90_inq_varid(nid,"Time",nvarid) |
---|
| 151 | ierr=nf90_get_var(nid,nvarid,time) |
---|
| 152 | IF (ierr .NE. nf90_noerr) THEN |
---|
| 153 | PRINT*, "Error: read_dust_scenario <Time> not found" |
---|
| 154 | write(*,*)trim(nf90_strerror(ierr)) |
---|
| 155 | stop |
---|
| 156 | ENDIF |
---|
| 157 | |
---|
| 158 | ierr=nf90_inq_varid(nid,"latitude",nvarid) |
---|
| 159 | ierr=nf90_get_var(nid,nvarid,lat) |
---|
| 160 | IF (ierr .NE. nf90_noerr) THEN |
---|
| 161 | PRINT*, "Error: read_dust_scenario <latitude> not found" |
---|
| 162 | write(*,*)trim(nf90_strerror(ierr)) |
---|
| 163 | stop |
---|
| 164 | ENDIF |
---|
| 165 | |
---|
| 166 | ierr=nf90_inq_varid(nid,"longitude",nvarid) |
---|
| 167 | ierr=nf90_get_var(nid,nvarid,lon) |
---|
| 168 | IF (ierr .NE. nf90_noerr) THEN |
---|
| 169 | PRINT*, "Error: read_dust_scenario <longitude> not found" |
---|
| 170 | write(*,*)trim(nf90_strerror(ierr)) |
---|
| 171 | stop |
---|
| 172 | ENDIF |
---|
| 173 | |
---|
| 174 | ierr=nf90_close(nid) |
---|
| 175 | |
---|
| 176 | endif ! of if (firstcall) |
---|
| 177 | |
---|
| 178 | do ig=1,ngrid |
---|
| 179 | |
---|
| 180 | ! Find the four nearest points, arranged as follows: |
---|
| 181 | ! 1 2 |
---|
| 182 | ! 3 4 |
---|
| 183 | |
---|
[1541] | 184 | latdeg=latitude(ig)*radeg ! latitude, in degrees |
---|
| 185 | londeg=longitude(ig)*radeg ! longitude, in degrees east |
---|
[607] | 186 | colat=90-latdeg ! colatitude, in degrees |
---|
| 187 | |
---|
| 188 | ! Ehouarn: rounding effects and/or specific compiler issues |
---|
| 189 | ! sometime cause londeg to be sligthly below -180 ... |
---|
| 190 | if (londeg.lt.-180) then |
---|
| 191 | ! check if it is by a large amount |
---|
| 192 | if ((londeg+180).lt.-1.e-3) then |
---|
| 193 | write(*,*) 'reattesassim: error!!' |
---|
| 194 | write(*,*) ' ig=',ig,' londeg=',londeg |
---|
| 195 | stop |
---|
| 196 | else |
---|
| 197 | londeg=-180 |
---|
| 198 | endif |
---|
| 199 | endif |
---|
| 200 | |
---|
| 201 | ! Find encompassing latitudes |
---|
| 202 | if (colat<(90-lat(1))) then ! between north pole and lat(1) |
---|
| 203 | ysup=1 |
---|
| 204 | yinf=1 |
---|
| 205 | else if (colat>=90-(lat(latlen))) then ! between south pole and lat(laten) |
---|
| 206 | ysup=latlen |
---|
| 207 | yinf=latlen |
---|
| 208 | else ! general case |
---|
| 209 | do iloop=2,latlen |
---|
| 210 | if(colat<(90-lat(iloop))) then |
---|
| 211 | ysup=iloop-1 |
---|
| 212 | yinf=iloop |
---|
| 213 | exit |
---|
| 214 | endif |
---|
| 215 | enddo |
---|
| 216 | endif |
---|
| 217 | |
---|
| 218 | latinf=lat(yinf) |
---|
| 219 | latsup=lat(ysup) |
---|
| 220 | |
---|
| 221 | |
---|
| 222 | ! Find encompassing longitudes |
---|
| 223 | ! Note: in input file, lon(1)=-180. |
---|
| 224 | if (londeg>lon(lonlen)) then |
---|
| 225 | xsup=1 |
---|
| 226 | xinf=lonlen |
---|
| 227 | loninf=lon(xsup) |
---|
| 228 | lonsup=180.0 ! since lon(1)=-180.0 |
---|
| 229 | else |
---|
| 230 | do iloop=2,lonlen |
---|
| 231 | if(londeg<lon(iloop)) then |
---|
| 232 | xsup=iloop |
---|
| 233 | xinf=iloop-1 |
---|
| 234 | exit |
---|
| 235 | endif |
---|
| 236 | enddo |
---|
| 237 | loninf=lon(xinf) |
---|
| 238 | lonsup=lon(xsup) |
---|
| 239 | endif |
---|
| 240 | |
---|
| 241 | if ((xsup.gt.lonlen).OR.(yinf.gt.latlen).OR.(xinf.lt.1)& |
---|
| 242 | .OR.(ysup.lt.1)) then |
---|
| 243 | write (*,*) "read_dust_scenario: SYSTEM ERROR on x or y in ts_gcm" |
---|
| 244 | write (*,*) "xinf: ",xinf |
---|
| 245 | write (*,*) "xsup: ",xsup |
---|
| 246 | write (*,*) "yinf: ",yinf |
---|
| 247 | write (*,*) "ysup: ",ysup |
---|
| 248 | stop |
---|
| 249 | endif |
---|
| 250 | |
---|
| 251 | ! loninf=lon(xinf) |
---|
| 252 | ! lonsup=lon(xsup) |
---|
| 253 | ! latinf=lat(yinf) |
---|
| 254 | ! latsup=lat(ysup) |
---|
| 255 | |
---|
| 256 | ! The four neighbouring points are arranged as follows: |
---|
| 257 | ! 1 2 |
---|
| 258 | ! 3 4 |
---|
| 259 | |
---|
| 260 | latp(1)=ysup |
---|
| 261 | latp(2)=ysup |
---|
| 262 | latp(3)=yinf |
---|
| 263 | latp(4)=yinf |
---|
| 264 | |
---|
| 265 | lonp(1)=xinf |
---|
| 266 | lonp(2)=xsup |
---|
| 267 | lonp(3)=xinf |
---|
| 268 | lonp(4)=xsup |
---|
| 269 | |
---|
| 270 | ! Linear interpolation on time, for all four neighbouring points |
---|
| 271 | |
---|
| 272 | if ((realday<time(1)).or.(realday>time(timelen))) then |
---|
| 273 | tinf=timelen |
---|
| 274 | tsup=1 |
---|
| 275 | timeflag=.true. |
---|
| 276 | else |
---|
| 277 | timeflag=.false. |
---|
| 278 | do iloop=2,timelen |
---|
[1544] | 279 | if (realday<=time(iloop)) then |
---|
[607] | 280 | tinf=iloop-1 |
---|
| 281 | tsup=iloop |
---|
| 282 | exit |
---|
| 283 | endif |
---|
| 284 | enddo |
---|
| 285 | endif |
---|
| 286 | |
---|
| 287 | ! Bilinear interpolation on the four nearest points |
---|
| 288 | |
---|
| 289 | if ((colat<(90-lat(1))).OR.(colat>(90-lat(latlen))).OR.(latsup==latinf)) then |
---|
| 290 | dlat=0 |
---|
| 291 | else |
---|
| 292 | dlat=((90-latinf)-colat)/((90-latinf)-(90-latsup)) |
---|
| 293 | endif |
---|
| 294 | |
---|
| 295 | if (lonsup==loninf) then |
---|
| 296 | dlon=0 |
---|
| 297 | else |
---|
| 298 | dlon=(londeg-loninf)/(lonsup-loninf) |
---|
| 299 | endif |
---|
| 300 | |
---|
| 301 | do iloop=1,4 |
---|
| 302 | taubuf(1)=tautes(lonp(iloop),latp(iloop),tinf) |
---|
| 303 | taubuf(2)=tautes(lonp(iloop),latp(iloop),tsup) |
---|
| 304 | if (timeflag) then |
---|
| 305 | if (realday>time(timelen)) then |
---|
| 306 | tau1(iloop)=taubuf(1)+(taubuf(2)-taubuf(1))*(realday-time(tinf))/(time(tsup)+(669-time(tinf))) |
---|
| 307 | else |
---|
| 308 | tau1(iloop)=taubuf(1)+(taubuf(2)-taubuf(1))*realday/(time(tsup)+(669-time(tinf))) |
---|
| 309 | endif |
---|
| 310 | else |
---|
| 311 | tau1(iloop)=taubuf(1)+(taubuf(2)-taubuf(1))*(realday-time(tinf))/(time(tsup)-time(tinf)) |
---|
| 312 | endif |
---|
| 313 | if (tau1(iloop)<0) then |
---|
| 314 | write (*,*) "read_dust_scenario: SYSTEM ERROR on tau" |
---|
| 315 | write (*,*) "utime ",realday |
---|
| 316 | write (*,*) "time(tinf) ",time(tinf) |
---|
| 317 | write (*,*) "time(tsup) ",time(tsup) |
---|
| 318 | write (*,*) "tau1 ",taubuf(1) |
---|
| 319 | write (*,*) "tau2 ",taubuf(2) |
---|
| 320 | write (*,*) "tau ",tau1(iloop) |
---|
| 321 | stop |
---|
| 322 | endif |
---|
| 323 | enddo |
---|
| 324 | |
---|
| 325 | if ((dlat>1).OR.(dlon>1) .OR. (dlat<0) .OR. (dlon<0)) then |
---|
| 326 | write (*,*) "read_dust_scenario: SYSTEM ERROR on dlat or dlon in ts_gcm" |
---|
| 327 | write (*,*) "dlat: ",dlat |
---|
| 328 | write (*,*) "lat: ",latdeg |
---|
| 329 | write (*,*) "dlon: ",dlon |
---|
| 330 | write (*,*) "lon: ",londeg |
---|
| 331 | stop |
---|
| 332 | endif |
---|
| 333 | |
---|
| 334 | tau= dlat*(dlon*(tau1(2)+tau1(3)-tau1(1)-tau1(4))+tau1(1)-tau1(3)) +dlon*(tau1(4)-tau1(3))+tau1(3) |
---|
| 335 | |
---|
| 336 | tauref(ig)=tau |
---|
| 337 | ! |
---|
| 338 | enddo ! of ig=1,ngrid |
---|
| 339 | |
---|
| 340 | if (filename(1:8)=="dust_tes") then |
---|
| 341 | ! when using old TES files, correction for: |
---|
| 342 | ! - the reference pressure was 700Pa (unlike 610Pa now) |
---|
| 343 | ! - the 1.3 conversion factor from IR absorbtion opacity to |
---|
| 344 | ! IR extinction opacity |
---|
| 345 | tauref(:)=tauref(:)*1.3*(610./700.) |
---|
| 346 | endif |
---|
| 347 | |
---|
| 348 | if (swrtype.eq.1) then ! Fouquart (NB: swrtype is set in callkeys.h) |
---|
| 349 | ! when using old radiative transfer (like in MCD 4.x) |
---|
| 350 | ! needed to decrease opacity (*0.825) to compensate overestimation of |
---|
| 351 | ! heating rates |
---|
| 352 | tauref(:)=tauref(:)*0.825/1.3 |
---|
| 353 | endif |
---|
| 354 | |
---|
| 355 | end |
---|