Changeset 3586 for trunk/LMDZ.MARS/libf/phymars
- Timestamp:
- Jan 20, 2025, 10:31:40 AM (8 days ago)
- Location:
- trunk/LMDZ.MARS/libf/phymars
- Files:
-
- 5 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.MARS/libf/phymars/albedocaps.F90
r3582 r3586 6 6 use ioipsl_getin_p_mod, only: getin_p 7 7 use geometry_mod, only: latitude ! grid point latitudes (rad) 8 use surfdat_h, only: TESicealbedo, TESice_Ncoef, TESice_Scoef, & 9 emisice, albedice, watercaptag, albedo_h2o_cap, & 10 emissiv, albedodat, albedo_perennialco2 8 use surfdat_h, only: TESicealbedo, TESice_Ncoef, TESice_Scoef, & 9 emisice, emissiv, albedice, albedodat, albedo_perennialco2 11 10 USE mod_phys_lmdz_transfert_para, ONLY: bcast 12 11 USE mod_phys_lmdz_para, ONLY: is_master … … 96 95 icap = 1 ! Northern hemisphere 97 96 endif 98 99 if (piceco2(ig) > 0) then ! CO2 frost 97 98 ! Initialization with bare ground 99 ! set emissivity of surface to be bare ground emissivity 100 emisref(ig) = emissiv 101 ! set the surface albedo to bare ground albedo 102 psolaralb(ig,:) = albedodat(ig) 103 104 ! CO2 perennial ice 105 if (paleoclimate .and. piceco2_peren(ig) > 0.) then 106 psolaralb(ig,:) = albedo_perennialco2(icap) 107 emisref(ig) = emisice(icap) 108 endif 109 ! CO2 frost 110 if (piceco2(ig) > 0.) then 100 111 ! set emissivity of surface to be the ice emissivity 101 112 emisref(ig) = emisice(icap) … … 105 116 psolaralb(ig,2) = psolaralb(ig,1) 106 117 else 107 psolaralb(ig,1) = albedice(icap) 108 psolaralb(ig,2) = albedice(icap) 118 psolaralb(ig,:) = albedice(icap) 109 119 endif 110 else if (paleoclimate .and. piceco2_peren(ig) > 0.) then ! CO2 perennial ice 111 psolaralb(ig,1) = albedo_perennialco2(icap) 112 psolaralb(ig,2) = albedo_perennialco2(icap) 113 emisref(ig) = emisice(icap) 114 else if (watercaptag(ig) .and. water) then ! H2O ice in polar caps 115 ! there is a water ice cap: set the surface albedo to the water ice one 116 ! to do : emissivity 117 emisref(ig) = 1 118 psolaralb(ig,:) = albedo_h2o_cap 119 else ! Bare ground 120 ! set emissivity of surface to be bare ground emissivity 121 emisref(ig) = emissiv 122 ! set the surface albedo to bare ground albedo 123 psolaralb(ig,:) = albedodat(ig) 124 endif ! of if (piceco2(ig).gt.0) 120 endif 125 121 enddo ! of ig=1,ngrid 126 122 -
trunk/LMDZ.MARS/libf/phymars/co2condens_mod.F
r3262 r3586 612 612 alb_tmp = psolaralb(:,:,islope) 613 613 emisref_tmp = 0. 614 perennial_co2ice_tmp = perennial_co2ice(:,islope) 614 perennial_co2ice_tmp = perennial_co2ice(:,islope) 615 615 CALL albedocaps(zls,ngrid,piceco2_tmp,perennial_co2ice_tmp, 616 616 & alb_tmp,emisref_tmp) 617 617 perennial_co2ice(:,islope) = perennial_co2ice_tmp 618 psolaralb(:,1,islope) = alb_tmp(:,1) 619 psolaralb(:,2,islope) = alb_tmp(:,2) 618 psolaralb(:,:,islope) = alb_tmp(:,:) 620 619 emisref(:,islope) = emisref_tmp 621 620 ENDDO -
trunk/LMDZ.MARS/libf/phymars/dyn1d/init_testphys1d_mod.F90
r3549 r3586 140 140 obliquit = 25.2 ! Obliquity (deg) ~25.2 141 141 eccentric = 0.0934 ! Eccentricity (0.0934) 142 Lsperi = 251. ! Solar longitude of perihelion 142 143 143 144 ! Planetary Boundary Layer and Turbulence parameters -
trunk/LMDZ.MARS/libf/phymars/dyn1d/testphys1d.F90
r3574 r3586 6 6 use phyredem, only: physdem0, physdem1 7 7 use watersat_mod, only: watersat 8 use tracer_mod, only: igcm_h2o_vap, igcm_h2o_ice, igcm_co2,noms8 use tracer_mod, only: igcm_h2o_vap, igcm_h2o_ice, noms 9 9 use comcstfi_h, only: pi, g, rcp, cpp 10 use time_phylmdz_mod, only: daysec10 use time_phylmdz_mod, only: daysec 11 11 use dimradmars_mod, only: tauvis, totcloudfrac, albedo 12 12 use dust_param_mod, only: tauscaling -
trunk/LMDZ.MARS/libf/phymars/physiq_mod.F
r3468 r3586 2366 2366 c ALWAYS PLACE these lines after co2condens !!! 2367 2367 c ------------------------------------------------------------- 2368 do ig =1,ngrid2368 do ig = 1,ngrid 2369 2369 do islope = 1,nslope 2370 if ((qsurf(ig,igcm_co2,islope).eq.0).and. 2371 & (qsurf(ig,igcm_h2o_ice,islope) 2372 & .gt.frost_albedo_threshold)) then 2373 if ((watercaptag(ig)).and.(cst_cap_albedo)) then 2374 albedo(ig,1,islope) = albedo_h2o_cap 2375 albedo(ig,2,islope) = albedo_h2o_cap 2376 else 2377 albedo(ig,1,islope) = albedo_h2o_frost 2378 albedo(ig,2,islope) = albedo_h2o_frost 2379 endif !((watercaptag(ig)).and.(cst_cap_albedo)) then 2380 c write(*,*) "frost thickness", qsurf(ig,igcm_h2o_ice) 2381 c write(*,*) "physiq.F frost :" 2382 c & ,latitude(ig)*180./pi, longitude(ig)*180./pi 2370 if (abs(qsurf(ig,igcm_co2,islope)) < 1.e-10) then ! No CO2 frost 2371 2372 if (qsurf(ig,igcm_h2o_ice,islope) > frost_albedo_threshold) 2373 & then ! There is H2O frost 2374 if (cst_cap_albedo .and. watercaptag(ig) .and. 2375 & abs(perennial_co2ice(ig,islope)) < 1.e-10) then ! Water cap remains unchanged by water frost deposition and no CO2 perennial ice 2376 albedo(ig,:,islope) = albedo_h2o_cap 2377 emis(ig,islope) = 1. 2378 else 2379 albedo(ig,:,islope) = albedo_h2o_frost 2380 emis(ig,islope) = 1. 2381 endif 2382 else ! No H2O frost 2383 if (abs(perennial_co2ice(ig,islope)) < 1.e-10 .and. 2384 & watercaptag(ig)) then ! No CO2 perennial ice but there is water cap 2385 albedo(ig,:,islope) = albedo_h2o_cap 2386 emis(ig,islope) = 1. 2387 endif 2388 endif 2389 2383 2390 endif 2384 2391 enddo ! islope
Note: See TracChangeset
for help on using the changeset viewer.