Changeset 2561 for trunk


Ignore:
Timestamp:
Sep 13, 2021, 12:19:38 PM (4 years ago)
Author:
jnaar
Message:

Water cycle flag update & addition :

  • "cap_albedo" renamed "cst_cap_albedo" (default false) : if true, water ice cap albedo remains constant even when frost with higher albedo condensates on it
  • "refill_watercap" added (default false) : turns h2o_ice_s into watercap when above a given threshold
  • "frost_metam_threshold" added (default 0.05 km.m-2) : threshold used by "refill_watercap"

JN

Location:
trunk/LMDZ.MARS
Files:
6 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.MARS/README

    r2559 r2561  
    34503450- incorporate auxilliary inistats and mkstats routines in wstats_mod.F90
    34513451- move flag "callstats" from callkeys.h to being a module variable of wstats_mod
     3452
     3453== 13/09/2021 == JN
     3454Water-cycle cleanup, clarifications and additions.
     3455- flag "cap_albedo" renamed "cst_cap_albedo" for clarity. When set to true, watercaptag grid cells albedo
     3456remains constant even when frost condensates on it. Only useful when albedo of frost and cap are different
     3457(default is false)
     3458- addition of two flags "refill_watercap" and "frost_metam_threshold". Refill_watercap is used to convert
     3459h2o_ice_s into watercap if the amount of h2o_ice_s is above the frost_metam_threshold (on watercaptag only).
     3460One could see it as water frost converted into perennial ice.
  • trunk/LMDZ.MARS/libf/phymars/callkeys.h

    r2559 r2561  
    1717     &   ,CLFvarying,satindexco2,rdstorm,slpwind,calllott_nonoro        &
    1818     &   ,latentheat_surfwater,gwd_convective_source,startphy_file      &
    19      &   ,hdo,hdofrac,cap_albedo,temp_dependant_m
     19     &   ,hdo,hdofrac,cst_cap_albedo,temp_dependant_m,refill_watercap
    2020     
    2121      COMMON/callkeys_i/iradia,iaervar,ilwd,ilwb,ilwn,ncouche           &
     
    6969      logical slpwind ! entrainment by slope wind parametrization
    7070      logical latentheat_surfwater ! latent heat release from ground water ice sublimation/condensation
    71       logical cap_albedo ! polar cap albedo remains unchanged by water frost deposition
     71      logical cst_cap_albedo ! polar cap albedo remains unchanged by water frost deposition
    7272      logical temp_dependant_m ! temperature-dependant water contact parameter
     73      logical refill_watercap ! h2o_ice_s is converted to watercap when above threshold
    7374      logical sedimentation
    7475      logical activice,tifeedback,supersat,caps
  • trunk/LMDZ.MARS/libf/phymars/conf_phys.F

    r2559 r2561  
    3737     &                       nuice_ref,nuiceco2_ref
    3838      use surfdat_h, only: albedo_h2o_cap,albedo_h2o_frost,
    39      &                     frost_albedo_threshold, inert_h2o_ice
     39     &                     frost_albedo_threshold, inert_h2o_ice,
     40     &                     frost_metam_threshold
    4041      use time_phylmdz_mod, only: ecritphy,day_step,iphysiq,ecritstart,
    4142     &                            daysec,dtphys
     
    805806         write(*,*)"Watercaptag albedo is unchanged by water frost",
    806807     &              " deposition (default is false)"
    807          cap_albedo=.false. ! default value
    808          call getin_p("cap_albedo",cap_albedo)
    809          write(*,*)"cap_albedo = ",cap_albedo
     808         cst_cap_albedo=.false. ! default value
     809         call getin_p("cst_cap_albedo",cst_cap_albedo)
     810         write(*,*)"cst_cap_albedo = ",cst_cap_albedo
     811
     812! Watercap evolution & refill (with discriminated albedos) (JN 2021)
     813         write(*,*)"Watercap is replenished by water frost",
     814     &              " accumulation (default is false)"
     815         refill_watercap=.false. ! default value
     816         call getin_p("refill_watercap",refill_watercap)
     817         write(*,*)"refill_watercap= ",refill_watercap
     818! frost thickness threshold for refill_watercap (ice metamorphism)
     819         write(*,*) "frost thickness threshold for metamorphism ?",
     820     &              "ie converted into watercap",
     821     &              "only if refill_watercap is .true."
     822         frost_metam_threshold=0.05 !  (i.e 0.05 kg.m-2)
     823         call getin_p("frost_metam_threshold",
     824     &    frost_metam_threshold)
     825         write(*,*) " frost_metam_threshold = ",
     826     &            frost_metam_threshold
     827
    810828
    811829! inert_h2o_ice
  • trunk/LMDZ.MARS/libf/phymars/nuclea.F

    r2522 r2561  
    5757
    5858      IF (temp_dependant_m) THEN
     59c     Simple linear parametrisation from Maattaanen 2014
     60c     Smectite sample
     61c     Maxed out at 0.97 for physical realism
    5962         mtetalocal = min(0.0044*temp + 0.1831,0.97)
    6063      ENDIF ! (temp_dependant_m) THEN
  • trunk/LMDZ.MARS/libf/phymars/physiq_mod.F

    r2559 r2561  
    4646      use surfdat_h, only: phisfi, albedodat, zmea, zstd, zsig, zgam,
    4747     &                     zthe, z0, albedo_h2o_cap,albedo_h2o_frost,
    48      &                     frost_albedo_threshold,
     48     &                     frost_albedo_threshold,frost_metam_threshold,
    4949     &                     tsurf, co2ice, emis,
    5050     &                     capcal, fluxgrd, qsurf,
     
    20782078           if ((co2ice(ig).eq.0).and.
    20792079     &        (qsurf(ig,igcm_h2o_ice).gt.frost_albedo_threshold)) then
    2080              if ((watercaptag(ig)).and.(cap_albedo)) then
     2080             if ((watercaptag(ig)).and.(cst_cap_albedo)) then
    20812081               albedo(ig,1) = albedo_h2o_cap
    20822082               albedo(ig,2) = albedo_h2o_cap
     
    20842084               albedo(ig,1) = albedo_h2o_frost
    20852085               albedo(ig,2) = albedo_h2o_frost
    2086              endif !((watercaptag(ig)).and.(cap_albedo)) then
     2086             endif !((watercaptag(ig)).and.(cst_cap_albedo)) then
    20872087c              write(*,*) "frost thickness", qsurf(ig,igcm_h2o_ice)
    20882088c              write(*,*) "physiq.F frost :"
     
    21512151         watercap(ig)=watercap(ig)+ptimestep*dwatercap(ig)
    21522152      ENDDO
     2153
     2154      IF (refill_watercap) THEN
     2155
     2156        DO ig=1,ngrid
     2157           if (watercaptag(ig).and.
     2158     &        (qsurf(ig,igcm_h2o_ice).gt.frost_metam_threshold)) then
     2159
     2160                watercap(ig)=watercap(ig)+qsurf(ig,igcm_h2o_ice)
     2161     &                  - frost_metam_threshold
     2162                qsurf(ig,igcm_h2o_ice) = frost_metam_threshold
     2163           endif ! (watercaptag(ig).and.
     2164        ENDDO
     2165
     2166      ENDIF ! (refill_watercap) THEN
     2167
    21532168
    21542169c-----------------------------------------------------------------------
  • trunk/LMDZ.MARS/libf/phymars/surfdat_h.F90

    r2508 r2561  
    1717  real,save :: inert_h2o_ice ! water ice thermal inertia
    1818  real,save :: frost_albedo_threshold ! water frost thickness on the ground (kg.m^-2, ie mm)
     19  real,save :: frost_metam_threshold ! water frost threshold before conversion to ice (kg.m^-2, ie mm)
    1920  real,save :: TESice_Ncoef ! coefficient for TES ice albedo in Northern hemisphere
    2021  real,save :: TESice_Scoef ! coefficient for TES ice albedo in Southern hemisphere
Note: See TracChangeset for help on using the changeset viewer.