source: trunk/LMDZ.MARS/libf/phymars/surfdat_h.F90 @ 1047

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

Mars GCM:

  • IMPORTANT CHANGE: Removed all reference/use of ngridmx (dimphys.h) in routines (necessary prerequisite to using parallel dynamics); in most cases this just means adding 'ngrid' as routine argument, and making local saved variables allocatable (and allocated at first call). In the process, had to convert many *.h files to equivalent modules: yomaer.h => yomaer_h.F90 , surfdat.h => surfdat_h.F90 , comsaison.h => comsaison_h.F90 , yomlw.h => yomlw_h.F90 , comdiurn.h => comdiurn_h.F90 , dimradmars.h => dimradmars_mod.F90 , comgeomfi.h => comgeomfi_h.F90, comsoil.h => comsoil_h.F90 , slope.h => slope_mod.F90
  • Also updated EOF routines, everything is now in eofdump_mod.F90
  • Removed unused routine lectfux.F (in dyn3d)

EM

File size: 1.6 KB
Line 
1module surfdat_h
2
3  ! arrays are allocated in inifis
4  real,save,allocatable :: albedodat(:) ! albedo of bare ground
5  real,save,allocatable :: phisfi(:) ! geopotential at ground level
6  real,save :: albedice(2) ! default albedo for ice (1: North H. 2: South H.)
7  real,save :: emisice(2) ! ice emissivity; 1:Northern hemisphere 2:Southern hemisphere
8  real,save :: emissiv ! emissivity of bare ground
9  logical,save :: TESicealbedo ! use TES ice cap albedoes (if set to .true.)
10  logical,save,allocatable :: watercaptag(:) ! flag for water ice surface
11     
12  logical,save :: temptag !temp tag for water caps
13     
14  real,save :: albedo_h2o_ice ! water ice albedo
15  real,save :: inert_h2o_ice ! water ice thermal inertia
16  real,save :: frost_albedo_threshold ! water frost thickness on the ground (kg.m^-2, ie mm)
17  real,save :: TESice_Ncoef ! coefficient for TES ice albedo in Northern hemisphere
18  real,save :: TESice_Scoef ! coefficient for TES ice albedo in Southern hemisphere
19  real,save :: iceradius(2) , dtemisice(2)
20  real,save,allocatable :: zmea(:),zstd(:),zsig(:),zgam(:),zthe(:)
21  real,save,allocatable :: z0(:) ! surface roughness length (m)
22  real,save :: z0_default ! default (constant over planet) surface roughness (m)
23
24contains
25
26  subroutine ini_surfdat_h(ngrid)
27 
28  implicit none
29  integer,intent(in) :: ngrid ! number of atmospheric columns
30 
31    allocate(albedodat(ngrid))
32    allocate(phisfi(ngrid))
33    allocate(watercaptag(ngrid))
34    allocate(zmea(ngrid))
35    allocate(zstd(ngrid))
36    allocate(zsig(ngrid))
37    allocate(zgam(ngrid))
38    allocate(zthe(ngrid))
39    allocate(z0(ngrid))
40 
41  end subroutine ini_surfdat_h
42
43end module surfdat_h
Note: See TracBrowser for help on using the repository browser.