Changeset 2999 for trunk/LMDZ.MARS/libf/phymars
- Timestamp:
- Jul 19, 2023, 11:40:38 AM (2 years ago)
- Location:
- trunk/LMDZ.MARS/libf/phymars
- Files:
-
- 8 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.MARS/libf/phymars/albedocaps.F90
r2609 r2999 1 subroutine albedocaps(zls,ngrid,piceco2,p solaralb,emisref)1 subroutine albedocaps(zls,ngrid,piceco2,piceco2_peren,psolaralb,emisref) 2 2 3 3 ! routine which changes the albedo (and emissivity) of the surface … … 11 11 USE mod_phys_lmdz_transfert_para, ONLY: bcast 12 12 USE mod_phys_lmdz_para, ONLY: is_master 13 USE paleoclimate_mod, ONLY: paleoclimate,albedo_perenialco2 14 13 15 implicit none 14 16 … … 19 21 integer,intent(in) :: ngrid 20 22 real,intent(in) :: piceco2(ngrid) ! amount of CO2 ice on the surface (kg/m2) 23 real,intent(inout) :: piceco2_peren(ngrid) ! amount of perenial co2 ice (kg/m^2) 21 24 real,intent(out) :: psolaralb(ngrid,2) ! albedo of the surface 22 25 real,intent(out) :: emisref(ngrid) ! emissivity of the surface … … 104 107 psolaralb(ig,1)=albedice(icap) 105 108 psolaralb(ig,2)=albedice(icap) 109 if(paleoclimate) then 110 if((piceco2_peren(ig).gt.0.).and.(piceco2(ig).lt.piceco2_peren(ig))) then 111 psolaralb(ig,1) = albedo_perenialco2 112 psolaralb(ig,2) = albedo_perenialco2 113 piceco2_peren(ig) = piceco2(ig) 114 endif 115 endif 106 116 endif 107 117 else if (watercaptag(ig) .and. water) then -
trunk/LMDZ.MARS/libf/phymars/co2condens_mod.F
r2977 r2999 11 11 $ pcapcal,pplay,pplev,ptsrf,pt, 12 12 $ pphi,pdt,pdu,pdv,pdtsrf,pu,pv,pq,pdq, 13 $ piceco2,psolaralb,pemisurf,rdust, 13 $ piceco2,perenial_co2ice, 14 $ psolaralb,pemisurf,rdust, 14 15 $ pdtc,pdtsrfc,pdpsrf,pduc,pdvc,pdqc, 15 16 $ fluxsurf_sw,zls, … … 35 36 #endif 36 37 use comslope_mod, ONLY: subslope_dist,def_slope_mean 38 USE paleoclimate_mod, ONLY: paleoclimate 39 37 40 IMPLICIT NONE 38 41 c======================================================================= … … 87 90 88 91 REAL,INTENT(INOUT) :: piceco2(ngrid,nslope) ! CO2 ice on the surface (kg.m-2) 92 REAL,INTENT(INOUT) :: perenial_co2ice(ngrid,nslope) ! Perenial CO2 ice on the surface (kg.m-2) 89 93 REAL,INTENT(INOUT) :: psolaralb(ngrid,2,nslope) ! albedo of the surface 90 94 REAL,INTENT(INOUT) :: pemisurf(ngrid,nslope) ! emissivity of the surface … … 180 184 REAL :: alb_tmp(ngrid,2) ! local 181 185 REAL :: zcondices_tmp(ngrid) ! local 182 REAL :: piceco2_tmp(ngrid) ! local 186 REAL :: piceco2_tmp(ngrid) ! local 187 REAL :: perenial_co2ice_tmp(ngrid) ! perenial ice on one subslope (kg/m^2) 183 188 REAL :: pemisurf_tmp(ngrid)! local 184 189 LOGICAL :: condsub_tmp(ngrid) !local … … 568 573 piceco2_tmp(:) = piceco2(:,islope) 569 574 alb_tmp(:,:) = psolaralb(:,:,islope) 570 emisref_tmp(:) = 0. 571 CALL albedocaps(zls,ngrid,piceco2_tmp,alb_tmp,emisref_tmp) 575 emisref_tmp(:) = 0. 576 perenial_co2ice_tmp(:) = perenial_co2ice(:,islope) 577 CALL albedocaps(zls,ngrid,piceco2_tmp,perenial_co2ice_tmp, 578 & alb_tmp,emisref_tmp) 579 perenial_co2ice(:,islope) = perenial_co2ice_tmp(:) 572 580 psolaralb(:,1,islope) = alb_tmp(:,1) 573 581 psolaralb(:,2,islope) = alb_tmp(:,2) … … 823 831 ! Extra special case for surface temperature tendency pdtsrfc: 824 832 ! we want to fix the south pole temperature to CO2 condensation temperature 825 if (caps.and.(obliquit.lt.27.) ) then833 if (caps.and.(obliquit.lt.27.).and.(.not.(paleoclimate))) then 826 834 ! check if last grid point is the south pole 827 835 if (abs(latitude(ngrid)-(-pi/2.)).lt.1.e-5) then -
trunk/LMDZ.MARS/libf/phymars/conf_phys.F
r2994 r2999 54 54 use aeropacity_mod, only: iddist, topdustref 55 55 USE mod_phys_lmdz_transfert_para, ONLY: bcast 56 USE paleoclimate_mod,ONLY: paleoclimate 56 USE paleoclimate_mod,ONLY: paleoclimate,albedo_perenialco2 57 57 IMPLICIT NONE 58 58 include "callkeys.h" … … 355 355 write(*,*)" paleoclimate = ",paleoclimate 356 356 357 write(*,*)"Albedo for perenial CO2 ice?" 358 albedo_perenialco2 = 0.85 ! default value 359 call getin_p("albedo_perenialco2",albedo_perenialco2) 360 write(*,*)"albedo_perenialco2 = ",albedo_perenialco2 357 361 358 362 ! TRACERS: -
trunk/LMDZ.MARS/libf/phymars/paleoclimate_mod.F90
r2994 r2999 11 11 IMPLICIT NONE 12 12 13 LOGICAL :: paleoclimate ! False by default, is activate for paleoclimates specific processes (e.g., lag layer)13 LOGICAL :: paleoclimate ! False by default, is activate for paleoclimates specific processes (e.g., lag layer) 14 14 real, save, allocatable :: lag_h2o_ice(:,:) ! Thickness of the lag before H2O ice [m] 15 15 real, save, allocatable :: lag_co2_ice(:,:) ! Thickness of the lag before CO2 ice [m] 16 16 real, save :: albedo_perenialco2 ! Albedo for perenial co2 ice [1] 17 17 CONTAINS 18 18 -
trunk/LMDZ.MARS/libf/phymars/phyetat0_mod.F90
r2994 r2999 10 10 subroutine phyetat0 (fichnom,tab0,Lmodif,nsoil,ngrid,nlay,nq, & 11 11 day_ini,time0,tsurf,tsoil,albedo,emis,q2,qsurf, & 12 tauscaling,totcloudfrac,wstar,watercap, def_slope, &13 def_slope _mean,subslope_dist)12 tauscaling,totcloudfrac,wstar,watercap,perenial_co2ice, & 13 def_slope,def_slope_mean,subslope_dist) 14 14 15 15 use tracer_mod, only: noms ! tracer names … … 28 28 USE comslope_mod, ONLY: nslope, major_slope 29 29 USE paleoclimate_mod, ONLY: paleoclimate, lag_h2o_ice, lag_co2_ice 30 USE comcstfi_h, only: pi 31 use geometry_mod, only: latitude 30 32 implicit none 31 33 … … 71 73 real,intent(out) :: wstar(ngrid) ! Max vertical velocity in thermals (m/s) 72 74 real,intent(out) :: watercap(ngrid,nslope) ! h2o_ice_cover 73 real, intent(out) :: def_slope(nslope+1) !boundaries for bining of the slopes 74 real, intent(out) :: def_slope_mean(nslope) 75 real, intent(out) :: subslope_dist(ngrid,nslope) !undermesh statistics 75 real,intent(out) :: perenial_co2ice(ngrid,nslope) ! perenial co2 ice(kg/m^2) 76 real,intent(out) :: def_slope(nslope+1) !boundaries for bining of the slopes 77 real,intent(out) :: def_slope_mean(nslope) 78 real,intent(out) :: subslope_dist(ngrid,nslope) !undermesh statistics 76 79 !====================================================================== 77 80 ! Local variables: … … 117 120 REAL :: sum_dist 118 121 REAL :: current_max !var to find max distrib slope 122 ! Variables for CO2 index 123 INTEGER :: igcm_co2_tmp 119 124 120 125 CHARACTER(len=5) :: modname="phyetat0" … … 755 760 756 761 if (paleoclimate) then 762 do iq=1,nq 763 txt=noms(iq) 764 if (txt.eq."co2") igcm_co2_tmp = iq 765 enddo 757 766 if (startphy_file) then 767 ! Depth of H2O lag 758 768 call get_field("lag_h2o_ice",lag_h2o_ice,found,indextime) 759 write(*,*) 'paleo found start?',found760 761 769 if (.not.found) then 762 770 write(*,*) "phyetat0: Failed loading <lag_h2o_ice> : ", & … … 765 773 endif 766 774 775 ! Depth of CO2 lag 767 776 call get_field("lag_co2_ice",lag_co2_ice,found,indextime) 768 777 if (.not.found) then … … 771 780 lag_co2_ice(:,:) = -1. 772 781 endif 773 else 782 783 ! Perenial CO2 ice 784 call get_field("perenial_co2ice",perenial_co2ice,found,indextime) 785 if (.not.found) then 786 write(*,*) "phyetat0: Failed loading <perenial_co2ice> : ", & 787 "<perenial_co2ice> is set as 10m at the South Pole" 788 perenial_co2ice(:,:) = 0. 789 if (abs(latitude(ngrid)-(-pi/2.)).lt.1.e-5) then 790 DO islope = 1,nslope 791 perenial_co2ice(ngrid,islope) = 10*1.6e3 ! 10m which is convert to kg/m^2 792 qsurf(ngrid,igcm_co2_tmp,islope) = qsurf(ngrid-1,igcm_co2_tmp,islope) + & ! perenial ice + frost 793 perenial_co2ice(ngrid,islope) 794 ENDDO 795 endif 796 endif ! notfound 797 else ! no startfiphyle 774 798 lag_h2o_ice(:,:) = -1. 775 799 lag_co2_ice(:,:) = -1. 800 perenial_co2ice(:,:) = 0. 801 if (abs(latitude(ngrid)-(-pi/2.)).lt.1.e-5) then 802 DO islope = 1,nslope 803 perenial_co2ice(ngrid,islope) = 10*1.6e3 ! 10m which is convert to kg/m^2 804 qsurf(ngrid,igcm_co2_tmp,islope) = qsurf(ngrid-1,igcm_co2_tmp,islope) + & ! perenial ice + frost 805 perenial_co2ice(ngrid,islope) 806 ENDDO 807 endif 776 808 endif !startphy_file 777 809 else 778 write(*,*) 'paleo found? nostart',found779 780 810 lag_h2o_ice(:,:) = -1. 781 811 lag_co2_ice(:,:) = -1. 812 perenial_co2ice(:,:) = 0. 782 813 endif !paleoclimate 783 814 -
trunk/LMDZ.MARS/libf/phymars/phyredem.F90
r2994 r2999 168 168 albedo,emis,q2,qsurf,& 169 169 tauscaling,totcloudfrac,wstar, & 170 watercap )170 watercap,perenial_co2ice) 171 171 ! write time-dependent variable to restart file 172 172 use iostart, only : open_restartphy, close_restartphy, & … … 178 178 use dust_param_mod, only: dustscaling_mode 179 179 use comsoil_h,only: flux_geo 180 use comslope_mod, ONLY: nslope 180 use comslope_mod, only: nslope 181 use paleoclimate_mod, only: paleoclimate 181 182 implicit none 182 183 … … 201 202 real,intent(in) :: wstar(ngrid) 202 203 real,intent(in) :: watercap(ngrid,nslope) 203 204 real,intent(in) :: perenial_co2ice(ngrid,nslope) 205 204 206 integer :: iq 205 207 character(len=30) :: txt ! to store some text … … 220 222 ! Water ice layer 221 223 call put_field("watercap","H2O ice cover",watercap,time) 222 224 225 ! Perenial CO2 ice layer 226 if(paleoclimate) then 227 call put_field("perenial_co2ice","CO2 ice cover",perenial_co2ice,time) 228 endif 223 229 ! Surface temperature 224 230 call put_field("tsurf","Surface temperature",tsurf,time) -
trunk/LMDZ.MARS/libf/phymars/physiq_mod.F
r2965 r2999 52 52 & tsurf, emis, 53 53 & capcal, fluxgrd, qsurf, 54 & hmons,summit,base,watercap,watercaptag 54 & hmons,summit,base,watercap,watercaptag, 55 & perenial_co2ice 55 56 use comsaison_h, only: dist_sol, declin, zls, 56 57 & mu0, fract, local_time … … 595 596 & tsurf,tsoil,albedo,emis, 596 597 & q2,qsurf,tauscaling,totcloudfrac,wstar, 597 & watercap,def_slope,def_slope_mean,subslope_dist) 598 & watercap,perenial_co2ice, 599 & def_slope,def_slope_mean,subslope_dist) 598 600 599 601 ! Sky view: … … 2199 2201 $ capcal,zplay,zplev,tsurf,pt, 2200 2202 $ pphi,pdt,pdu,pdv,zdtsurf,pu,pv,pq,pdq, 2201 $ qsurf(:,igcm_co2,:),albedo,emis,rdust, 2203 $ qsurf(:,igcm_co2,:),perenial_co2ice, 2204 $ albedo,emis,rdust, 2202 2205 $ zdtc,zdtsurfc,pdpsrf,zduc,zdvc,zdqc, 2203 2206 $ fluxsurf_dn_sw,zls, … … 2614 2617 . tsurf,tsoil,inertiesoil,albedo, 2615 2618 . emis,q2,qsurf,tauscaling,totcloudfrac,wstar, 2616 . watercap )2619 . watercap,perenial_co2ice) 2617 2620 ENDIF ! of IF (write_restart) 2618 2621 -
trunk/LMDZ.MARS/libf/phymars/surfdat_h.F90
r2909 r2999 50 50 REAL,ALLOCATABLE,SAVE :: qsurf(:,:,:) ! tracer on surface (e.g. kg.m-2) 51 51 REAL,SAVE,ALLOCATABLE :: watercap(:,:) ! Surface water ice (kg.m-2) 52 53 !$OMP THREADPRIVATE(tsurf,emis,capcal,fluxgrd,qsurf,watercap )52 REAL,SAVE,ALLOCATABLE :: perenial_co2ice(:,:) ! Perenial CO2 ice (kg.m-2) 53 !$OMP THREADPRIVATE(tsurf,emis,capcal,fluxgrd,qsurf,watercap,perenial_co2ice) 54 54 55 55 contains … … 74 74 allocate(tsurf(ngrid,nslope)) 75 75 allocate(watercap(ngrid,nslope)) 76 allocate(perenial_co2ice(ngrid,nslope)) 76 77 allocate(emis(ngrid,nslope)) 77 78 allocate(capcal(ngrid,nslope)) … … 91 92 implicit none 92 93 93 if (allocated(albedodat)) deallocate(albedodat) 94 if (allocated(phisfi)) deallocate(phisfi) 95 if (allocated(watercaptag)) deallocate(watercaptag) 96 if (allocated(dryness)) deallocate(dryness) 97 if (allocated(zmea)) deallocate(zmea) 98 if (allocated(zstd)) deallocate(zstd) 99 if (allocated(zsig)) deallocate(zsig) 100 if (allocated(zgam)) deallocate(zgam) 101 if (allocated(zthe)) deallocate(zthe) 102 if (allocated(z0)) deallocate(z0) 103 if (allocated(qsurf)) deallocate(qsurf) 104 if (allocated(tsurf)) deallocate(tsurf) 105 if (allocated(watercap)) deallocate(watercap) 106 if (allocated(emis)) deallocate(emis) 107 if (allocated(capcal)) deallocate(capcal) 108 if (allocated(fluxgrd)) deallocate(fluxgrd) 109 if (allocated(hmons)) deallocate(hmons) 110 if (allocated(summit)) deallocate(summit) 111 if (allocated(base)) deallocate(base) 112 if (allocated(alpha_hmons)) deallocate(alpha_hmons) 113 if (allocated(hsummit)) deallocate(hsummit) 114 if (allocated(contains_mons)) deallocate(contains_mons) 94 if (allocated(albedodat)) deallocate(albedodat) 95 if (allocated(phisfi)) deallocate(phisfi) 96 if (allocated(watercaptag)) deallocate(watercaptag) 97 if (allocated(dryness)) deallocate(dryness) 98 if (allocated(zmea)) deallocate(zmea) 99 if (allocated(zstd)) deallocate(zstd) 100 if (allocated(zsig)) deallocate(zsig) 101 if (allocated(zgam)) deallocate(zgam) 102 if (allocated(zthe)) deallocate(zthe) 103 if (allocated(z0)) deallocate(z0) 104 if (allocated(qsurf)) deallocate(qsurf) 105 if (allocated(tsurf)) deallocate(tsurf) 106 if (allocated(watercap)) deallocate(watercap) 107 if (allocated(perenial_co2ice)) deallocate(perenial_co2ice) 108 if (allocated(emis)) deallocate(emis) 109 if (allocated(capcal)) deallocate(capcal) 110 if (allocated(fluxgrd)) deallocate(fluxgrd) 111 if (allocated(hmons)) deallocate(hmons) 112 if (allocated(summit)) deallocate(summit) 113 if (allocated(base)) deallocate(base) 114 if (allocated(alpha_hmons)) deallocate(alpha_hmons) 115 if (allocated(hsummit)) deallocate(hsummit) 116 if (allocated(contains_mons)) deallocate(contains_mons) 115 117 116 118 end subroutine end_surfdat_h … … 125 127 allocate(tsurf(ngrid,nslope)) 126 128 allocate(watercap(ngrid,nslope)) 129 allocate(perenial_co2ice(ngrid,nslope)) 127 130 allocate(emis(ngrid,nslope)) 128 131 allocate(capcal(ngrid,nslope)) … … 135 138 implicit none 136 139 137 if (allocated(qsurf)) deallocate(qsurf) 138 if (allocated(tsurf)) deallocate(tsurf) 139 if (allocated(watercap)) deallocate(watercap) 140 if (allocated(emis)) deallocate(emis) 141 if (allocated(capcal)) deallocate(capcal) 142 if (allocated(fluxgrd)) deallocate(fluxgrd) 140 if (allocated(qsurf)) deallocate(qsurf) 141 if (allocated(tsurf)) deallocate(tsurf) 142 if (allocated(watercap)) deallocate(watercap) 143 if (allocated(perenial_co2ice)) deallocate(perenial_co2ice) 144 if (allocated(emis)) deallocate(emis) 145 if (allocated(capcal)) deallocate(capcal) 146 if (allocated(fluxgrd)) deallocate(fluxgrd) 143 147 144 148 end subroutine end_surfdat_h_slope_var
Note: See TracChangeset
for help on using the changeset viewer.