[38] | 1 | subroutine albedocaps(zls,ngrid,piceco2,psolaralb,emisref) |
---|
| 2 | |
---|
| 3 | ! routine which changes the albedo (and emissivity) of the surface |
---|
| 4 | ! depending on the presence of CO2 ice on the surface |
---|
| 5 | |
---|
| 6 | ! to use the 'getin' routine |
---|
| 7 | use ioipsl_getincom |
---|
| 8 | |
---|
| 9 | implicit none |
---|
| 10 | |
---|
| 11 | #include"dimensions.h" |
---|
| 12 | #include"dimphys.h" |
---|
| 13 | #include"surfdat.h" |
---|
[283] | 14 | #include"callkeys.h" |
---|
[801] | 15 | #ifdef MESOSCALE |
---|
| 16 | #include"comgeomfi.h" |
---|
| 17 | #endif |
---|
[38] | 18 | |
---|
| 19 | ! arguments: |
---|
| 20 | real,intent(in) :: zls ! solar longitude (rad) |
---|
| 21 | integer,intent(in) :: ngrid |
---|
| 22 | real,intent(in) :: piceco2(ngrid) ! amount of CO2 ice on the surface (kg/m2) |
---|
| 23 | real,intent(out) :: psolaralb(ngrid,2) ! albedo of the surface |
---|
| 24 | real,intent(out) :: emisref(ngrid) ! emissivity of the surface |
---|
| 25 | |
---|
| 26 | |
---|
| 27 | ! local variables: |
---|
| 28 | logical,save :: firstcall=.true. |
---|
| 29 | integer :: ig,icap |
---|
| 30 | |
---|
| 31 | ! 1. Initializations |
---|
| 32 | if (firstcall) then |
---|
| 33 | ! find out if user wants to use TES cap albedoes or not |
---|
| 34 | TESicealbedo=.false. ! default value |
---|
| 35 | write(*,*)" albedocaps: Use TES Cap albedoes ?" |
---|
| 36 | call getin("TESicealbedo",TESicealbedo) |
---|
| 37 | write(*,*)" albedocaps: TESicealbedo = ",TESicealbedo |
---|
| 38 | |
---|
| 39 | ! if using TES albedoes, load coeffcients |
---|
| 40 | if (TESicealbedo) then |
---|
| 41 | write(*,*)" albedocaps: Coefficient for Northern Cap ?" |
---|
| 42 | TESice_Ncoef=1.0 ! default value |
---|
| 43 | call getin("TESice_Ncoef",TESice_Ncoef) |
---|
| 44 | write(*,*)" albedocaps: TESice_Ncoef = ",TESice_Ncoef |
---|
| 45 | |
---|
| 46 | write(*,*)" albedocaps: Coefficient for Southern Cap ?" |
---|
| 47 | TESice_Scoef=1.0 ! default value |
---|
| 48 | call getin("TESice_Scoef",TESice_Scoef) |
---|
| 49 | write(*,*)" albedocaps: TESice_Scoef = ",TESice_Scoef |
---|
| 50 | endif |
---|
| 51 | |
---|
| 52 | firstcall=.false. |
---|
| 53 | endif ! of if (firstcall) |
---|
| 54 | |
---|
| 55 | do ig=1,ngrid |
---|
[801] | 56 | #ifndef MESOSCALE |
---|
[38] | 57 | if (ig.gt.ngrid/2+1) then |
---|
[801] | 58 | #else |
---|
| 59 | if (lati(ig)*180./acos(-1.).lt.0.) then |
---|
| 60 | #endif |
---|
[38] | 61 | icap=2 ! Southern hemisphere |
---|
| 62 | else |
---|
| 63 | icap=1 ! Northern hemisphere |
---|
| 64 | endif |
---|
| 65 | |
---|
| 66 | if (piceco2(ig).gt.0) then |
---|
| 67 | ! set emissivity of surface to be the ice emissivity |
---|
| 68 | emisref(ig)=emisice(icap) |
---|
| 69 | ! set the surface albedo to be the ice albedo |
---|
| 70 | if (TESicealbedo) then |
---|
| 71 | ! write(*,*) "albedocaps: call TES_icecap_albedo" |
---|
| 72 | ! write(*,*) "albedocaps: zls=",zls," ig=",ig |
---|
[801] | 73 | call TES_icecap_albedo(zls,ig,psolaralb(ig,1),icap) |
---|
[38] | 74 | ! write(*,*) "albedocaps: psolaralb(ig,1)=",psolaralb(ig,1) |
---|
| 75 | psolaralb(ig,2)=psolaralb(ig,1) |
---|
| 76 | else |
---|
| 77 | psolaralb(ig,1)=albedice(icap) |
---|
| 78 | psolaralb(ig,2)=albedice(icap) |
---|
| 79 | endif |
---|
[283] | 80 | else if (watercaptag(ig) .and. water) then |
---|
| 81 | ! there is a water ice cap: set the surface albedo to the water ice one |
---|
| 82 | ! to do : emissivity |
---|
| 83 | !write(*,*) "watercaptag in albedocaps:",ig |
---|
| 84 | emisref(ig) = 1 |
---|
| 85 | psolaralb(ig,1)=albedo_h2o_ice |
---|
| 86 | psolaralb(ig,2)=albedo_h2o_ice |
---|
[38] | 87 | else |
---|
| 88 | ! set emissivity of surface to be bare ground emissivity |
---|
| 89 | emisref(ig)=emissiv |
---|
| 90 | ! set the surface albedo to bare ground albedo |
---|
| 91 | psolaralb(ig,1)=albedodat(ig) |
---|
| 92 | psolaralb(ig,2)=albedodat(ig) |
---|
| 93 | endif ! of if (piceco2(ig).gt.0) |
---|
| 94 | enddo ! of ig=1,ngrid |
---|
| 95 | end subroutine albedocaps |
---|
| 96 | |
---|
| 97 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
[801] | 98 | subroutine TES_icecap_albedo(zls,ig,alb,icap) |
---|
[38] | 99 | |
---|
| 100 | implicit none |
---|
| 101 | #include"dimensions.h" |
---|
| 102 | #include"dimphys.h" |
---|
| 103 | #include"surfdat.h" |
---|
| 104 | #include"comgeomfi.h" |
---|
| 105 | #include"netcdf.inc" |
---|
| 106 | #include"datafile.h" |
---|
| 107 | |
---|
| 108 | ! arguments: |
---|
| 109 | real,intent(in) :: zls ! solar longitude (rad) |
---|
| 110 | integer,intent(in) :: ig ! grid point index |
---|
| 111 | real,intent(out) :: alb ! (interpolated) TES ice albedo at that grid point |
---|
[801] | 112 | integer :: icap ! =1: Northern hemisphere =2: Southern hemisphere |
---|
[38] | 113 | |
---|
| 114 | ! local variables: |
---|
| 115 | logical,save :: firstcall=.true. |
---|
| 116 | real,save :: zls_old ! value of zls from a previous call |
---|
| 117 | integer,save :: tinf,tsup ! encompassing time indexes of TES data |
---|
| 118 | real,save :: reltime ! relative position in-between time indexes (in [0;1]) |
---|
| 119 | integer :: latinf,latsup ! encompassing latitude indexes of TES data |
---|
| 120 | real :: rellat ! relative position in-between latitude indexes (in[0;1]) |
---|
| 121 | integer :: loninf,lonsup ! encompassing longitude indexes of TES data |
---|
| 122 | real :: rellon !relative position in-between longitude indexes (in[0;1]) |
---|
| 123 | real,save :: pi,radeg ! to convert radians to degrees |
---|
| 124 | real :: zlsd ! solar longitude, in degrees |
---|
| 125 | real :: latd ! latitude, in degrees |
---|
| 126 | real :: lond ! longitude, in degrees |
---|
| 127 | integer :: i |
---|
| 128 | |
---|
| 129 | ! TES datasets: (hard coded fixed length/sizes; for now) |
---|
| 130 | integer,parameter :: TESlonsize=72 |
---|
| 131 | real,parameter :: TESdeltalon=5.0 ! step in longitude in TES files |
---|
| 132 | ! longitudes, in TES files, in degrees, from TESlon(1)=-177.5 to TESlon(72)=177.5 |
---|
| 133 | real,save :: TESlon(TESlonsize) |
---|
| 134 | integer,parameter :: TESlatsize=30 |
---|
| 135 | real,parameter :: TESdeltalat=2.0 ! step in latitude in TES files |
---|
| 136 | ! latitudes (north hemisphere file), in degrees, from TESlatn(1)=31, |
---|
| 137 | ! to TESlatn(30)=89 ; TESlatn(8)=45 |
---|
| 138 | real,parameter :: TESlatnmin=45. ! minimum TES latitude (North hemisphere) |
---|
| 139 | real,parameter :: TESlatsmax=-45. ! maximum TES latitude (South hemisphere) |
---|
| 140 | real,save :: TESlatn(TESlatsize) |
---|
| 141 | ! latitudes (south hemisphere file), in degrees, from TESlats(1)=-89, |
---|
| 142 | ! to TESlats(30)=-31 ; TESlats(23)=-45 |
---|
| 143 | real,save :: TESlats(TESlatsize) |
---|
| 144 | integer,parameter :: TESlssize=72 |
---|
| 145 | real,parameter :: TESdeltals=5.0 ! step in solar longitude in TES files |
---|
| 146 | ! Solar longitude in TES files, TESls(1)=2.5 to TESls(72)=357.5 |
---|
| 147 | real,save :: TESls(TESlssize) |
---|
| 148 | ! TES North albedo (=-1 for missing values) |
---|
| 149 | real,save :: TESalbn(TESlonsize,TESlatsize,TESlssize) |
---|
| 150 | ! TES South albedo (=-1 for missing values) |
---|
| 151 | real,save :: TESalbs(TESlonsize,TESlatsize,TESlssize) |
---|
| 152 | ! encompassing nodes arranged as follow : 4 3 |
---|
| 153 | real :: val(4) ! 1 2 |
---|
| 154 | |
---|
| 155 | !NetCDF variables: |
---|
| 156 | integer :: ierr ! NetCDF status |
---|
| 157 | integer :: nid ! NetCDF file ID |
---|
| 158 | integer :: nvarid ! NetCDF variable ID |
---|
| 159 | |
---|
| 160 | ! 0. Preliminary stuff |
---|
| 161 | if (firstcall) then |
---|
| 162 | ! Load TES albedoes for Northern Hemisphere |
---|
| 163 | ! Note: datafile() is defined in "datafile.h" |
---|
| 164 | ierr=NF_OPEN(trim(datafile)//"/npsc_albedo.nc",NF_NOWRITE,nid) |
---|
| 165 | IF (ierr.NE.NF_NOERR) THEN |
---|
| 166 | write(*,*)'Problem opening npsc_albedo.nc (phymars/albedocaps.F90)' |
---|
| 167 | write(*,*)'It should be in :',trim(datafile),'/' |
---|
[707] | 168 | write(*,*)'1) You can change this directory address in callfis.def with' |
---|
| 169 | write(*,*)' datadir=/path/to/datafiles' |
---|
[38] | 170 | write(*,*)'2) If necessary, npsc_albedo.nc (and other datafiles)' |
---|
| 171 | write(*,*)' can be obtained online on:' |
---|
| 172 | write(*,*)' http://www.lmd.jussieu.fr/~forget/datagcm/datafile' |
---|
| 173 | CALL ABORT |
---|
| 174 | ENDIF |
---|
| 175 | |
---|
| 176 | ierr=NF_INQ_VARID(nid,"longitude",nvarid) |
---|
| 177 | if (ierr.ne.NF_NOERR) then |
---|
| 178 | write(*,*) "Failed to find longitude in file!" |
---|
| 179 | else |
---|
| 180 | ierr=NF_GET_VAR_REAL(nid,nvarid,TESlon) |
---|
| 181 | endif |
---|
| 182 | |
---|
| 183 | ierr=NF_INQ_VARID(nid,"latitude",nvarid) |
---|
| 184 | if (ierr.ne.NF_NOERR) then |
---|
| 185 | write(*,*) "Failed to find latitude in file!" |
---|
| 186 | else |
---|
| 187 | ierr=NF_GET_VAR_REAL(nid,nvarid,TESlatn) |
---|
| 188 | endif |
---|
| 189 | |
---|
| 190 | ierr=NF_INQ_VARID(nid,"time",nvarid) |
---|
| 191 | if (ierr.ne.NF_NOERR) then |
---|
| 192 | write(*,*) "Failed to find time in file!" |
---|
| 193 | else |
---|
| 194 | ierr=NF_GET_VAR_REAL(nid,nvarid,TESls) |
---|
| 195 | endif |
---|
| 196 | |
---|
| 197 | ierr=NF_INQ_VARID(nid,"albedo",nvarid) |
---|
| 198 | if (ierr.ne.NF_NOERR) then |
---|
| 199 | write(*,*) "Failed to find albedo in file!" |
---|
| 200 | else |
---|
| 201 | ierr=NF_GET_VAR_REAL(nid,nvarid,TESalbn) |
---|
| 202 | endif |
---|
| 203 | |
---|
| 204 | ! Load albedoes for Southern Hemisphere |
---|
| 205 | ierr=NF_OPEN(trim(datafile)//"/spsc_albedo.nc",NF_NOWRITE,nid) |
---|
| 206 | IF (ierr.NE.NF_NOERR) THEN |
---|
| 207 | write(*,*)'Problem opening spsc_albedo.nc (phymars/albedocaps.F90)' |
---|
| 208 | write(*,*)'It should be in :',trim(datafile),'/' |
---|
[707] | 209 | write(*,*)'1) You can change this directory address in callfis.def with' |
---|
| 210 | write(*,*)' datadir=/path/to/datafiles' |
---|
[38] | 211 | write(*,*)'2) If necessary, spsc_albedo.nc (and other datafiles)' |
---|
| 212 | write(*,*)' can be obtained online on:' |
---|
| 213 | write(*,*)' http://www.lmd.jussieu.fr/~forget/datagcm/datafile' |
---|
| 214 | CALL ABORT |
---|
| 215 | ENDIF |
---|
| 216 | |
---|
| 217 | ierr=NF_INQ_VARID(nid,"latitude",nvarid) |
---|
| 218 | if (ierr.ne.NF_NOERR) then |
---|
| 219 | write(*,*) "Failed to find latitude in file!" |
---|
| 220 | else |
---|
| 221 | ierr=NF_GET_VAR_REAL(nid,nvarid,TESlats) |
---|
| 222 | endif |
---|
| 223 | |
---|
| 224 | ierr=NF_INQ_VARID(nid,"albedo",nvarid) |
---|
| 225 | if (ierr.ne.NF_NOERR) then |
---|
| 226 | write(*,*) "Failed to find albedo in file!" |
---|
| 227 | else |
---|
| 228 | ierr=NF_GET_VAR_REAL(nid,nvarid,TESalbs) |
---|
| 229 | endif |
---|
| 230 | |
---|
| 231 | ! constants: |
---|
| 232 | pi=acos(-1.) |
---|
| 233 | radeg=180/pi |
---|
| 234 | |
---|
| 235 | zls_old=-999 ! dummy initialization |
---|
| 236 | |
---|
| 237 | firstcall=.false. |
---|
| 238 | endif ! of if firstcall |
---|
| 239 | |
---|
[801] | 240 | ! 1. Identify encompassing latitudes |
---|
[38] | 241 | |
---|
| 242 | ! Check that latitude is such that there is TES data to use |
---|
| 243 | ! (ie: latitude 45 deg and poleward) otherwise use 'default' albedoes |
---|
| 244 | latd=lati(ig)*radeg ! latitude, in degrees |
---|
| 245 | if (icap.eq.1) then |
---|
| 246 | ! North hemisphere |
---|
| 247 | if (latd.lt.TESlatnmin) then |
---|
| 248 | alb=albedice(1) |
---|
| 249 | ! the job is done; quit this routine |
---|
| 250 | return |
---|
| 251 | else |
---|
| 252 | ! find encompassing latitudes |
---|
| 253 | if (latd.ge.TESlatn(TESlatsize)) then |
---|
| 254 | latinf=TESlatsize |
---|
| 255 | latsup=TESlatsize |
---|
| 256 | rellat=0. |
---|
| 257 | else |
---|
| 258 | do i=1,TESlatsize-1 |
---|
| 259 | if ((latd.ge.TESlatn(i)).and.(latd.lt.TESlatn(i+1))) then |
---|
| 260 | latinf=i |
---|
| 261 | latsup=i+1 |
---|
| 262 | rellat=(latd-TESlatn(i))/TESdeltalat |
---|
| 263 | exit ! found encompassing indexes; quit loop |
---|
| 264 | endif |
---|
| 265 | enddo |
---|
| 266 | endif |
---|
| 267 | endif ! of if (latd.lt.TESlatnmin) |
---|
| 268 | else ! icap=2 |
---|
| 269 | ! South hemisphere |
---|
| 270 | if (latd.gt.TESlatsmax) then |
---|
| 271 | alb=albedice(2) |
---|
| 272 | ! the job is done; quit this routine |
---|
| 273 | return |
---|
| 274 | else |
---|
| 275 | ! find encompassing latitudes |
---|
| 276 | if (latd.lt.TESlats(1)) then |
---|
| 277 | latinf=1 |
---|
| 278 | latsup=1 |
---|
| 279 | rellat=0. |
---|
| 280 | else |
---|
| 281 | do i=1,TESlatsize-1 |
---|
| 282 | if ((latd.ge.TESlats(i)).and.(latd.lt.TESlats(i+1))) then |
---|
| 283 | latinf=i |
---|
| 284 | latsup=i+1 |
---|
| 285 | rellat=(latd-TESlats(i))/TESdeltalat |
---|
| 286 | exit ! found encompassing indexes; quit loop |
---|
| 287 | endif |
---|
| 288 | enddo |
---|
| 289 | endif |
---|
| 290 | endif ! of if (latd.gt.-45.) |
---|
| 291 | endif ! of if (icap.eq.1) |
---|
| 292 | |
---|
| 293 | ! 2. Identify encompassing time indexes |
---|
| 294 | if (zls.ne.zls_old) then |
---|
| 295 | zlsd=zls*radeg ! solar longitude, in degrees |
---|
| 296 | |
---|
| 297 | if (zlsd.lt.TESls(1)) then |
---|
| 298 | tinf=TESlssize |
---|
| 299 | tsup=1 |
---|
| 300 | reltime=0.5+zlsd/TESdeltals |
---|
| 301 | else |
---|
| 302 | if (zlsd.ge.TESls(TESlssize)) then |
---|
| 303 | tinf=TESlssize |
---|
| 304 | tsup=1 |
---|
| 305 | reltime=(360.-zlsd)/TESdeltals |
---|
| 306 | else |
---|
| 307 | ! look for encompassing indexes |
---|
| 308 | do i=1,TESlssize-1 |
---|
| 309 | if ((zlsd.ge.TESls(i)).and.(zlsd.lt.TESls(i+1))) then |
---|
| 310 | tinf=i |
---|
| 311 | tsup=i+1 |
---|
| 312 | reltime=(zlsd-TESls(i))/TESdeltals |
---|
| 313 | exit ! quit loop, we found the indexes |
---|
| 314 | endif |
---|
| 315 | enddo |
---|
| 316 | endif |
---|
| 317 | endif ! of if (zlsd.lt.TESls(1)) |
---|
| 318 | |
---|
| 319 | zls_old=zls ! store current zls |
---|
| 320 | endif ! of if (zls.ne.zls_old) |
---|
| 321 | |
---|
| 322 | ! 3. Identify encompassing longitudes |
---|
| 323 | lond=long(ig)*radeg ! east longitude, in degrees |
---|
| 324 | if (lond.lt.TESlon(1)) then |
---|
| 325 | loninf=TESlonsize |
---|
| 326 | lonsup=1 |
---|
| 327 | rellon=0.5+(180.+lond)/TESdeltalon |
---|
| 328 | else |
---|
| 329 | if (lond.ge.TESlon(TESlonsize)) then |
---|
| 330 | loninf=TESlonsize |
---|
| 331 | lonsup=1 |
---|
| 332 | rellon=(180-lond)/TESdeltalon |
---|
| 333 | else |
---|
| 334 | do i=1,TESlonsize-1 |
---|
| 335 | if ((lond.ge.TESlon(i)).and.(lond.lt.TESlon(i+1))) then |
---|
| 336 | loninf=i |
---|
| 337 | lonsup=i+1 |
---|
| 338 | rellon=(lond-TESlon(i))/TESdeltalon |
---|
| 339 | exit ! quit loop, we found the indexes |
---|
| 340 | endif |
---|
| 341 | enddo |
---|
| 342 | endif ! of if (lond.ge.TESlon(TESlonsize)) |
---|
| 343 | endif ! of if (lond.lt.TESlon(1)) |
---|
| 344 | |
---|
| 345 | ! 4. Use linear interpolation in time to build encompassing nodal values |
---|
| 346 | ! encompassing nodes are arranged as follow : 4 3 |
---|
| 347 | ! 1 2 |
---|
| 348 | if (icap.eq.1) then |
---|
| 349 | ! Northern hemisphere |
---|
| 350 | val(1)=(1.-reltime)*TESalbn(loninf,latinf,tinf) & |
---|
| 351 | +reltime*TESalbn(loninf,latinf,tsup) |
---|
| 352 | val(2)=(1.-reltime)*TESalbn(lonsup,latinf,tinf) & |
---|
| 353 | +reltime*TESalbn(lonsup,latinf,tsup) |
---|
| 354 | val(3)=(1.-reltime)*TESalbn(lonsup,latsup,tinf) & |
---|
| 355 | +reltime*TESalbn(lonsup,latsup,tsup) |
---|
| 356 | val(4)=(1.-reltime)*TESalbn(loninf,latsup,tinf) & |
---|
| 357 | +reltime*TESalbn(loninf,latsup,tsup) |
---|
| 358 | else |
---|
| 359 | ! Southern hemisphere |
---|
| 360 | val(1)=(1.-reltime)*TESalbs(loninf,latinf,tinf) & |
---|
| 361 | +reltime*TESalbs(loninf,latinf,tsup) |
---|
| 362 | val(2)=(1.-reltime)*TESalbs(lonsup,latinf,tinf) & |
---|
| 363 | +reltime*TESalbs(lonsup,latinf,tsup) |
---|
| 364 | val(3)=(1.-reltime)*TESalbs(lonsup,latsup,tinf) & |
---|
| 365 | +reltime*TESalbs(lonsup,latsup,tsup) |
---|
| 366 | val(4)=(1.-reltime)*TESalbs(loninf,latsup,tinf) & |
---|
| 367 | +reltime*TESalbs(loninf,latsup,tsup) |
---|
| 368 | endif ! of if (icap.eq.1) |
---|
| 369 | |
---|
| 370 | ! 5. Use bilinear interpolation to compute albedo |
---|
| 371 | alb=(1.-rellon)*(1.-rellat)*val(1) & |
---|
| 372 | +rellon*(1.-rellat)*val(2) & |
---|
| 373 | +rellon*rellat*val(3) & |
---|
| 374 | +(1.-rellon)*rellat*val(4) |
---|
| 375 | |
---|
| 376 | ! 6. Apply coefficient to interpolated TES albedo |
---|
| 377 | if (icap.eq.1) then |
---|
| 378 | alb=alb*TESice_Ncoef |
---|
| 379 | else |
---|
| 380 | alb=alb*TESice_Scoef |
---|
| 381 | endif ! of if (icap.eq.1) |
---|
| 382 | |
---|
[707] | 383 | ! Make sure that returned albedo is never greater than 0.90 |
---|
| 384 | if (alb.gt.0.90) alb=0.90 |
---|
| 385 | |
---|
[38] | 386 | end subroutine TES_icecap_albedo |
---|
| 387 | |
---|