Ignore:
Timestamp:
Jul 19, 2023, 11:40:38 AM (2 years ago)
Author:
llange
Message:

Mars PCM
Include perenial_co2ice (equivalent of watercap) to distinguich between CO2 frost and perenial CO2 ice for paleoclimate studies.
When no frost is present and we dig into perenial ice, the surface albedo is changed. The albedo for seasonal ice is set to 0.65, and the perenial ice albedo can be fixed in the callphys.def. I recommand values between 0.8 and 0.9.
To use this, paleoclimate must be set to True and TESalbedo to false in the callphys.def. Else, it runs as usual with TES albedo
LL

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,psolaralb,emisref)
     1subroutine albedocaps(zls,ngrid,piceco2,piceco2_peren,psolaralb,emisref)
    22
    33! routine which changes the albedo (and emissivity) of the surface
     
    1111USE mod_phys_lmdz_transfert_para, ONLY: bcast
    1212USE mod_phys_lmdz_para, ONLY: is_master
     13USE paleoclimate_mod, ONLY: paleoclimate,albedo_perenialco2
     14
    1315implicit none
    1416
     
    1921integer,intent(in) :: ngrid
    2022real,intent(in) :: piceco2(ngrid) ! amount of CO2 ice on the surface (kg/m2)
     23real,intent(inout) :: piceco2_peren(ngrid) ! amount of perenial co2 ice (kg/m^2)
    2124real,intent(out) :: psolaralb(ngrid,2) ! albedo of the surface
    2225real,intent(out) :: emisref(ngrid) ! emissivity of the surface
     
    104107      psolaralb(ig,1)=albedice(icap)
    105108      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
    106116    endif
    107117  else if (watercaptag(ig) .and. water) then
  • trunk/LMDZ.MARS/libf/phymars/co2condens_mod.F

    r2977 r2999  
    1111     $                  pcapcal,pplay,pplev,ptsrf,pt,
    1212     $                  pphi,pdt,pdu,pdv,pdtsrf,pu,pv,pq,pdq,
    13      $                  piceco2,psolaralb,pemisurf,rdust,
     13     $                  piceco2,perenial_co2ice,
     14     $                  psolaralb,pemisurf,rdust,
    1415     $                  pdtc,pdtsrfc,pdpsrf,pduc,pdvc,pdqc,
    1516     $                  fluxsurf_sw,zls,
     
    3536#endif
    3637      use comslope_mod, ONLY: subslope_dist,def_slope_mean
     38      USE paleoclimate_mod, ONLY: paleoclimate
     39
    3740       IMPLICIT NONE
    3841c=======================================================================
     
    8790
    8891      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)
    8993      REAL,INTENT(INOUT) :: psolaralb(ngrid,2,nslope) ! albedo of the surface
    9094      REAL,INTENT(INOUT) :: pemisurf(ngrid,nslope) ! emissivity of the surface
     
    180184      REAL   :: alb_tmp(ngrid,2) ! local
    181185      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)
    183188      REAL   :: pemisurf_tmp(ngrid)! local
    184189      LOGICAL :: condsub_tmp(ngrid) !local
     
    568573        piceco2_tmp(:) = piceco2(:,islope)
    569574        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(:)
    572580        psolaralb(:,1,islope) =  alb_tmp(:,1)
    573581        psolaralb(:,2,islope) =  alb_tmp(:,2)
     
    823831! Extra special case for surface temperature tendency pdtsrfc:
    824832! we want to fix the south pole temperature to CO2 condensation temperature
    825          if (caps.and.(obliquit.lt.27.)) then
     833         if (caps.and.(obliquit.lt.27.).and.(.not.(paleoclimate))) then
    826834           ! check if last grid point is the south pole
    827835           if (abs(latitude(ngrid)-(-pi/2.)).lt.1.e-5) then
  • trunk/LMDZ.MARS/libf/phymars/conf_phys.F

    r2994 r2999  
    5454      use aeropacity_mod, only: iddist, topdustref
    5555      USE mod_phys_lmdz_transfert_para, ONLY: bcast
    56       USE paleoclimate_mod,ONLY: paleoclimate
     56      USE paleoclimate_mod,ONLY: paleoclimate,albedo_perenialco2
    5757      IMPLICIT NONE
    5858      include "callkeys.h"
     
    355355         write(*,*)" paleoclimate = ",paleoclimate
    356356
     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
    357361
    358362! TRACERS:
  • trunk/LMDZ.MARS/libf/phymars/paleoclimate_mod.F90

    r2994 r2999  
    1111    IMPLICIT NONE
    1212
    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)
    1414    real, save, allocatable :: lag_h2o_ice(:,:)  ! Thickness of the lag before H2O ice [m]
    1515    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]
    1717    CONTAINS
    1818
  • trunk/LMDZ.MARS/libf/phymars/phyetat0_mod.F90

    r2994 r2999  
    1010subroutine phyetat0 (fichnom,tab0,Lmodif,nsoil,ngrid,nlay,nq, &
    1111                     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)
    1414
    1515  use tracer_mod, only: noms ! tracer names
     
    2828  USE comslope_mod, ONLY: nslope, major_slope
    2929  USE paleoclimate_mod, ONLY: paleoclimate, lag_h2o_ice, lag_co2_ice
     30  USE comcstfi_h, only: pi
     31  use geometry_mod, only: latitude
    3032  implicit none
    3133 
     
    7173  real,intent(out) :: wstar(ngrid) ! Max vertical velocity in thermals (m/s)
    7274  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
    7679!======================================================================
    7780!  Local variables:
     
    117120      REAL :: sum_dist
    118121      REAL :: current_max !var to find max distrib slope
     122! Variables for CO2 index
     123      INTEGER :: igcm_co2_tmp
    119124
    120125      CHARACTER(len=5) :: modname="phyetat0"
     
    755760
    756761if (paleoclimate) then
     762  do iq=1,nq
     763   txt=noms(iq)
     764   if (txt.eq."co2") igcm_co2_tmp = iq
     765  enddo
    757766  if (startphy_file) then
     767! Depth of H2O lag
    758768   call get_field("lag_h2o_ice",lag_h2o_ice,found,indextime)
    759    write(*,*) 'paleo found start?',found
    760    
    761769   if (.not.found) then
    762770     write(*,*) "phyetat0: Failed loading <lag_h2o_ice> : ", &
     
    765773   endif
    766774
     775!  Depth of CO2 lag
    767776   call get_field("lag_co2_ice",lag_co2_ice,found,indextime)
    768777   if (.not.found) then
     
    771780     lag_co2_ice(:,:) = -1.
    772781   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
    774798     lag_h2o_ice(:,:) = -1.
    775799     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
    776808  endif !startphy_file
    777809else
    778    write(*,*) 'paleo found? nostart',found
    779    
    780810  lag_h2o_ice(:,:) = -1.
    781811  lag_co2_ice(:,:) = -1.
     812  perenial_co2ice(:,:) = 0.
    782813endif !paleoclimate
    783814
  • trunk/LMDZ.MARS/libf/phymars/phyredem.F90

    r2994 r2999  
    168168                    albedo,emis,q2,qsurf,&
    169169                    tauscaling,totcloudfrac,wstar, &
    170                     watercap)
     170                    watercap,perenial_co2ice)
    171171  ! write time-dependent variable to restart file
    172172  use iostart, only : open_restartphy, close_restartphy, &
     
    178178  use dust_param_mod, only: dustscaling_mode
    179179  use comsoil_h,only: flux_geo
    180   use comslope_mod, ONLY: nslope
     180  use comslope_mod, only: nslope
     181  use paleoclimate_mod, only: paleoclimate
    181182  implicit none
    182183 
     
    201202  real,intent(in) :: wstar(ngrid)
    202203  real,intent(in) :: watercap(ngrid,nslope)
    203  
     204  real,intent(in) :: perenial_co2ice(ngrid,nslope)
     205
    204206  integer :: iq
    205207  character(len=30) :: txt ! to store some text
     
    220222  ! Water ice layer
    221223  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
    223229  ! Surface temperature
    224230  call put_field("tsurf","Surface temperature",tsurf,time)
  • trunk/LMDZ.MARS/libf/phymars/physiq_mod.F

    r2965 r2999  
    5252     &                     tsurf, emis,
    5353     &                     capcal, fluxgrd, qsurf,
    54      &                     hmons,summit,base,watercap,watercaptag
     54     &                     hmons,summit,base,watercap,watercaptag,
     55     &                     perenial_co2ice
    5556      use comsaison_h, only: dist_sol, declin, zls,
    5657     &                       mu0, fract, local_time
     
    595596     &         tsurf,tsoil,albedo,emis,
    596597     &         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)
    598600
    599601!     Sky view:
     
    21992201     $              capcal,zplay,zplev,tsurf,pt,
    22002202     $              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,
    22022205     $              zdtc,zdtsurfc,pdpsrf,zduc,zdvc,zdqc,
    22032206     $              fluxsurf_dn_sw,zls,
     
    26142617     .                tsurf,tsoil,inertiesoil,albedo,
    26152618     .                emis,q2,qsurf,tauscaling,totcloudfrac,wstar,
    2616      .                watercap)
     2619     .                watercap,perenial_co2ice)
    26172620          ENDIF ! of IF (write_restart)
    26182621
  • trunk/LMDZ.MARS/libf/phymars/surfdat_h.F90

    r2909 r2999  
    5050  REAL,ALLOCATABLE,SAVE :: qsurf(:,:,:) ! tracer on surface (e.g. kg.m-2)
    5151  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)
    5454
    5555contains
     
    7474    allocate(tsurf(ngrid,nslope))
    7575    allocate(watercap(ngrid,nslope))
     76    allocate(perenial_co2ice(ngrid,nslope))
    7677    allocate(emis(ngrid,nslope))
    7778    allocate(capcal(ngrid,nslope))
     
    9192  implicit none
    9293
    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)
    115117   
    116118  end subroutine end_surfdat_h
     
    125127    allocate(tsurf(ngrid,nslope))
    126128    allocate(watercap(ngrid,nslope))
     129    allocate(perenial_co2ice(ngrid,nslope))
    127130    allocate(emis(ngrid,nslope))
    128131    allocate(capcal(ngrid,nslope))
     
    135138  implicit none
    136139
    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)
    143147   
    144148  end subroutine end_surfdat_h_slope_var
Note: See TracChangeset for help on using the changeset viewer.