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

Last change on this file since 1242 was 1130, checked in by emillour, 11 years ago

Mars GCM:
Series of changes to enable running in parallel (using LMDZ.COMMON dynamics);
Current LMDZ.MARS can still notheless be compiled and run in serial mode
"as previously".
Summary of main changes:

  • Main programs (newstart, start2archive, xvik) that used to be in dyn3d have been moved to phymars.
  • dyn3d/control.h is now module control_mod.F90
  • rearanged input/outputs routines everywhere to handle serial/MPI cases. physdem.F => phyredem.F90 , phyetat0.F => phyetat0.F90 ; all read/write routines for startfi files are gathered in module iostart.F90
  • added parallelism related routines init_phys_lmdz.F90, comgeomphy.F90, dimphy.F90, iniphysiq.F90, mod_grid_phy_lmdz.F90, mod_phys_lmdz_mpi_data.F90, mod_phys_lmdz_mpi_transfert.F90, mod_phys_lmdz_omp_data.F90, mod_phys_lmdz_omp_transfert.F90, mod_phys_lmdz_para.F90, mod_phys_lmdz_transfert_para.F90 in phymars and mod_const_mpi.F90 in dyn3d (for compliance with parallel case)
  • created generic routines 'planetwide_maxval' and 'planetwide_minval', in module "planetwide_mod", that enable obtaining the min and max of a field over the whole planet.

EM

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