Ignore:
Timestamp:
May 23, 2024, 3:05:55 PM (7 months ago)
Author:
llange
Message:

Mars PCM
Small correction of how qsurf of co2 is initialized in the south pole.
Some correction for the albedo of co2 ice when using paleoclimate: a uniform value (albedo_co2_cap) is used for both hemisphere for CO2 frost, and a different (higher) value is used for perennial ice.
LL

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

Legend:

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

    r3159 r3343  
    1111USE mod_phys_lmdz_transfert_para, ONLY: bcast
    1212USE mod_phys_lmdz_para, ONLY: is_master
    13 USE paleoclimate_mod, ONLY: paleoclimate,albedo_perennialco2
     13USE paleoclimate_mod, ONLY: paleoclimate,albedo_perennialco2,albedo_co2_cap
    1414
    1515implicit none
     
    105105      psolaralb(ig,2)=psolaralb(ig,1)
    106106    else
    107       psolaralb(ig,:)=albedice(icap)
    108       if (paleoclimate .and. piceco2_peren(ig) > 0. .and. abs(piceco2(ig)) < 1.e-10) psolaralb(ig,:) = albedo_perennialco2
    109     endif
     107      if(paleoclimate) then
     108            psolaralb(ig,1) = albedo_co2_cap
     109            psolaralb(ig,2) = albedo_co2_cap
     110      else
     111            psolaralb(ig,1) = albedice(icap)
     112            psolaralb(ig,2) = albedice(icap)
     113      endif
     114    endif
     115  else if(paleoclimate .and. piceco2_peren(ig) > 0.) then
     116     psolaralb(ig,1) = albedo_perennialco2
     117     psolaralb(ig,2) = albedo_perennialco2
    110118  else if (watercaptag(ig) .and. water) then
    111119    ! there is a water ice cap: set the surface albedo to the water ice one
  • trunk/LMDZ.MARS/libf/phymars/conf_phys.F

    r3333 r3343  
    3838      USE mod_phys_lmdz_transfert_para, ONLY: bcast
    3939      USE paleoclimate_mod,ONLY: paleoclimate,albedo_perennialco2,
     40     &                           albedo_co2_cap,
    4041     &                           lag_layer, include_waterbuoyancy
    4142      use microphys_h, only: mteta
     
    375376         call getin_p("albedo_perennialco2",albedo_perennialco2)
    376377         write(*,*)"albedo_perennialco2 = ",albedo_perennialco2
     378
     379         write(*,*)"Albedo for seasonal CO2 ice?"
     380         albedo_co2_cap = 0.6 ! default value
     381         call getin_p("albedo_co2_cap",albedo_co2_cap)
     382         write(*,*)"albedo_co2_cap = ",albedo_co2_cap
    377383
    378384         write(*,*)"Include water buoyancy effect??"
  • trunk/LMDZ.MARS/libf/phymars/paleoclimate_mod.F90

    r3336 r3343  
    2121    logical, save                              :: lag_layer             ! Does lag layer is present?
    2222    logical, save                              :: include_waterbuoyancy ! Include the effect of water buoyancy when computing the sublimation of water ice ?
    23 !$OMP THREADPRIVATE(h2o_ice_depth,lag_co2_ice,d_coef,albedo_perennialco2,lag_layer,include_waterbuoyancy)
     23    real,    save                              :: albedo_co2_cap
     24!$OMP THREADPRIVATE(h2o_ice_depth,lag_co2_ice,d_coef,albedo_perennialco2,lag_layer,include_waterbuoyancy,albedo_co2_cap)
    2425
    2526!=======================================================================
  • trunk/LMDZ.MARS/libf/phymars/phyetat0_mod.F90

    r3316 r3343  
    842842        if (abs(latitude(ngrid) - (-pi/2.)) < 1.e-5) then
    843843            perennial_co2ice(ngrid,:) = 10*1.6e3 ! 10m which is convert to kg/m^2
    844             qsurf(ngrid,igcm_co2_tmp,:) = qsurf(ngrid - 1,igcm_co2_tmp,:) + perennial_co2ice(ngrid,:) ! perennial ice + frost
     844            qsurf(ngrid,igcm_co2_tmp,:) = qsurf(ngrid - 1,igcm_co2_tmp,:)
    845845        endif
    846846    endif ! not found
     
    852852    if (abs(latitude(ngrid) - (-pi/2.)) < 1.e-5) then
    853853        perennial_co2ice(ngrid,:) = 10*1.6e3 ! 10m which is convert to kg/m^2
    854         qsurf(ngrid,igcm_co2_tmp,:) = qsurf(ngrid - 1,igcm_co2_tmp,:) + perennial_co2ice(ngrid,:) ! perennial ice + frost
     854        qsurf(ngrid,igcm_co2_tmp,:) = qsurf(ngrid - 1,igcm_co2_tmp,:)
    855855    endif
    856856  endif !startphy_file
Note: See TracChangeset for help on using the changeset viewer.