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