Ignore:
Timestamp:
May 1, 2021, 5:36:56 PM (4 years ago)
Author:
jnaar
Message:

Added a water-cycle related flag : "cap_albedo" (default is false). When true, the albedo on watercaptag is unchanged by seasonal water frost deposition.
JN

Location:
trunk/LMDZ.MARS/libf/phymars
Files:
3 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.MARS/libf/phymars/callkeys.h

    r2447 r2512  
    1717     &   ,CLFvarying,satindexco2,rdstorm,slpwind,calllott_nonoro        &
    1818     &   ,latentheat_surfwater,gwd_convective_source,startphy_file      &
    19      &   ,hdo,hdofrac
     19     &   ,hdo,hdofrac,cap_albedo
    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
    7172      logical sedimentation
    7273      logical activice,tifeedback,supersat,caps
  • trunk/LMDZ.MARS/libf/phymars/conf_phys.F

    r2508 r2512  
    781781! albedo_h2o_frost. Retrocompatible with old
    782782! callphys.def with albedo_h2o_ice
    783          write(*,*) "water ice albedo ?"
     783         write(*,*) "water ice albedo ? Old settings use ",
     784     &              "albedo_h2o_ice, new settings use ",
     785     &              "albedo_h2o_cap and albedo_h2o_frost "
    784786         albedo_h2o_cap=0.35
    785787         albedo_h2o_frost=0.35
     
    790792         call getin_p("albedo_h2o_frost",albedo_h2o_frost)
    791793         write(*,*) " albedo_h2o_frost = ",albedo_h2o_frost
     794
     795! Northern polar cap albedo (JN 2021)
     796         write(*,*)"Watercaptag albedo is unchanged by water frost",
     797     &              " deposition (default is false)"
     798         cap_albedo=.false. ! default value
     799         call getin_p("cap_albedo",cap_albedo)
     800         write(*,*)"cap_albedo = ",cap_albedo
     801
    792802! inert_h2o_ice
    793803         write(*,*) "water ice thermal inertia ?"
  • trunk/LMDZ.MARS/libf/phymars/physiq_mod.F

    r2511 r2512  
    4949     &                     tsurf, co2ice, emis,
    5050     &                     capcal, fluxgrd, qsurf,
    51      &                     hmons,summit,base,watercap
     51     &                     hmons,summit,base,watercap,watercaptag
    5252      use comsaison_h, only: dist_sol, declin, mu0, fract, local_time
    5353      use slope_mod, only: theta_sl, psi_sl
     
    20722072           if ((co2ice(ig).eq.0).and.
    20732073     &        (qsurf(ig,igcm_h2o_ice).gt.frost_albedo_threshold)) then
    2074               albedo(ig,1) = albedo_h2o_frost
    2075               albedo(ig,2) = albedo_h2o_frost
     2074             if ((watercaptag(ig)).and.(cap_albedo)) then
     2075               albedo(ig,1) = albedo_h2o_cap
     2076               albedo(ig,2) = albedo_h2o_cap
     2077             else
     2078               albedo(ig,1) = albedo_h2o_frost
     2079               albedo(ig,2) = albedo_h2o_frost
     2080             endif !((watercaptag(ig)).and.(cap_albedo)) then
    20762081c              write(*,*) "frost thickness", qsurf(ig,igcm_h2o_ice)
    20772082c              write(*,*) "physiq.F frost :"
Note: See TracChangeset for help on using the changeset viewer.