Changeset 3159 for trunk/LMDZ.MARS


Ignore:
Timestamp:
Dec 14, 2023, 6:56:40 PM (11 months ago)
Author:
jbclement
Message:

PEM:

  • A number of corrections related to r3149 for CO2 ice management:

    'co2_ice' was not transferred to 'perennial_co2ice' at the end of the PEM run;
    In the Mars PCM, 'perennial_co2ice' is now correctly handled regarding the frost in case of CO2 sublimation/condensation.

  • Addition of (CO2 and H2O) ice metamorphism: if the frost is above a given value, then the excess frost is transferred to the perennial ice.
  • Thresholds value for ice management can now be set in the "run_PEM.def".

JBC

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

Legend:

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

    r3130 r3159  
    105105      psolaralb(ig,2)=psolaralb(ig,1)
    106106    else
    107       psolaralb(ig,1)=albedice(icap)
    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_perennialco2
    112              psolaralb(ig,2) = albedo_perennialco2
    113              piceco2_peren(ig) = piceco2(ig)
    114          endif
    115       endif
     107      psolaralb(ig,:)=albedice(icap)
     108      if (paleoclimate .and. piceco2_peren(ig) > 0. .and. abs(piceco2(ig)) < 1.e-10) psolaralb(ig,:) = albedo_perennialco2
    116109    endif
    117110  else if (watercaptag(ig) .and. water) then
  • trunk/LMDZ.MARS/libf/phymars/co2condens_mod.F

    r3142 r3159  
    55      logical, save :: scavco2cond = .false. ! flag for using scavenging_by_co2
    66!$OMP THREADPRIVATE(scavco2cond)
    7      
     7
    88      CONTAINS
    99
     
    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   :: perennial_co2ice_tmp(ngrid)    ! temporary amount of perennial 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]
    191191      REAL :: condens_layer_tmp(ngrid,nlayer) ! temporary condensation rate in layer l (co2 cloud) for a specific sub-grid surface [kg/m^2/s]
    192192      INTEGER :: islope                       ! index for loop variables
    193       REAL :: pdqsc_tmp(ngrid,nq)         ! tendency on surface tracers (grid-mesh average) after scavenging
     193      REAL :: pdqsc_tmp(ngrid,nq)             ! tendency on surface tracers (grid-mesh average) after scavenging
    194194c----------------------------------------------------------------------
    195195
     
    227227      ENDIF ! of IF (firstcall)
    228228      zcpi=1./cpp
     229      if (paleoclimate) piceco2 = piceco2 + perennial_co2ice
    229230
    230231c
     
    564565c  ********************************************************************
    565566
     567!     Redistribute piceco2 into piceco2 (frost) and perennial_co2ice
     568!     --------------------------------------------------------------
     569      if (paleoclimate) then
     570          where (piceco2 > perennial_co2ice) ! Perennial co2 ice has not been affected
     571              ! It means:
     572              !     - In case of sublimation, only frost is lost
     573              !     - In case of condensation, only frost accumulates new ice
     574              piceco2 = piceco2 - perennial_co2ice
     575          else where ! Perennial co2 ice has been affected
     576              ! It means that frost disappeared with sublimation and perennial ice is being lost
     577              perennial_co2ice = piceco2
     578              piceco2 = 0.
     579          end where
     580      endif
     581
    566582!     Check that amont of CO2 ice is not problematic
    567       DO ig=1,ngrid
     583!     ----------------------------------------------
     584      DO ig = 1,ngrid
    568585         DO islope = 1,nslope
    569            if(.not.piceco2(ig,islope).ge.0.) THEN
    570               if(piceco2(ig,islope).le.-5.e-8) print*,
    571      $         'WARNING co2condens piceco2(',ig,')=', piceco2(ig,islope)
    572                piceco2(ig,islope)=0.
     586           if (piceco2(ig,islope) < 0.) then
     587              if (piceco2(ig,islope) <= -5.e-8) print*,
     588     $         'WARNING co2condens piceco2(',ig,',',islope,
     589     $         ') =',piceco2(ig,islope)
     590               piceco2(ig,islope) = 0.
    573591           endif
    574592        ENDDO
    575593      ENDDO
     594      if (paleoclimate) then
     595      DO ig = 1,ngrid
     596         DO islope = 1,nslope
     597           if (perennial_co2ice(ig,islope) < 0.) then
     598              if (perennial_co2ice(ig,islope) <= -5.e-8) print*,
     599     $         'WARNING co2condens perennial_co2ice(',ig,',',islope,
     600     $         ') =',perennial_co2ice(ig,islope)
     601               perennial_co2ice(ig,islope) = 0.
     602           endif
     603        ENDDO
     604      ENDDO
     605      endif
    576606     
    577607!     Set albedo and emissivity of the surface
    578608!     ----------------------------------------
    579609      DO islope = 1,nslope
    580         piceco2_tmp(:) = piceco2(:,islope)
    581         alb_tmp(:,:) = psolaralb(:,:,islope)
    582         emisref_tmp(:) = 0.
    583         perennial_co2ice_tmp(:) =  perennial_co2ice(:,islope)
     610        piceco2_tmp = piceco2(:,islope)
     611        alb_tmp = psolaralb(:,:,islope)
     612        emisref_tmp = 0.
     613        perennial_co2ice_tmp =  perennial_co2ice(:,islope)
    584614        CALL albedocaps(zls,ngrid,piceco2_tmp,perennial_co2ice_tmp,
    585615     &                  alb_tmp,emisref_tmp)
    586         perennial_co2ice(:,islope) = perennial_co2ice_tmp(:)
     616        perennial_co2ice(:,islope) = perennial_co2ice_tmp
    587617        psolaralb(:,1,islope) =  alb_tmp(:,1)
    588618        psolaralb(:,2,islope) =  alb_tmp(:,2)
    589         emisref(:,islope) = emisref_tmp(:)
     619        emisref(:,islope) = emisref_tmp
    590620      ENDDO
    591621
    592 ! set pemisurf() to emissiv when there is bare surface (needed for co2snow)
    593       DO ig=1,ngrid
    594         DO islope = 1,nslope
    595           if (piceco2(ig,islope).eq.0) then
    596             pemisurf(ig,islope)=emissiv
    597           endif
    598         ENDDO
    599       ENDDO
     622!     Set pemisurf() to emissiv when there is bare surface (needed for co2snow)
     623!     -------------------------------------------------------------------------
     624      where (piceco2 == 0.) pemisurf = emissiv
     625
    600626
    601627!         firstcall2=.false.
  • trunk/LMDZ.MARS/libf/phymars/physiq_mod.F

    r3157 r3159  
    24732473      IF (refill_watercap) THEN
    24742474
    2475         DO ig=1,ngrid
     2475        DO ig = 1,ngrid
    24762476         DO islope = 1,nslope
    2477            if (watercaptag(ig).and. (qsurf(ig,igcm_h2o_ice,islope)
    2478      &        .gt.frost_metam_threshold)) then
    2479 
    2480                 watercap(ig,islope)=watercap(ig,islope)
    2481      &          +qsurf(ig,igcm_h2o_ice,islope)
    2482      &         - frost_metam_threshold
     2477           if (watercaptag(ig) .and. (qsurf(ig,igcm_h2o_ice,islope)
     2478     &        > frost_metam_threshold)) then
     2479
     2480                watercap(ig,islope) = watercap(ig,islope)
     2481     &          + qsurf(ig,igcm_h2o_ice,islope)
     2482     &          - frost_metam_threshold
    24832483                qsurf(ig,igcm_h2o_ice,islope) = frost_metam_threshold
    2484            endif ! (watercaptag(ig).and.
     2484           endif ! watercaptag
    24852485          ENDDO
    24862486        ENDDO
    24872487
    2488       ENDIF ! (refill_watercap) THEN
     2488      ENDIF ! refill_watercap
    24892489
    24902490
Note: See TracChangeset for help on using the changeset viewer.