[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 | |
---|
[2304] | 6 | use ioipsl_getin_p_mod, only: getin_p |
---|
[1543] | 7 | use geometry_mod, only: latitude ! grid point latitudes (rad) |
---|
[1047] | 8 | use surfdat_h, only: TESicealbedo, TESice_Ncoef, TESice_Scoef, & |
---|
[2508] | 9 | emisice, albedice, watercaptag, albedo_h2o_cap, & |
---|
[1047] | 10 | emissiv, albedodat |
---|
[2609] | 11 | USE mod_phys_lmdz_transfert_para, ONLY: bcast |
---|
| 12 | USE mod_phys_lmdz_para, ONLY: is_master |
---|
[38] | 13 | implicit none |
---|
| 14 | |
---|
[1528] | 15 | include"callkeys.h" |
---|
[38] | 16 | |
---|
| 17 | ! arguments: |
---|
| 18 | real,intent(in) :: zls ! solar longitude (rad) |
---|
| 19 | integer,intent(in) :: ngrid |
---|
| 20 | real,intent(in) :: piceco2(ngrid) ! amount of CO2 ice on the surface (kg/m2) |
---|
| 21 | real,intent(out) :: psolaralb(ngrid,2) ! albedo of the surface |
---|
| 22 | real,intent(out) :: emisref(ngrid) ! emissivity of the surface |
---|
| 23 | |
---|
| 24 | |
---|
| 25 | ! local variables: |
---|
| 26 | logical,save :: firstcall=.true. |
---|
[2609] | 27 | !$OMP THREADPRIVATE(firstcall) |
---|
[38] | 28 | integer :: ig,icap |
---|
| 29 | |
---|
[2609] | 30 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 31 | |
---|
| 32 | ! local variables: |
---|
| 33 | real,save :: zls_old ! value of zls from a previous call |
---|
| 34 | real,save :: pi,radeg ! to convert radians to degrees |
---|
| 35 | |
---|
| 36 | !$OMP THREADPRIVATE(zls_old,pi,radeg) |
---|
| 37 | |
---|
| 38 | ! TES datasets: (hard coded fixed length/sizes; for now) |
---|
| 39 | integer,parameter :: TESlonsize=72 |
---|
| 40 | ! longitudes, in TES files, in degrees, from TESlon(1)=-177.5 to TESlon(72)=177.5 |
---|
| 41 | real,save :: TESlon(TESlonsize) |
---|
| 42 | integer,parameter :: TESlatsize=30 |
---|
| 43 | ! latitudes (north hemisphere file), in degrees, from TESlatn(1)=31, |
---|
| 44 | ! to TESlatn(30)=89 ; TESlatn(8)=45 |
---|
| 45 | real,save :: TESlatn(TESlatsize) |
---|
| 46 | ! latitudes (south hemisphere file), in degrees, from TESlats(1)=-89, |
---|
| 47 | ! to TESlats(30)=-31 ; TESlats(23)=-45 |
---|
| 48 | real,save :: TESlats(TESlatsize) |
---|
| 49 | integer,parameter :: TESlssize=72 |
---|
| 50 | ! Solar longitude in TES files, TESls(1)=2.5 to TESls(72)=357.5 |
---|
| 51 | real,save :: TESls(TESlssize) |
---|
| 52 | ! TES North albedo (=-1 for missing values) |
---|
| 53 | real,save :: TESalbn(TESlonsize,TESlatsize,TESlssize) |
---|
| 54 | ! TES South albedo (=-1 for missing values) |
---|
| 55 | real,save :: TESalbs(TESlonsize,TESlatsize,TESlssize) |
---|
| 56 | |
---|
| 57 | !$OMP THREADPRIVATE(TESlon,TESlatn,TESlats,TESls,TESalbn,TESalbs) |
---|
| 58 | |
---|
| 59 | |
---|
| 60 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 61 | |
---|
[38] | 62 | ! 1. Initializations |
---|
[1779] | 63 | ! AS: OK firstcall absolute |
---|
[38] | 64 | if (firstcall) then |
---|
| 65 | ! find out if user wants to use TES cap albedoes or not |
---|
| 66 | TESicealbedo=.false. ! default value |
---|
| 67 | write(*,*)" albedocaps: Use TES Cap albedoes ?" |
---|
[2304] | 68 | call getin_p("TESicealbedo",TESicealbedo) |
---|
[38] | 69 | write(*,*)" albedocaps: TESicealbedo = ",TESicealbedo |
---|
| 70 | |
---|
| 71 | ! if using TES albedoes, load coeffcients |
---|
| 72 | if (TESicealbedo) then |
---|
| 73 | write(*,*)" albedocaps: Coefficient for Northern Cap ?" |
---|
| 74 | TESice_Ncoef=1.0 ! default value |
---|
[2304] | 75 | call getin_p("TESice_Ncoef",TESice_Ncoef) |
---|
[38] | 76 | write(*,*)" albedocaps: TESice_Ncoef = ",TESice_Ncoef |
---|
| 77 | |
---|
| 78 | write(*,*)" albedocaps: Coefficient for Southern Cap ?" |
---|
| 79 | TESice_Scoef=1.0 ! default value |
---|
[2304] | 80 | call getin_p("TESice_Scoef",TESice_Scoef) |
---|
[38] | 81 | write(*,*)" albedocaps: TESice_Scoef = ",TESice_Scoef |
---|
| 82 | endif |
---|
| 83 | |
---|
[2609] | 84 | call read_TES_icecap_albedo(zls_old,pi,radeg,TESlon,TESlatn,TESlats,TESls,TESalbn,TESalbs) |
---|
| 85 | |
---|
[38] | 86 | firstcall=.false. |
---|
| 87 | endif ! of if (firstcall) |
---|
[2609] | 88 | |
---|
[38] | 89 | do ig=1,ngrid |
---|
[1541] | 90 | if (latitude(ig).lt.0.) then |
---|
[38] | 91 | icap=2 ! Southern hemisphere |
---|
| 92 | else |
---|
| 93 | icap=1 ! Northern hemisphere |
---|
| 94 | endif |
---|
[2609] | 95 | |
---|
[38] | 96 | if (piceco2(ig).gt.0) then |
---|
| 97 | ! set emissivity of surface to be the ice emissivity |
---|
| 98 | emisref(ig)=emisice(icap) |
---|
| 99 | ! set the surface albedo to be the ice albedo |
---|
| 100 | if (TESicealbedo) then |
---|
[2609] | 101 | call TES_icecap_albedo(zls,ig,psolaralb(ig,1),icap,zls_old,pi,radeg,TESlon,TESlatn,TESlats,TESls,TESalbn,TESalbs) |
---|
[38] | 102 | psolaralb(ig,2)=psolaralb(ig,1) |
---|
| 103 | else |
---|
| 104 | psolaralb(ig,1)=albedice(icap) |
---|
| 105 | psolaralb(ig,2)=albedice(icap) |
---|
| 106 | endif |
---|
[283] | 107 | else if (watercaptag(ig) .and. water) then |
---|
| 108 | ! there is a water ice cap: set the surface albedo to the water ice one |
---|
| 109 | ! to do : emissivity |
---|
| 110 | emisref(ig) = 1 |
---|
[2508] | 111 | psolaralb(ig,1)=albedo_h2o_cap |
---|
| 112 | psolaralb(ig,2)=albedo_h2o_cap |
---|
[38] | 113 | else |
---|
| 114 | ! set emissivity of surface to be bare ground emissivity |
---|
| 115 | emisref(ig)=emissiv |
---|
| 116 | ! set the surface albedo to bare ground albedo |
---|
| 117 | psolaralb(ig,1)=albedodat(ig) |
---|
| 118 | psolaralb(ig,2)=albedodat(ig) |
---|
| 119 | endif ! of if (piceco2(ig).gt.0) |
---|
| 120 | enddo ! of ig=1,ngrid |
---|
[2609] | 121 | |
---|
[38] | 122 | end subroutine albedocaps |
---|
| 123 | |
---|
| 124 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 125 | |
---|
[2609] | 126 | subroutine read_TES_icecap_albedo(zls_old,pi,radeg,TESlon,TESlatn,TESlats,TESls,TESalbn,TESalbs) |
---|
| 127 | |
---|
[1918] | 128 | use datafile_mod, only: datadir |
---|
[1130] | 129 | use netcdf, only: nf90_open, NF90_NOWRITE, NF90_NOERR, & |
---|
| 130 | nf90_strerror, nf90_inq_varid, nf90_get_var, nf90_close |
---|
[2609] | 131 | USE mod_phys_lmdz_para, ONLY: is_master |
---|
| 132 | USE mod_phys_lmdz_transfert_para, ONLY: bcast |
---|
[1130] | 133 | |
---|
[38] | 134 | implicit none |
---|
| 135 | |
---|
| 136 | ! arguments: |
---|
| 137 | |
---|
| 138 | ! local variables: |
---|
[2609] | 139 | real:: zls_old ! value of zls from a previous call |
---|
| 140 | real:: pi,radeg ! to convert radians to degrees |
---|
[2304] | 141 | character(len=20),parameter :: modname="TES_icecap_albedo" |
---|
[38] | 142 | |
---|
| 143 | ! TES datasets: (hard coded fixed length/sizes; for now) |
---|
[2609] | 144 | integer,parameter:: TESlonsize=72 |
---|
[38] | 145 | ! longitudes, in TES files, in degrees, from TESlon(1)=-177.5 to TESlon(72)=177.5 |
---|
[2609] | 146 | real:: TESlon(TESlonsize) |
---|
| 147 | integer,parameter:: TESlatsize=30 |
---|
[38] | 148 | ! latitudes (north hemisphere file), in degrees, from TESlatn(1)=31, |
---|
| 149 | ! to TESlatn(30)=89 ; TESlatn(8)=45 |
---|
[2609] | 150 | real:: TESlatn(TESlatsize) |
---|
[38] | 151 | ! latitudes (south hemisphere file), in degrees, from TESlats(1)=-89, |
---|
| 152 | ! to TESlats(30)=-31 ; TESlats(23)=-45 |
---|
[2609] | 153 | real:: TESlats(TESlatsize) |
---|
| 154 | integer,parameter:: TESlssize=72 |
---|
[38] | 155 | ! Solar longitude in TES files, TESls(1)=2.5 to TESls(72)=357.5 |
---|
[2609] | 156 | real:: TESls(TESlssize) |
---|
[38] | 157 | ! TES North albedo (=-1 for missing values) |
---|
[2609] | 158 | real:: TESalbn(TESlonsize,TESlatsize,TESlssize) |
---|
[38] | 159 | ! TES South albedo (=-1 for missing values) |
---|
[2609] | 160 | real:: TESalbs(TESlonsize,TESlatsize,TESlssize) |
---|
[38] | 161 | |
---|
| 162 | !NetCDF variables: |
---|
| 163 | integer :: ierr ! NetCDF status |
---|
| 164 | integer :: nid ! NetCDF file ID |
---|
| 165 | integer :: nvarid ! NetCDF variable ID |
---|
| 166 | |
---|
| 167 | ! 0. Preliminary stuff |
---|
[2609] | 168 | |
---|
| 169 | if(is_master) then |
---|
| 170 | |
---|
[38] | 171 | ! Load TES albedoes for Northern Hemisphere |
---|
[1918] | 172 | ierr=nf90_open(trim(datadir)//"/npsc_albedo.nc",NF90_NOWRITE,nid) |
---|
[1130] | 173 | IF (ierr.NE.NF90_NOERR) THEN |
---|
[38] | 174 | write(*,*)'Problem opening npsc_albedo.nc (phymars/albedocaps.F90)' |
---|
[1918] | 175 | write(*,*)'It should be in :',trim(datadir),'/' |
---|
[707] | 176 | write(*,*)'1) You can change this directory address in callfis.def with' |
---|
| 177 | write(*,*)' datadir=/path/to/datafiles' |
---|
[38] | 178 | write(*,*)'2) If necessary, npsc_albedo.nc (and other datafiles)' |
---|
| 179 | write(*,*)' can be obtained online on:' |
---|
[1381] | 180 | write(*,*)' http://www.lmd.jussieu.fr/~lmdz/planets/mars/datadir' |
---|
[2304] | 181 | CALL abort_physic(modname,"missing input file",1) |
---|
[1130] | 182 | ELSE |
---|
[1918] | 183 | write(*,*) "albedocaps: using file ",trim(datadir)//"/npsc_albedo.nc" |
---|
[38] | 184 | ENDIF |
---|
| 185 | |
---|
[1130] | 186 | ierr=nf90_inq_varid(nid,"longitude",nvarid) |
---|
| 187 | if (ierr.ne.NF90_NOERR) then |
---|
[38] | 188 | write(*,*) "Failed to find longitude in file!" |
---|
[1130] | 189 | write(*,*)trim(nf90_strerror(ierr)) |
---|
[2304] | 190 | call abort_physic(modname,"failed finding longitude",1) |
---|
[38] | 191 | else |
---|
[1130] | 192 | ierr=nf90_get_var(nid,nvarid,TESlon) |
---|
| 193 | if (ierr.ne.NF90_NOERR) then |
---|
| 194 | write(*,*) "Failed loading longitude data from file!" |
---|
| 195 | write(*,*)trim(nf90_strerror(ierr)) |
---|
[2304] | 196 | call abort_physic(modname,"failed loading longitude",1) |
---|
[1130] | 197 | endif |
---|
[38] | 198 | endif |
---|
| 199 | |
---|
[1130] | 200 | ierr=nf90_inq_varid(nid,"latitude",nvarid) |
---|
| 201 | if (ierr.ne.NF90_NOERR) then |
---|
[38] | 202 | write(*,*) "Failed to find latitude in file!" |
---|
[1130] | 203 | write(*,*)trim(nf90_strerror(ierr)) |
---|
[2304] | 204 | call abort_physic(modname,"failed finding latitude",1) |
---|
[38] | 205 | else |
---|
[1130] | 206 | ierr=nf90_get_var(nid,nvarid,TESlatn) |
---|
| 207 | if (ierr.ne.NF90_NOERR) then |
---|
| 208 | write(*,*) "Failed loading latitude data from file!" |
---|
| 209 | write(*,*)trim(nf90_strerror(ierr)) |
---|
[2304] | 210 | call abort_physic(modname,"failed loading latitude",1) |
---|
[1130] | 211 | endif |
---|
[38] | 212 | endif |
---|
| 213 | |
---|
[1130] | 214 | ierr=nf90_inq_varid(nid,"time",nvarid) |
---|
| 215 | if (ierr.ne.NF90_NOERR) then |
---|
[38] | 216 | write(*,*) "Failed to find time in file!" |
---|
[1130] | 217 | write(*,*)trim(nf90_strerror(ierr)) |
---|
[2304] | 218 | call abort_physic(modname,"failed finding time",1) |
---|
[38] | 219 | else |
---|
[1130] | 220 | ierr=nf90_get_var(nid,nvarid,TESls) |
---|
| 221 | if (ierr.ne.NF90_NOERR) then |
---|
| 222 | write(*,*) "Failed loading time data from file!" |
---|
| 223 | write(*,*)trim(nf90_strerror(ierr)) |
---|
[2304] | 224 | call abort_physic(modname,"failed loading time",1) |
---|
[1130] | 225 | endif |
---|
[38] | 226 | endif |
---|
| 227 | |
---|
[1130] | 228 | ierr=nf90_inq_varid(nid,"albedo",nvarid) |
---|
| 229 | if (ierr.ne.NF90_NOERR) then |
---|
[38] | 230 | write(*,*) "Failed to find albedo in file!" |
---|
[1130] | 231 | write(*,*)trim(nf90_strerror(ierr)) |
---|
[2304] | 232 | call abort_physic(modname,"failed finding albedo",1) |
---|
[38] | 233 | else |
---|
[1130] | 234 | ierr=nf90_get_var(nid,nvarid,TESalbn) |
---|
| 235 | if (ierr.ne.NF90_NOERR) then |
---|
| 236 | write(*,*) "Failed loading albedo data from file!" |
---|
| 237 | write(*,*)trim(nf90_strerror(ierr)) |
---|
[2304] | 238 | call abort_physic(modname,"failed loading albedo",1) |
---|
[1130] | 239 | endif |
---|
[38] | 240 | endif |
---|
| 241 | |
---|
[1130] | 242 | ierr=nf90_close(nid) |
---|
| 243 | |
---|
[38] | 244 | ! Load albedoes for Southern Hemisphere |
---|
[1918] | 245 | ierr=nf90_open(trim(datadir)//"/spsc_albedo.nc",NF90_NOWRITE,nid) |
---|
[1130] | 246 | IF (ierr.NE.NF90_NOERR) THEN |
---|
[38] | 247 | write(*,*)'Problem opening spsc_albedo.nc (phymars/albedocaps.F90)' |
---|
[1918] | 248 | write(*,*)'It should be in :',trim(datadir),'/' |
---|
[707] | 249 | write(*,*)'1) You can change this directory address in callfis.def with' |
---|
| 250 | write(*,*)' datadir=/path/to/datafiles' |
---|
[38] | 251 | write(*,*)'2) If necessary, spsc_albedo.nc (and other datafiles)' |
---|
| 252 | write(*,*)' can be obtained online on:' |
---|
[1381] | 253 | write(*,*)' http://www.lmd.jussieu.fr/~lmdz/planets/mars/datadir' |
---|
[2304] | 254 | CALL abort_physic(modname,"missing input file",1) |
---|
[1130] | 255 | ELSE |
---|
[1918] | 256 | write(*,*) "albedocaps: using file ",trim(datadir)//"/spsc_albedo.nc" |
---|
[38] | 257 | ENDIF |
---|
| 258 | |
---|
[1130] | 259 | ierr=nf90_inq_varid(nid,"latitude",nvarid) |
---|
| 260 | if (ierr.ne.NF90_NOERR) then |
---|
[38] | 261 | write(*,*) "Failed to find latitude in file!" |
---|
[1130] | 262 | write(*,*)trim(nf90_strerror(ierr)) |
---|
[2304] | 263 | call abort_physic(modname,"failed finding latitude",1) |
---|
[38] | 264 | else |
---|
[1130] | 265 | ierr=nf90_get_var(nid,nvarid,TESlats) |
---|
| 266 | if (ierr.ne.NF90_NOERR) then |
---|
| 267 | write(*,*) "Failed loading latitude data from file!" |
---|
| 268 | write(*,*)trim(nf90_strerror(ierr)) |
---|
[2304] | 269 | call abort_physic(modname,"failed loading latitude",1) |
---|
[1130] | 270 | endif |
---|
[38] | 271 | endif |
---|
| 272 | |
---|
[1130] | 273 | ierr=nf90_inq_varid(nid,"albedo",nvarid) |
---|
| 274 | if (ierr.ne.NF90_NOERR) then |
---|
[38] | 275 | write(*,*) "Failed to find albedo in file!" |
---|
[1130] | 276 | write(*,*)trim(nf90_strerror(ierr)) |
---|
[2304] | 277 | call abort_physic(modname,"failed finding albedo",1) |
---|
[38] | 278 | else |
---|
[1130] | 279 | ierr=nf90_get_var(nid,nvarid,TESalbs) |
---|
| 280 | if (ierr.ne.NF90_NOERR) then |
---|
| 281 | write(*,*) "Failed loading albedo data from file!" |
---|
| 282 | write(*,*)trim(nf90_strerror(ierr)) |
---|
[2304] | 283 | call abort_physic(modname,"failed loading albedo",1) |
---|
[1130] | 284 | endif |
---|
[38] | 285 | endif |
---|
| 286 | |
---|
[1130] | 287 | ierr=nf90_close(nid) |
---|
| 288 | |
---|
[38] | 289 | ! constants: |
---|
| 290 | pi=acos(-1.) |
---|
| 291 | radeg=180/pi |
---|
| 292 | |
---|
| 293 | zls_old=-999 ! dummy initialization |
---|
[2609] | 294 | |
---|
| 295 | endif !is_master |
---|
| 296 | |
---|
| 297 | call bcast(TESlon) |
---|
| 298 | call bcast(TESlatn) |
---|
| 299 | call bcast(TESlats) |
---|
| 300 | call bcast(TESls) |
---|
| 301 | call bcast(TESalbn) |
---|
| 302 | call bcast(TESalbs) |
---|
| 303 | call bcast(pi) |
---|
| 304 | call bcast(zls_old) |
---|
| 305 | call bcast(radeg) |
---|
| 306 | |
---|
| 307 | end subroutine read_TES_icecap_albedo |
---|
[38] | 308 | |
---|
[2609] | 309 | subroutine TES_icecap_albedo(zls,ig,alb,icap,zls_old,pi,radeg,TESlon,TESlatn,TESlats,TESls,TESalbn,TESalbs) |
---|
[38] | 310 | |
---|
[2609] | 311 | use geometry_mod, only: latitude, longitude ! in radians |
---|
| 312 | use surfdat_h, only: albedice, TESice_Ncoef, TESice_Scoef |
---|
| 313 | use datafile_mod, only: datadir |
---|
| 314 | use netcdf, only: nf90_open, NF90_NOWRITE, NF90_NOERR, & |
---|
| 315 | nf90_strerror, nf90_inq_varid, nf90_get_var, nf90_close |
---|
| 316 | USE mod_phys_lmdz_transfert_para, ONLY: bcast |
---|
| 317 | |
---|
| 318 | implicit none |
---|
| 319 | |
---|
| 320 | ! arguments: |
---|
| 321 | real,intent(in) :: zls ! solar longitude (rad) |
---|
| 322 | integer,intent(in) :: ig ! grid point index |
---|
| 323 | real,intent(out) :: alb ! (interpolated) TES ice albedo at that grid point |
---|
| 324 | integer :: icap ! =1: Northern hemisphere =2: Southern hemisphere |
---|
| 325 | |
---|
| 326 | ! local variables: |
---|
| 327 | integer,save :: tinf,tsup ! encompassing time indexes of TES data |
---|
| 328 | real,save :: reltime ! relative position in-between time indexes (in [0;1]) |
---|
| 329 | integer :: latinf,latsup ! encompassing latitude indexes of TES data |
---|
| 330 | real :: rellat ! relative position in-between latitude indexes (in[0;1]) |
---|
| 331 | integer :: loninf,lonsup ! encompassing longitude indexes of TES data |
---|
| 332 | real :: rellon !relative position in-between longitude indexes (in[0;1]) |
---|
| 333 | real :: zlsd ! solar longitude, in degrees |
---|
| 334 | real :: latd ! latitude, in degrees |
---|
| 335 | real :: lond ! longitude, in degrees |
---|
| 336 | integer :: i |
---|
| 337 | |
---|
| 338 | !$OMP THREADPRIVATE(tinf,tsup,reltime) |
---|
| 339 | |
---|
| 340 | ! TES datasets: (hard coded fixed length/sizes; for now) |
---|
| 341 | real,parameter :: TESdeltalon=5.0 ! step in longitude in TES files |
---|
| 342 | real,parameter :: TESdeltalat=2.0 ! step in latitude in TES files |
---|
| 343 | real,parameter :: TESlatnmin=45. ! minimum TES latitude (North hemisphere) |
---|
| 344 | real,parameter :: TESlatsmax=-45. ! maximum TES latitude (South hemisphere) |
---|
| 345 | ! to TESlats(30)=-31 ; TESlats(23)=-45 |
---|
| 346 | real,parameter :: TESdeltals=5.0 ! step in solar longitude in TES files |
---|
| 347 | |
---|
| 348 | real:: zls_old ! value of zls from a previous call |
---|
| 349 | real:: pi,radeg ! to convert radians to degrees |
---|
| 350 | |
---|
| 351 | ! TES datasets: (hard coded fixed length/sizes; for now) |
---|
| 352 | integer,parameter:: TESlonsize=72 |
---|
| 353 | ! longitudes, in TES files, in degrees, from TESlon(1)=-177.5 to TESlon(72)=177.5 |
---|
| 354 | real:: TESlon(TESlonsize) |
---|
| 355 | integer,parameter:: TESlatsize=30 |
---|
| 356 | ! latitudes (north hemisphere file), in degrees, from TESlatn(1)=31, |
---|
| 357 | ! to TESlatn(30)=89 ; TESlatn(8)=45 |
---|
| 358 | real:: TESlatn(TESlatsize) |
---|
| 359 | ! latitudes (south hemisphere file), in degrees, from TESlats(1)=-89, |
---|
| 360 | ! to TESlats(30)=-31 ; TESlats(23)=-45 |
---|
| 361 | real:: TESlats(TESlatsize) |
---|
| 362 | integer,parameter:: TESlssize=72 |
---|
| 363 | ! Solar longitude in TES files, TESls(1)=2.5 to TESls(72)=357.5 |
---|
| 364 | real:: TESls(TESlssize) |
---|
| 365 | ! TES North albedo (=-1 for missing values) |
---|
| 366 | real:: TESalbn(TESlonsize,TESlatsize,TESlssize) |
---|
| 367 | ! TES South albedo (=-1 for missing values) |
---|
| 368 | real:: TESalbs(TESlonsize,TESlatsize,TESlssize) |
---|
| 369 | ! encompassing nodes arranged as follow : 4 3 |
---|
| 370 | real :: val(4) ! 1 2 |
---|
| 371 | |
---|
[801] | 372 | ! 1. Identify encompassing latitudes |
---|
[38] | 373 | |
---|
| 374 | ! Check that latitude is such that there is TES data to use |
---|
| 375 | ! (ie: latitude 45 deg and poleward) otherwise use 'default' albedoes |
---|
[2609] | 376 | |
---|
[1541] | 377 | latd=latitude(ig)*radeg ! latitude, in degrees |
---|
[38] | 378 | if (icap.eq.1) then |
---|
| 379 | ! North hemisphere |
---|
| 380 | if (latd.lt.TESlatnmin) then |
---|
| 381 | alb=albedice(1) |
---|
| 382 | ! the job is done; quit this routine |
---|
| 383 | return |
---|
| 384 | else |
---|
| 385 | ! find encompassing latitudes |
---|
| 386 | if (latd.ge.TESlatn(TESlatsize)) then |
---|
| 387 | latinf=TESlatsize |
---|
| 388 | latsup=TESlatsize |
---|
| 389 | rellat=0. |
---|
| 390 | else |
---|
| 391 | do i=1,TESlatsize-1 |
---|
| 392 | if ((latd.ge.TESlatn(i)).and.(latd.lt.TESlatn(i+1))) then |
---|
| 393 | latinf=i |
---|
| 394 | latsup=i+1 |
---|
| 395 | rellat=(latd-TESlatn(i))/TESdeltalat |
---|
| 396 | exit ! found encompassing indexes; quit loop |
---|
| 397 | endif |
---|
| 398 | enddo |
---|
| 399 | endif |
---|
| 400 | endif ! of if (latd.lt.TESlatnmin) |
---|
| 401 | else ! icap=2 |
---|
| 402 | ! South hemisphere |
---|
| 403 | if (latd.gt.TESlatsmax) then |
---|
| 404 | alb=albedice(2) |
---|
| 405 | ! the job is done; quit this routine |
---|
| 406 | return |
---|
| 407 | else |
---|
| 408 | ! find encompassing latitudes |
---|
| 409 | if (latd.lt.TESlats(1)) then |
---|
| 410 | latinf=1 |
---|
| 411 | latsup=1 |
---|
| 412 | rellat=0. |
---|
| 413 | else |
---|
| 414 | do i=1,TESlatsize-1 |
---|
| 415 | if ((latd.ge.TESlats(i)).and.(latd.lt.TESlats(i+1))) then |
---|
| 416 | latinf=i |
---|
| 417 | latsup=i+1 |
---|
| 418 | rellat=(latd-TESlats(i))/TESdeltalat |
---|
| 419 | exit ! found encompassing indexes; quit loop |
---|
| 420 | endif |
---|
| 421 | enddo |
---|
| 422 | endif |
---|
| 423 | endif ! of if (latd.gt.-45.) |
---|
| 424 | endif ! of if (icap.eq.1) |
---|
| 425 | |
---|
| 426 | ! 2. Identify encompassing time indexes |
---|
| 427 | if (zls.ne.zls_old) then |
---|
| 428 | zlsd=zls*radeg ! solar longitude, in degrees |
---|
| 429 | |
---|
| 430 | if (zlsd.lt.TESls(1)) then |
---|
| 431 | tinf=TESlssize |
---|
| 432 | tsup=1 |
---|
| 433 | reltime=0.5+zlsd/TESdeltals |
---|
| 434 | else |
---|
| 435 | if (zlsd.ge.TESls(TESlssize)) then |
---|
| 436 | tinf=TESlssize |
---|
| 437 | tsup=1 |
---|
| 438 | reltime=(360.-zlsd)/TESdeltals |
---|
| 439 | else |
---|
| 440 | ! look for encompassing indexes |
---|
| 441 | do i=1,TESlssize-1 |
---|
| 442 | if ((zlsd.ge.TESls(i)).and.(zlsd.lt.TESls(i+1))) then |
---|
| 443 | tinf=i |
---|
| 444 | tsup=i+1 |
---|
| 445 | reltime=(zlsd-TESls(i))/TESdeltals |
---|
| 446 | exit ! quit loop, we found the indexes |
---|
| 447 | endif |
---|
| 448 | enddo |
---|
| 449 | endif |
---|
| 450 | endif ! of if (zlsd.lt.TESls(1)) |
---|
| 451 | |
---|
| 452 | zls_old=zls ! store current zls |
---|
| 453 | endif ! of if (zls.ne.zls_old) |
---|
| 454 | |
---|
| 455 | ! 3. Identify encompassing longitudes |
---|
[1541] | 456 | lond=longitude(ig)*radeg ! east longitude, in degrees |
---|
[38] | 457 | if (lond.lt.TESlon(1)) then |
---|
| 458 | loninf=TESlonsize |
---|
| 459 | lonsup=1 |
---|
| 460 | rellon=0.5+(180.+lond)/TESdeltalon |
---|
| 461 | else |
---|
| 462 | if (lond.ge.TESlon(TESlonsize)) then |
---|
| 463 | loninf=TESlonsize |
---|
| 464 | lonsup=1 |
---|
| 465 | rellon=(180-lond)/TESdeltalon |
---|
| 466 | else |
---|
| 467 | do i=1,TESlonsize-1 |
---|
| 468 | if ((lond.ge.TESlon(i)).and.(lond.lt.TESlon(i+1))) then |
---|
| 469 | loninf=i |
---|
| 470 | lonsup=i+1 |
---|
| 471 | rellon=(lond-TESlon(i))/TESdeltalon |
---|
| 472 | exit ! quit loop, we found the indexes |
---|
| 473 | endif |
---|
| 474 | enddo |
---|
| 475 | endif ! of if (lond.ge.TESlon(TESlonsize)) |
---|
| 476 | endif ! of if (lond.lt.TESlon(1)) |
---|
| 477 | |
---|
| 478 | ! 4. Use linear interpolation in time to build encompassing nodal values |
---|
| 479 | ! encompassing nodes are arranged as follow : 4 3 |
---|
| 480 | ! 1 2 |
---|
| 481 | if (icap.eq.1) then |
---|
| 482 | ! Northern hemisphere |
---|
| 483 | val(1)=(1.-reltime)*TESalbn(loninf,latinf,tinf) & |
---|
| 484 | +reltime*TESalbn(loninf,latinf,tsup) |
---|
| 485 | val(2)=(1.-reltime)*TESalbn(lonsup,latinf,tinf) & |
---|
| 486 | +reltime*TESalbn(lonsup,latinf,tsup) |
---|
| 487 | val(3)=(1.-reltime)*TESalbn(lonsup,latsup,tinf) & |
---|
| 488 | +reltime*TESalbn(lonsup,latsup,tsup) |
---|
| 489 | val(4)=(1.-reltime)*TESalbn(loninf,latsup,tinf) & |
---|
| 490 | +reltime*TESalbn(loninf,latsup,tsup) |
---|
| 491 | else |
---|
| 492 | ! Southern hemisphere |
---|
| 493 | val(1)=(1.-reltime)*TESalbs(loninf,latinf,tinf) & |
---|
| 494 | +reltime*TESalbs(loninf,latinf,tsup) |
---|
| 495 | val(2)=(1.-reltime)*TESalbs(lonsup,latinf,tinf) & |
---|
| 496 | +reltime*TESalbs(lonsup,latinf,tsup) |
---|
| 497 | val(3)=(1.-reltime)*TESalbs(lonsup,latsup,tinf) & |
---|
| 498 | +reltime*TESalbs(lonsup,latsup,tsup) |
---|
| 499 | val(4)=(1.-reltime)*TESalbs(loninf,latsup,tinf) & |
---|
| 500 | +reltime*TESalbs(loninf,latsup,tsup) |
---|
| 501 | endif ! of if (icap.eq.1) |
---|
| 502 | |
---|
| 503 | ! 5. Use bilinear interpolation to compute albedo |
---|
| 504 | alb=(1.-rellon)*(1.-rellat)*val(1) & |
---|
| 505 | +rellon*(1.-rellat)*val(2) & |
---|
| 506 | +rellon*rellat*val(3) & |
---|
| 507 | +(1.-rellon)*rellat*val(4) |
---|
| 508 | |
---|
| 509 | ! 6. Apply coefficient to interpolated TES albedo |
---|
| 510 | if (icap.eq.1) then |
---|
| 511 | alb=alb*TESice_Ncoef |
---|
| 512 | else |
---|
| 513 | alb=alb*TESice_Scoef |
---|
| 514 | endif ! of if (icap.eq.1) |
---|
| 515 | |
---|
[707] | 516 | ! Make sure that returned albedo is never greater than 0.90 |
---|
| 517 | if (alb.gt.0.90) alb=0.90 |
---|
| 518 | |
---|
[38] | 519 | end subroutine TES_icecap_albedo |
---|
| 520 | |
---|