Changeset 2512 for trunk/LMDZ.MARS/libf/phymars
- Timestamp:
- May 1, 2021, 5:36:56 PM (4 years ago)
- Location:
- trunk/LMDZ.MARS/libf/phymars
- Files:
-
- 3 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.MARS/libf/phymars/callkeys.h
r2447 r2512 17 17 & ,CLFvarying,satindexco2,rdstorm,slpwind,calllott_nonoro & 18 18 & ,latentheat_surfwater,gwd_convective_source,startphy_file & 19 & ,hdo,hdofrac 19 & ,hdo,hdofrac,cap_albedo 20 20 21 21 COMMON/callkeys_i/iradia,iaervar,ilwd,ilwb,ilwn,ncouche & … … 69 69 logical slpwind ! entrainment by slope wind parametrization 70 70 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 72 logical sedimentation 72 73 logical activice,tifeedback,supersat,caps -
trunk/LMDZ.MARS/libf/phymars/conf_phys.F
r2508 r2512 781 781 ! albedo_h2o_frost. Retrocompatible with old 782 782 ! 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 " 784 786 albedo_h2o_cap=0.35 785 787 albedo_h2o_frost=0.35 … … 790 792 call getin_p("albedo_h2o_frost",albedo_h2o_frost) 791 793 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 792 802 ! inert_h2o_ice 793 803 write(*,*) "water ice thermal inertia ?" -
trunk/LMDZ.MARS/libf/phymars/physiq_mod.F
r2511 r2512 49 49 & tsurf, co2ice, emis, 50 50 & capcal, fluxgrd, qsurf, 51 & hmons,summit,base,watercap 51 & hmons,summit,base,watercap,watercaptag 52 52 use comsaison_h, only: dist_sol, declin, mu0, fract, local_time 53 53 use slope_mod, only: theta_sl, psi_sl … … 2072 2072 if ((co2ice(ig).eq.0).and. 2073 2073 & (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 2076 2081 c write(*,*) "frost thickness", qsurf(ig,igcm_h2o_ice) 2077 2082 c write(*,*) "physiq.F frost :"
Note: See TracChangeset
for help on using the changeset viewer.