Changeset 2561 for trunk/LMDZ.MARS/libf/phymars
- Timestamp:
- Sep 13, 2021, 12:19:38 PM (3 years ago)
- Location:
- trunk/LMDZ.MARS/libf/phymars
- Files:
-
- 5 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.MARS/libf/phymars/callkeys.h
r2559 r2561 17 17 & ,CLFvarying,satindexco2,rdstorm,slpwind,calllott_nonoro & 18 18 & ,latentheat_surfwater,gwd_convective_source,startphy_file & 19 & ,hdo,hdofrac,c ap_albedo,temp_dependant_m19 & ,hdo,hdofrac,cst_cap_albedo,temp_dependant_m,refill_watercap 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 c ap_albedo ! polar cap albedo remains unchanged by water frost deposition71 logical cst_cap_albedo ! polar cap albedo remains unchanged by water frost deposition 72 72 logical temp_dependant_m ! temperature-dependant water contact parameter 73 logical refill_watercap ! h2o_ice_s is converted to watercap when above threshold 73 74 logical sedimentation 74 75 logical activice,tifeedback,supersat,caps -
trunk/LMDZ.MARS/libf/phymars/conf_phys.F
r2559 r2561 37 37 & nuice_ref,nuiceco2_ref 38 38 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 40 41 use time_phylmdz_mod, only: ecritphy,day_step,iphysiq,ecritstart, 41 42 & daysec,dtphys … … 805 806 write(*,*)"Watercaptag albedo is unchanged by water frost", 806 807 & " 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 810 828 811 829 ! inert_h2o_ice -
trunk/LMDZ.MARS/libf/phymars/nuclea.F
r2522 r2561 57 57 58 58 IF (temp_dependant_m) THEN 59 c Simple linear parametrisation from Maattaanen 2014 60 c Smectite sample 61 c Maxed out at 0.97 for physical realism 59 62 mtetalocal = min(0.0044*temp + 0.1831,0.97) 60 63 ENDIF ! (temp_dependant_m) THEN -
trunk/LMDZ.MARS/libf/phymars/physiq_mod.F
r2559 r2561 46 46 use surfdat_h, only: phisfi, albedodat, zmea, zstd, zsig, zgam, 47 47 & zthe, z0, albedo_h2o_cap,albedo_h2o_frost, 48 & frost_albedo_threshold, 48 & frost_albedo_threshold,frost_metam_threshold, 49 49 & tsurf, co2ice, emis, 50 50 & capcal, fluxgrd, qsurf, … … 2078 2078 if ((co2ice(ig).eq.0).and. 2079 2079 & (qsurf(ig,igcm_h2o_ice).gt.frost_albedo_threshold)) then 2080 if ((watercaptag(ig)).and.(c ap_albedo)) then2080 if ((watercaptag(ig)).and.(cst_cap_albedo)) then 2081 2081 albedo(ig,1) = albedo_h2o_cap 2082 2082 albedo(ig,2) = albedo_h2o_cap … … 2084 2084 albedo(ig,1) = albedo_h2o_frost 2085 2085 albedo(ig,2) = albedo_h2o_frost 2086 endif !((watercaptag(ig)).and.(c ap_albedo)) then2086 endif !((watercaptag(ig)).and.(cst_cap_albedo)) then 2087 2087 c write(*,*) "frost thickness", qsurf(ig,igcm_h2o_ice) 2088 2088 c write(*,*) "physiq.F frost :" … … 2151 2151 watercap(ig)=watercap(ig)+ptimestep*dwatercap(ig) 2152 2152 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 2153 2168 2154 2169 c----------------------------------------------------------------------- -
trunk/LMDZ.MARS/libf/phymars/surfdat_h.F90
r2508 r2561 17 17 real,save :: inert_h2o_ice ! water ice thermal inertia 18 18 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) 19 20 real,save :: TESice_Ncoef ! coefficient for TES ice albedo in Northern hemisphere 20 21 real,save :: TESice_Scoef ! coefficient for TES ice albedo in Southern hemisphere
Note: See TracChangeset
for help on using the changeset viewer.