Changeset 2508


Ignore:
Timestamp:
Apr 28, 2021, 5:38:26 PM (4 years ago)
Author:
jnaar
Message:

The water ice albedo is now conceptually decoupled between old perennial ice
and seasonal frost. For now, "old ice" is only on watercaptag = .true.
The variable albedo_h2o_ice has been decomposed into albedo_h2o_cap and
albedo_h2o_frost. Default values are both 0.35 and retrocompatible with
old callphys.def and albedo_h2o_ice if needed.
JN

Location:
trunk/LMDZ.MARS
Files:
6 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.MARS/README

    r2507 r2508  
    33443344of the simulation and the remaining time is in Time (often=0), if we write
    33453345multiple restart nothing changes
     3346
     3347== 22/04/2021 == JN
     3348 - The variable albedo_h2o_ice is now decomposed into albedo_h2o_cap
     3349and albedo_h2o_frost, in order to discriminate between perennial ice
     3350and seasonnal frost. Retrocompatible with old callphys.def and
     3351albedo_h2o_ice. Current default values are 0.35 for both albedos.
     3352
  • trunk/LMDZ.MARS/libf/phymars/albedocaps.F90

    r2304 r2508  
    77use geometry_mod, only: latitude ! grid point latitudes (rad)
    88use surfdat_h, only: TESicealbedo, TESice_Ncoef, TESice_Scoef, &
    9                      emisice, albedice, watercaptag, albedo_h2o_ice, &
     9                     emisice, albedice, watercaptag, albedo_h2o_cap, &
    1010                     emissiv, albedodat
    1111implicit none
     
    7272    ! to do : emissivity
    7373    emisref(ig) = 1
    74     psolaralb(ig,1)=albedo_h2o_ice
    75     psolaralb(ig,2)=albedo_h2o_ice
     74    psolaralb(ig,1)=albedo_h2o_cap
     75    psolaralb(ig,2)=albedo_h2o_cap
    7676  else
    7777    ! set emissivity of surface to be bare ground emissivity
  • trunk/LMDZ.MARS/libf/phymars/conf_phys.F

    r2447 r2508  
    3636      use tracer_mod, only : nuice_sed, ccn_factor, nuiceco2_sed,
    3737     &                       nuice_ref,nuiceco2_ref
    38       use surfdat_h, only: albedo_h2o_ice, inert_h2o_ice,
    39      &                     frost_albedo_threshold
     38      use surfdat_h, only: albedo_h2o_cap,albedo_h2o_frost,
     39     &                     frost_albedo_threshold, inert_h2o_ice
    4040      use time_phylmdz_mod, only: ecritphy,day_step,iphysiq,ecritstart,
    4141     &                            daysec,dtphys
     
    778778         write(*,*) " caps = ",caps
    779779
    780 ! albedo_h2o_ice
     780! JN : now separated between albedo_h2o_cap and
     781! albedo_h2o_frost. Retrocompatible with old
     782! callphys.def with albedo_h2o_ice
    781783         write(*,*) "water ice albedo ?"
    782          albedo_h2o_ice=0.45
    783          call getin_p("albedo_h2o_ice",albedo_h2o_ice)
    784          write(*,*) " albedo_h2o_ice = ",albedo_h2o_ice
     784         albedo_h2o_cap=0.35
     785         albedo_h2o_frost=0.35
     786         call getin_p("albedo_h2o_ice",albedo_h2o_cap)
     787         albedo_h2o_frost=albedo_h2o_cap
     788         call getin_p("albedo_h2o_cap",albedo_h2o_cap)
     789         write(*,*) " albedo_h2o_cap = ",albedo_h2o_cap
     790         call getin_p("albedo_h2o_frost",albedo_h2o_frost)
     791         write(*,*) " albedo_h2o_frost = ",albedo_h2o_frost
    785792! inert_h2o_ice
    786793         write(*,*) "water ice thermal inertia ?"
  • trunk/LMDZ.MARS/libf/phymars/physiq_mod.F

    r2507 r2508  
    4545      use comgeomfi_h, only: sinlon, coslon, sinlat, coslat
    4646      use surfdat_h, only: phisfi, albedodat, zmea, zstd, zsig, zgam,
    47      &                     zthe, z0, albedo_h2o_ice,
     47     &                     zthe, z0, albedo_h2o_cap,albedo_h2o_frost,
    4848     &                     frost_albedo_threshold,
    4949     &                     tsurf, co2ice, emis,
     
    680680
    681681        IF (tracer.AND.water.AND.(ngrid.NE.1)) THEN
    682           write(*,*)"physiq: water_param Surface water ice albedo:",
    683      .                  albedo_h2o_ice
     682          write(*,*)"physiq: water_param Surface water frost albedo:",
     683     .                  albedo_h2o_frost
     684          write(*,*)"physiq: water_param Surface watercap albedo:",
     685     .                  albedo_h2o_cap
    684686        ENDIF
    685687
     
    20702072           if ((co2ice(ig).eq.0).and.
    20712073     &        (qsurf(ig,igcm_h2o_ice).gt.frost_albedo_threshold)) then
    2072               albedo(ig,1) = albedo_h2o_ice
    2073               albedo(ig,2) = albedo_h2o_ice
     2074              albedo(ig,1) = albedo_h2o_frost
     2075              albedo(ig,2) = albedo_h2o_frost
    20742076c              write(*,*) "frost thickness", qsurf(ig,igcm_h2o_ice)
    20752077c              write(*,*) "physiq.F frost :"
  • trunk/LMDZ.MARS/libf/phymars/surfdat_h.F90

    r2260 r2508  
    1313  logical,save :: temptag !temp tag for water caps
    1414     
    15   real,save :: albedo_h2o_ice ! water ice albedo
     15  real,save :: albedo_h2o_cap ! water cap albedo
     16  real,save :: albedo_h2o_frost ! water frost albedo
    1617  real,save :: inert_h2o_ice ! water ice thermal inertia
    1718  real,save :: frost_albedo_threshold ! water frost thickness on the ground (kg.m^-2, ie mm)
  • trunk/LMDZ.MARS/libf/phymars/surfini.F

    r2502 r2508  
    77     &                     cell_area ! for watercaptag diagnosis
    88      use surfdat_h, only: watercaptag, frost_albedo_threshold,
    9      &                     albedo_h2o_ice, inert_h2o_ice, albedodat,
     9     &                     albedo_h2o_cap, inert_h2o_ice, albedodat,
    1010     &                     albedice, dryness
    1111#ifndef MESOSCALE
     
    8484                 watercaptag(ig)  = .true.
    8585                 dryness(ig)      = 1.
    86                  albedodat(ig)    = albedo_h2o_ice  !! pour output
     86                 albedodat(ig)    = albedo_h2o_cap  !! pour output
    8787         else
    8888                 watercaptag(ig)  = .false.
Note: See TracChangeset for help on using the changeset viewer.