source: trunk/LMDZ.MARS/libf/phymars/albedocaps.F90 @ 174

Last change on this file since 174 was 38, checked in by emillour, 14 years ago

Ajout du modè Martien (mon LMDZ.MARS.BETA, du 28/01/2011) dans le rértoire mars, pour pouvoir suivre plus facilement les modifs.
EM

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