Ignore:
Timestamp:
May 16, 2017, 9:57:31 AM (8 years ago)
Author:
aslmd
Message:

Added surfalbedo and surfemis keywords to be used with startphy_file = .false. Also made default values in phyetat0_mod unambiguously float

Location:
trunk/LMDZ.GENERIC/libf/phystd
Files:
3 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.GENERIC/libf/phystd/callcorrk.F90

    r1699 r1709  
    2020      USE tracer_h
    2121      use comcstfi_mod, only: pi, mugaz, cpp
    22       use callkeys_mod, only: varactive,diurnal,tracer,water,nosurf,varfixed,satval,        &
     22      use callkeys_mod, only: varactive,diurnal,tracer,water,varfixed,satval,        &
    2323                              kastprof,strictboundcorrk,specOLR,CLFvarying
    2424
     
    494494            albv(nw)=albedo(ig,nw)
    495495         ENDDO
    496 
    497          if (nosurf) then ! Case with no surface.
    498             DO nw=1,L_NSPECTV
    499                if(albv(nw).gt.0.0) then
    500                   print*,'For open lower boundary in callcorrk must'
    501                   print*,'have spectral surface band albedos all set to zero!'
    502                   call abort
    503                endif
    504             ENDDO         
    505          endif
    506496
    507497      if ((ngrid.eq.1).and.(global1d)) then ! Fixed zenith angle 'szangle' in 1D simulations w/ globally-averaged sunlight.
  • trunk/LMDZ.GENERIC/libf/phystd/callkeys_mod.F90

    r1677 r1709  
    107107      real,save :: MassPlanet
    108108!$OMP THREADPRIVATE(flatten,Rmean,J2,MassPlanet)
     109      real,save :: surfalbedo
     110      real,save :: surfemis
     111!$OMP THREADPRIVATE(surfalbedo,surfemis)
    109112
    110113      logical,save :: iscallphys=.false.!existence of callphys.def
  • trunk/LMDZ.GENERIC/libf/phystd/phyetat0_mod.F90

    r1669 r1709  
    1010                     emis,q2,qsurf,cloudfrac,totcloudfrac,hice, &
    1111                     rnat,pctsrf_sic,tslab,tsea_ice,sea_ice)
     12
     13! to use  'getin_p'
     14      use ioipsl_getin_p_mod, only: getin_p
    1215
    1316  use tabfi_mod, only: tabfi
     
    1821                     inquire_dimension, inquire_dimension_length
    1922  use slab_ice_h, only: noceanmx
    20   use callkeys_mod, only: CLFvarying
     23  use callkeys_mod, only: CLFvarying,surfalbedo,surfemis
    2124
    2225  implicit none
     
    115118  endif
    116119else
    117   phisfi(:)=0
     120  phisfi(:)=0.
    118121endif ! of if (startphy_file)
    119122write(*,*) "phyetat0: surface geopotential <phisfi> range:", &
     
    127130  endif
    128131else
    129   albedodat(:)=0.5 ! would be better to read value from def file...
     132  ! If no startfi file, use parameter surfalbedo in def file
     133  surfalbedo=0.5
     134  call getin_p("surfalbedo",surfalbedo)
     135  print*,"surfalbedo",surfalbedo
     136  albedodat(:)=surfalbedo
    130137endif ! of if (startphy_file)
    131138write(*,*) "phyetat0: Bare ground albedo <albedodat> range:", &
     
    139146  endif
    140147else
    141   zmea(:)=0
     148  zmea(:)=0.
    142149endif ! of if (startphy_file)
    143150write(*,*) "phyetat0: <ZMEA> range:", &
     
    151158  endif
    152159else
    153   zstd(:)=0
     160  zstd(:)=0.
    154161endif ! of if (startphy_file)
    155162write(*,*) "phyetat0: <ZSTD> range:", &
     
    163170  endif
    164171else
    165   zsig(:)=0
     172  zsig(:)=0.
    166173endif ! of if (startphy_file)
    167174write(*,*) "phyetat0: <ZSIG> range:", &
     
    175182  endif
    176183else
    177   zgam(:)=0
     184  zgam(:)=0.
    178185endif ! of if (startphy_file)
    179186write(*,*) "phyetat0: <ZGAM> range:", &
     
    187194  endif
    188195else
    189   zthe(:)=0
     196  zthe(:)=0.
    190197endif ! of if (startphy_file)
    191198write(*,*) "phyetat0: <ZTHE> range:", &
     
    199206  endif
    200207else
    201   tsurf(:)=0 ! will be updated afterwards in physiq !
     208  tsurf(:)=0. ! will be updated afterwards in physiq !
    202209endif ! of if (startphy_file)
    203210write(*,*) "phyetat0: Surface temperature <tsurf> range:", &
     
    211218  endif
    212219else
    213   emis(:)=1 ! would be better to read value from def file...
     220  ! If no startfi file, use parameter surfemis in def file
     221  surfemis=1.0
     222  call getin_p("surfemis",surfemis)
     223  print*,"surfemis",surfemis
     224  emis(:)=surfemis
    214225endif ! of if (startphy_file)
    215226write(*,*) "phyetat0: Surface emissivity <emis> range:", &
     
    224235    endif
    225236  else
    226     cloudfrac(:,:)=0
     237    cloudfrac(:,:)=0.0
    227238  endif ! of if (startphy_file)
    228239  write(*,*) "phyetat0: Cloud fraction <cloudfrac> range:", &
     
    240251    endif
    241252  else
    242     totcloudfrac(:)=0
     253    totcloudfrac(:)=0.0
    243254  endif ! of if (startphy_file)
    244255  write(*,*) "phyetat0: Total cloud fraction <totcloudfrac> range:", &
     
    254265    write(*,*) "phyetat0: Failed loading <hice>"
    255266    !  call abort
    256     hice(:)=0
    257   endif
    258 else
    259   hice(:)=0
     267    hice(:)=0.
     268  endif
     269else
     270  hice(:)=0.
    260271endif ! of if (startphy_file)
    261272write(*,*) "phyetat0: Height of oceanic ice <hice> range:", &
     
    279290  endif ! of if (.not.found)
    280291else
    281   rnat(:)=1
     292  rnat(:)=1.
    282293endif ! of if (startphy_file)
    283294write(*,*) "phyetat0: Nature of surface <rnat> range:", &
     
    292303  endif
    293304else
    294   pctsrf_sic(:)=0
     305  pctsrf_sic(:)=0.
    295306endif ! of if (startphy_file)
    296307write(*,*) "phyetat0: Pourcentage of sea ice cover <pctsrf_sic> range:", &
     
    335346  endif
    336347else
    337   tsea_ice(1:ngrid)=0
     348  tsea_ice(1:ngrid)=0.
    338349endif ! of if (startphy_file)
    339350write(*,*) "phyetat0: Oceanic ice quantity <sea_ice> range:", &
     
    348359  endif
    349360else
    350   q2(:,:)=0
     361  q2(:,:)=0.
    351362endif ! of if (startphy_file)
    352363write(*,*) "phyetat0: PBL wind variance <q2> range:", &
     
    365376      endif
    366377    else
    367       qsurf(:,iq)=0
     378      qsurf(:,iq)=0.
    368379    endif ! of if (startphy_file)
    369380    write(*,*) "phyetat0: Surface tracer <",trim(txt),"> range:", &
Note: See TracChangeset for help on using the changeset viewer.