Ignore:
Timestamp:
Nov 20, 2023, 1:25:31 PM (20 months ago)
Author:
jbclement
Message:

PEM:
The perennial co2 ice is now taken into account with co2 frost (qsurf) to compute the tendency and to make the update + Rework of how co2 frost is converted to perennial co2 ice at the end of the PEM run + Correction of the value of 'threshold_co2_frost2perennial' to correspond to 10 m + Perennial co2 ice is now handled outside 'paleoclimate' in "phyetat0_mod.F90" of the Mars PCM + Some cleanings.

/!\ Commit for the PEM management of co2 ice before a rework of ice management in the PEM!
JBC

Location:
trunk/LMDZ.MARS/libf/phymars
Files:
9 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.MARS/libf/phymars/albedocaps.F90

    r2999 r3130  
    1111USE mod_phys_lmdz_transfert_para, ONLY: bcast
    1212USE mod_phys_lmdz_para, ONLY: is_master
    13 USE paleoclimate_mod, ONLY: paleoclimate,albedo_perenialco2
     13USE paleoclimate_mod, ONLY: paleoclimate,albedo_perennialco2
    1414
    1515implicit none
     
    2121integer,intent(in) :: ngrid
    2222real,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)
     23real,intent(inout) :: piceco2_peren(ngrid) ! amount of perennial co2 ice (kg/m^2)
    2424real,intent(out) :: psolaralb(ngrid,2) ! albedo of the surface
    2525real,intent(out) :: emisref(ngrid) ! emissivity of the surface
     
    109109      if(paleoclimate) then
    110110         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
     111             psolaralb(ig,1) = albedo_perennialco2
     112             psolaralb(ig,2) = albedo_perennialco2
    113113             piceco2_peren(ig) = piceco2(ig)
    114114         endif
  • trunk/LMDZ.MARS/libf/phymars/co2condens_mod.F

    r3091 r3130  
    1111     $                  pcapcal,pplay,pplev,ptsrf,pt,
    1212     $                  pphi,pdt,pdu,pdv,pdtsrf,pu,pv,pq,pdq,
    13      $                  piceco2,perenial_co2ice,
     13     $                  piceco2,perennial_co2ice,
    1414     $                  psolaralb,pemisurf,rdust,
    1515     $                  pdtc,pdtsrfc,pdpsrf,pduc,pdvc,pdqc,
     
    9090
    9191      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)
     92      REAL,INTENT(INOUT) :: perennial_co2ice(ngrid,nslope) ! perennial CO2 ice on the surface (kg.m-2)
    9393      REAL,INTENT(INOUT) :: psolaralb(ngrid,2,nslope) ! albedo of the surface
    9494      REAL,INTENT(INOUT) :: pemisurf(ngrid,nslope) ! emissivity of the surface
     
    186186      REAL   :: zcondices_tmp(ngrid)          ! temporary condensation rate [kg/m^2/s]
    187187      REAL   :: piceco2_tmp(ngrid)            ! temporary amount of CO2 frost on the surface [kg/m^2]
    188       REAL   :: perenial_co2ice_tmp(ngrid)    ! temporary amount of perenial CO2 frost on the surface [kg/m^2]
     188      REAL   :: perennial_co2ice_tmp(ngrid)    ! temporary amount of perennial CO2 frost on the surface [kg/m^2]
    189189      LOGICAL :: condsub_tmp(ngrid)           ! Boolean to check if CO2 ice is condensing / sublimating on the sub grid surface [1]
    190190      REAL :: zfallice_tmp(ngrid,nlayer+1)    ! temporary amount of ice falling from layer l for a specific sub-grid surface [kg/m^2/s]
     
    581581        alb_tmp(:,:) = psolaralb(:,:,islope)
    582582        emisref_tmp(:) = 0.
    583         perenial_co2ice_tmp(:) =  perenial_co2ice(:,islope)
    584         CALL albedocaps(zls,ngrid,piceco2_tmp,perenial_co2ice_tmp,
     583        perennial_co2ice_tmp(:) =  perennial_co2ice(:,islope)
     584        CALL albedocaps(zls,ngrid,piceco2_tmp,perennial_co2ice_tmp,
    585585     &                  alb_tmp,emisref_tmp)
    586         perenial_co2ice(:,islope) = perenial_co2ice_tmp(:)
     586        perennial_co2ice(:,islope) = perennial_co2ice_tmp(:)
    587587        psolaralb(:,1,islope) =  alb_tmp(:,1)
    588588        psolaralb(:,2,islope) =  alb_tmp(:,2)
  • trunk/LMDZ.MARS/libf/phymars/conf_phys.F

    r3128 r3130  
    3737      use aeropacity_mod, only: iddist, topdustref
    3838      USE mod_phys_lmdz_transfert_para, ONLY: bcast
    39       USE paleoclimate_mod,ONLY: paleoclimate,albedo_perenialco2,
     39      USE paleoclimate_mod,ONLY: paleoclimate,albedo_perennialco2,
    4040     &                           lag_layer
    4141      use microphys_h, only: mteta
     
    351351
    352352
    353          write(*,*)"Albedo for perenial CO2 ice?"
    354          albedo_perenialco2 = 0.85 ! default value
    355          call getin_p("albedo_perenialco2",albedo_perenialco2)
    356          write(*,*)"albedo_perenialco2 = ",albedo_perenialco2
     353         write(*,*)"Albedo for perennial CO2 ice?"
     354         albedo_perennialco2 = 0.85 ! default value
     355         call getin_p("albedo_perennialco2",albedo_perennialco2)
     356         write(*,*)"albedo_perennialco2 = ",albedo_perennialco2
    357357
    358358! TRACERS:
  • trunk/LMDZ.MARS/libf/phymars/dyn1d/testphys1d.F90

    r3129 r3130  
    22
    33use comsoil_h,           only: inertiedat, inertiesoil, nsoilmx, tsoil, nqsoil, qsoil
    4 use surfdat_h,           only: albedodat, perenial_co2ice, watercap, tsurf, emis, qsurf
     4use surfdat_h,           only: albedodat, perennial_co2ice, watercap, tsurf, emis, qsurf
    55use comslope_mod,        only: def_slope, subslope_dist
    66use phyredem,            only: physdem0, physdem1
     
    161161    call physdem1("startfi.nc",nsoilmx,ngrid,llm,nq,nqsoil,dttestphys,time,      &
    162162                  tsurf,tsoil,inertiesoil,albedo,emis,q2,qsurf,qsoil,tauscaling, &
    163                   totcloudfrac,wstar,watercap,perenial_co2ice)
     163                  totcloudfrac,wstar,watercap,perennial_co2ice)
    164164endif !(.not. therestartfi)
    165165
  • trunk/LMDZ.MARS/libf/phymars/paleoclimate_mod.F90

    r3111 r3130  
    1717    real, save, allocatable :: h2o_ice_depth(:,:)  ! Thickness of the lag before H2O ice [m]
    1818    real, save, allocatable :: lag_co2_ice(:,:)  ! Thickness of the lag before CO2 ice [m]
    19     real, save :: albedo_perenialco2             ! Albedo for perenial co2 ice [1]
     19    real, save :: albedo_perennialco2             ! Albedo for perennial co2 ice [1]
    2020    real, save, allocatable :: d_coef(:,:)  ! Diffusion coeficent
    2121    LOGICAL,SAVE :: lag_layer ! does lag layer is present?
    22 !$OMP THREADPRIVATE(h2o_ice_depth,d_coef,lag_co2_ice,albedo_perenialco2)
     22!$OMP THREADPRIVATE(h2o_ice_depth,d_coef,lag_co2_ice,albedo_perennialco2)
    2323
    2424    CONTAINS
  • trunk/LMDZ.MARS/libf/phymars/phyetat0_mod.F90

    r3113 r3130  
    1010subroutine phyetat0 (fichnom,tab0,Lmodif,nsoil,ngrid,nlay,nq,nqsoil, &
    1111                     day_ini,time0,tsurf,tsoil,albedo,emis,q2,qsurf,qsoil, &
    12                      tauscaling,totcloudfrac,wstar,watercap,perenial_co2ice, &
     12                     tauscaling,totcloudfrac,wstar,watercap,perennial_co2ice, &
    1313                     def_slope,def_slope_mean,subslope_dist)
    1414
     
    7575  real,intent(out) :: wstar(ngrid) ! Max vertical velocity in thermals (m/s)
    7676  real,intent(out) :: watercap(ngrid,nslope) ! h2o_ice_cover
    77   real,intent(out) :: perenial_co2ice(ngrid,nslope) ! perenial co2 ice(kg/m^2)
     77  real,intent(out) :: perennial_co2ice(ngrid,nslope) ! perennial co2 ice(kg/m^2)
    7878  real,intent(out) :: def_slope(nslope+1) !boundaries for bining of the slopes
    7979  real,intent(out) :: def_slope_mean(nslope)
     
    753753   endif
    754754
    755 ! Diffuesion coeficent
     755! Diffusion coeficent
    756756   call get_field("d_coef",d_coef,found,indextime)
    757757   if (.not.found) then
     
    761761   endif
    762762
    763 
    764 !  Depth of CO2 lag
     763! Depth of CO2 lag
    765764   call get_field("lag_co2_ice",lag_co2_ice,found,indextime)
    766765   if (.not.found) then
     
    769768     lag_co2_ice(:,:) = -1.
    770769   endif
    771 
    772 ! Perenial CO2 ice   
    773    call get_field("perenial_co2ice",perenial_co2ice,found,indextime)
    774    if (.not.found) then
    775      write(*,*) "phyetat0: Failed loading <perenial_co2ice> : ", &
    776                 "<perenial_co2ice> is set as 10m at the South Pole"
    777      perenial_co2ice(:,:) = 0.
    778      if (abs(latitude(ngrid)-(-pi/2.)).lt.1.e-5) then
    779        DO islope = 1,nslope
    780         perenial_co2ice(ngrid,islope) = 10*1.6e3 ! 10m which is convert to kg/m^2
    781         qsurf(ngrid,igcm_co2_tmp,islope) = qsurf(ngrid-1,igcm_co2_tmp,islope) + & ! perenial ice + frost
    782                                         perenial_co2ice(ngrid,islope)
    783        ENDDO
    784      endif
    785    endif ! notfound
    786770  else ! no startfiphyle
    787771     h2o_ice_depth(:,:) = -1.
    788772     lag_co2_ice(:,:) = -1.
    789      perenial_co2ice(:,:) = 0.
    790773     d_coef(:,:)= 4.e-4
    791      if (abs(latitude(ngrid)-(-pi/2.)).lt.1.e-5) then
    792        DO islope = 1,nslope
    793         perenial_co2ice(ngrid,islope) = 10*1.6e3 ! 10m which is convert to kg/m^2
    794         qsurf(ngrid,igcm_co2_tmp,islope) = qsurf(ngrid-1,igcm_co2_tmp,islope) + & ! perenial ice + frost
    795                                         perenial_co2ice(ngrid,islope)
    796        ENDDO
    797      endif
    798774  endif !startphy_file
    799775else
    800776  h2o_ice_depth(:,:) = -1.
    801777  lag_co2_ice(:,:) = -1.
    802   perenial_co2ice(:,:) = 0.
    803778  d_coef(:,:)= 4.e-4
    804779endif !paleoclimate
     780
     781! Perennial CO2 ice
     782perennial_co2ice(:,:) = 0.
     783if (startphy_file) then
     784    call get_field("perennial_co2ice",perennial_co2ice,found,indextime)
     785    if (.not. found) then
     786        write(*,*) "phyetat0: Failed loading <perennial_co2ice> : ", &
     787                   "<perennial_co2ice> is set as 10m at the South Pole"
     788        if (abs(latitude(ngrid) - (-pi/2.)) < 1.e-5) then
     789            do islope = 1,nslope
     790                perennial_co2ice(ngrid,islope) = 10*1.6e3 ! 10m which is convert to kg/m^2
     791               qsurf(ngrid,igcm_co2_tmp,islope) = qsurf(ngrid - 1,igcm_co2_tmp,islope) + perennial_co2ice(ngrid,islope) ! perennial ice + frost
     792            enddo
     793        endif
     794    endif ! not found
     795else
     796    if (abs(latitude(ngrid) - (-pi/2.)) < 1.e-5) then
     797        do islope = 1,nslope
     798            perennial_co2ice(ngrid,islope) = 10*1.6e3 ! 10m which is convert to kg/m^2
     799            qsurf(ngrid,igcm_co2_tmp,islope) = qsurf(ngrid - 1,igcm_co2_tmp,islope) + perennial_co2ice(ngrid,islope) ! perennial ice + frost
     800        enddo
     801    endif
     802endif !startphy_file
    805803
    806804! close file:
  • trunk/LMDZ.MARS/libf/phymars/phyredem.F90

    r3113 r3130  
    168168                    albedo,emis,q2,qsurf,qsoil,&
    169169                    tauscaling,totcloudfrac,wstar, &
    170                     watercap,perenial_co2ice)
     170                    watercap,perennial_co2ice)
    171171  ! write time-dependent variable to restart file
    172172  use iostart, only : open_restartphy, close_restartphy, &
     
    205205  real,intent(in) :: wstar(ngrid)
    206206  real,intent(in) :: watercap(ngrid,nslope)
    207   real,intent(in) :: perenial_co2ice(ngrid,nslope)
     207  real,intent(in) :: perennial_co2ice(ngrid,nslope)
    208208
    209209  integer :: iq
     
    226226  call put_field("watercap","H2O ice cover",watercap,time)
    227227 
    228  ! Perenial CO2 ice layer
    229   if (paleoclimate) call put_field("perenial_co2ice","CO2 ice cover",perenial_co2ice,time)
     228  ! Perennial CO2 ice layer
     229  call put_field("perennial_co2ice","CO2 ice cover",perennial_co2ice,time)
    230230
    231231  ! Surface temperature
  • trunk/LMDZ.MARS/libf/phymars/physiq_mod.F

    r3122 r3130  
    5959     &                     capcal, fluxgrd, qsurf,
    6060     &                     hmons,summit,base,watercap,watercaptag,
    61      &                     perenial_co2ice
     61     &                     perennial_co2ice
    6262      use comsaison_h, only: dist_sol, declin, zls,
    6363     &                       mu0, fract, local_time
     
    604604     &         tsurf,tsoil,albedo,emis,
    605605     &         q2,qsurf,qsoil,tauscaling,totcloudfrac,wstar,
    606      &         watercap,perenial_co2ice,
     606     &         watercap,perennial_co2ice,
    607607     &         def_slope,def_slope_mean,subslope_dist)
    608608
     
    22092209     $              capcal,zplay,zplev,tsurf,pt,
    22102210     $              pphi,pdt,pdu,pdv,zdtsurf,pu,pv,pq,pdq,
    2211      $              qsurf(:,igcm_co2,:),perenial_co2ice,
     2211     $              qsurf(:,igcm_co2,:),perennial_co2ice,
    22122212     $              albedo,emis,rdust,
    22132213     $              zdtc,zdtsurfc,pdpsrf,zduc,zdvc,zdqc,
     
    26252625     .                tsurf,tsoil,inertiesoil,albedo,
    26262626     .                emis,q2,qsurf,qsoil,tauscaling,totcloudfrac,
    2627      .                wstar,watercap,perenial_co2ice)
     2627     .                wstar,watercap,perennial_co2ice)
    26282628          ENDIF ! of IF (write_restart)
    26292629
     
    32053205     &         "Perennial water ice thickness"
    32063206     &         ,"kg.m-2",watercap(:,islope))
     3207         ENDDO
     3208         call write_output("perennial_co2ice",
     3209     &         "Perennial co2 ice thickness","kg.m-2",
     3210     &         perennial_co2ice(:,iflat))
     3211         do islope=1,nslope
     3212          write(str2(1:2),'(i2.2)') islope
     3213         call write_output("perennial_co2ice_slope"//str2,
     3214     &         "Perennial c ice thickness"
     3215     &         ,"kg.m-2",perennial_co2ice(:,islope))
    32073216         ENDDO
    32083217         call write_output("temp_layer1","temperature in layer 1",
  • trunk/LMDZ.MARS/libf/phymars/surfdat_h.F90

    r3111 r3130  
    5353  REAL,ALLOCATABLE,SAVE :: qsurf(:,:,:) ! tracer on surface (e.g. kg.m-2)
    5454  REAL,SAVE,ALLOCATABLE :: watercap(:,:) ! Surface water ice (kg.m-2)
    55   REAL,SAVE,ALLOCATABLE :: perenial_co2ice(:,:) ! Perenial CO2 ice (kg.m-2)
    56 !$OMP THREADPRIVATE(tsurf,emis,capcal,fluxgrd,qsurf,watercap,perenial_co2ice)
     55  REAL,SAVE,ALLOCATABLE :: perennial_co2ice(:,:) ! Perennial CO2 ice (kg.m-2)
     56!$OMP THREADPRIVATE(tsurf,emis,capcal,fluxgrd,qsurf,watercap,perennial_co2ice)
    5757
    5858contains
     
    7777    allocate(tsurf(ngrid,nslope))
    7878    allocate(watercap(ngrid,nslope))
    79     allocate(perenial_co2ice(ngrid,nslope))
     79    allocate(perennial_co2ice(ngrid,nslope))
    8080    allocate(emis(ngrid,nslope))
    8181    allocate(capcal(ngrid,nslope))
     
    108108    if (allocated(tsurf))            deallocate(tsurf)
    109109    if (allocated(watercap))         deallocate(watercap)
    110     if (allocated(perenial_co2ice))  deallocate(perenial_co2ice)
     110    if (allocated(perennial_co2ice))  deallocate(perennial_co2ice)
    111111    if (allocated(emis))             deallocate(emis)
    112112    if (allocated(capcal))           deallocate(capcal)
     
    130130    allocate(tsurf(ngrid,nslope))
    131131    allocate(watercap(ngrid,nslope))
    132     allocate(perenial_co2ice(ngrid,nslope))
     132    allocate(perennial_co2ice(ngrid,nslope))
    133133    allocate(emis(ngrid,nslope))
    134134    allocate(capcal(ngrid,nslope))
     
    144144    if (allocated(tsurf))            deallocate(tsurf)
    145145    if (allocated(watercap))         deallocate(watercap)
    146     if (allocated(perenial_co2ice))  deallocate(perenial_co2ice)
     146    if (allocated(perennial_co2ice)) deallocate(perennial_co2ice)
    147147    if (allocated(emis))             deallocate(emis)
    148148    if (allocated(capcal))           deallocate(capcal)
Note: See TracChangeset for help on using the changeset viewer.