Changeset 2634


Ignore:
Timestamp:
Mar 3, 2022, 10:17:13 PM (3 years ago)
Author:
abierjon
Message:

Mars GCM:
Changes on dust_rad_adjust for GCM6:

  • put an upper limit on dust_rad_adjust to avoid skyrocketing values when there are strong diurnal variations of opacity
  • move the condition to compute dust_rad_adjust only once per time step before the whole call to compute_dust_rad_adjust, instead of repetitive local checks
  • some cleaning

AB

Location:
trunk/LMDZ.MARS
Files:
7 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.MARS/README

    r2630 r2634  
    36303630startfi.nc in the output file start_archive.nc, and newstart reads the variable
    36313631albedo in start_archive.nc before outputing it in restartfi.nc
     3632
     3633== 03/03/2022 == AB
     3634Changes on dust_rad_adjust for GCM6:
     3635- put an upper limit on dust_rad_adjust to avoid skyrocketing values when there are strong diurnal
     3636  variations of opacity
     3637- move the condition to compute dust_rad_adjust only once per time step before the whole call to
     3638  compute_dust_rad_adjust, instead of repetitive local checks
     3639- some cleaning
  • trunk/LMDZ.MARS/libf/phymars/aeropacity_mod.F

    r2628 r2634  
    109109                             ! false to compute RT in the topdust
    110110      REAL,INTENT(OUT) :: tauscaling(ngrid) ! Scaling factor for qdust and Ndust
    111       REAL,INTENT(OUT) :: dust_rad_adjust(ngrid) ! Radiative adjustment
     111      REAL,INTENT(INOUT) :: dust_rad_adjust(ngrid) ! Radiative adjustment
    112112                          ! factor for dust
    113113      REAL,INTENT(IN) :: totcloudfrac(ngrid) ! total water ice cloud fraction
     
    597597c       and once in the part of the mesh without the sub-grid mountain (nohmons=true)
    598598        aerosol(1:ngrid,1:nlayer,iaer) = 0.
    599         IF (nohmons) THEN  ! considering part of the mesh without top flows
     599        IF (nohmons) THEN  ! considering part of the mesh without top dust flows
    600600          aerosol(1:ngrid,1:nlayer,iaer)=1.e-25
    601601        ELSE  ! part of the mesh with concentrated dust flows
  • trunk/LMDZ.MARS/libf/phymars/callradite_mod.F

    r2628 r2634  
    177177      REAL,INTENT(INOUT) :: tauscaling(ngrid) ! Conversion factor for
    178178                               ! qdust and Ndust
    179       REAL,INTENT(OUT) :: dust_rad_adjust(ngrid) ! Radiative adjustment
     179      REAL,INTENT(INOUT) :: dust_rad_adjust(ngrid) ! Radiative adjustment
    180180                          ! factor for dust
    181181      REAL,INTENT(IN) :: albedo(ngrid,2),emis(ngrid)
  • trunk/LMDZ.MARS/libf/phymars/dust_rad_adjust_mod.F90

    r2616 r2634  
    2323  integer,intent(in) :: ngrid ! number of atmospheric columns
    2424  integer,intent(in) :: nlayer ! number of atmospheric levels
    25   real,intent(in) :: zday ! tim (in sols and fraction thereof)
     25  real,intent(in) :: zday ! time (in sols and fraction thereof)
    2626  real,intent(in) :: pplev(ngrid,nlayer+1) ! pressure (Pa) at layer boundaries
    2727  real,intent(in) :: taudust(ngrid) ! visible dust columns opacity in the GCM
     
    3131  real,allocatable,save :: local_time(:) ! LT at current physics time step
    3232  real,allocatable,save :: local_time_prevdt(:) ! LT at previous physics time step
    33   real :: zday_prevdt !value of zday at previous physics time step
     33  real :: zday_prevdt ! value of zday at previous physics time step
    3434  real,save :: zday_scenario ! to fetch dod values from the scenario
    3535  real,save :: zday_scenario_next ! to fetch dod values from the scenario the next day
    36   logical,save :: firstcall=.true.
     36  logical,save :: firstcall=.true. ! is it the first call of the run?
    3737  integer :: ig
    3838!  real,allocatable,save :: tau_pref_scenario(:)
    3939  real,allocatable,save :: tau_pref_scenario_next(:)
    4040  real :: weight ! interpolation weight
    41   real,save :: zday_prev_call=-666. ! stored value of zday from previous call
    4241 
    43 !$OMP THREADPRIVATE(  local_time,local_time_prevdt,zday_scenario,zday_scenario_next)
    44 !$OMP THREADPRIVATE(firstcall,tau_pref_scenario_next,zday_prev_call)
     42!$OMP THREADPRIVATE(local_time,local_time_prevdt,zday_scenario,zday_scenario_next)
     43!$OMP THREADPRIVATE(firstcall,tau_pref_scenario_next)
    4544
    4645 
     
    6059  endif ! of if firstcall
    6160 
    62   ! 1. Compute local times (in sol fraction), if not already done
    63   if (zday/=zday_prev_call) then
    64     local_time(1:ngrid)=modulo(1.+(zday-INT(zday)) + &
    65                          (longitude_deg(1:ngrid)/15)/24,1.)
    66     zday_prevdt=zday-dtphys/daysec
    67     local_time_prevdt(1:ngrid)=modulo(1.+(zday_prevdt-INT(zday_prevdt)) + &
    68                          (longitude_deg(1:ngrid)/15)/24,1.)
    69    
    70     zday_scenario=zday-modulo(zday,1.) ! integer value of the day: the scenario
    71     ! opacity is assumed to be measured at 2pm but stored at nidnight
    72     zday_scenario_next=zday_scenario+1
    73   endif ! of if (zday/=zday_prev_call)
     61  ! 1. Compute local times (in sol fraction)
     62  local_time(1:ngrid)=modulo(1.+(zday-INT(zday)) + &
     63                       (longitude_deg(1:ngrid)/15)/24,1.)
     64  zday_prevdt=zday-dtphys/daysec
     65  local_time_prevdt(1:ngrid)=modulo(1.+(zday_prevdt-INT(zday_prevdt)) + &
     66                       (longitude_deg(1:ngrid)/15)/24,1.)
     67
     68  zday_scenario=zday-modulo(zday,1.) ! integer value of the day: the scenario
     69  ! opacity is assumed to be measured at 2pm but stored at nidnight
     70  zday_scenario_next=zday_scenario+1
    7471 
    7572  ! 2. Load dust opacities for zday_scenario and zday_scenario_next
    76   !    if not already done
    77   if (zday/=zday_prev_call) then
    78 !    call read_dust_scenario(ngrid,nlayer,zday_scenario,pplev,     &
    79 !                                         tau_pref_scenario)
    80     call read_dust_scenario(ngrid,nlayer,zday_scenario_next,pplev, &
    81                                          tau_pref_scenario_next)
    82   endif ! of if (zday/=zday_prev_call)
     73!  call read_dust_scenario(ngrid,nlayer,zday_scenario,pplev,     &
     74!                                       tau_pref_scenario)
     75  call read_dust_scenario(ngrid,nlayer,zday_scenario_next,pplev, &
     76                                       tau_pref_scenario_next)
    8377
    8478  ! 3. Update dust_rad_adjust_* for grid points which just reached 2pm
    85   ! but only when this routine is called for the first time
    86   ! during this time step
    87   if (zday/=zday_prev_call) then
    88    do ig=1,ngrid
    89     if ((local_time(ig).ge.t_scenario_sol).and. &
    90                  (local_time_prevdt(ig).lt.(t_scenario_sol))) then
    91       ! store previous "next" as "prev" (NB we could also decide to recompute
    92       ! it using the current taudust...)
    93       dust_rad_adjust_prev(ig)=dust_rad_adjust_next(ig)
    94       ! compute new target based on current dust opacity
    95       dust_rad_adjust_next(ig)=tau_pref_scenario_next(ig)* &
    96                                pplev(ig,1)/odpref/taudust(ig)
    97     endif
    98    enddo
    99   endif ! of if (zday/=zday_prev_call)
     79  do ig=1,ngrid
     80   if ((local_time(ig).ge.t_scenario_sol).and. &
     81                (local_time_prevdt(ig).lt.(t_scenario_sol))) then
     82     ! store previous "next" as "prev" (NB we could also decide to recompute
     83     ! it using the current taudust...)
     84     dust_rad_adjust_prev(ig)=dust_rad_adjust_next(ig)
     85     ! compute new target based on current dust opacity
     86     dust_rad_adjust_next(ig)=tau_pref_scenario_next(ig)* &
     87                              pplev(ig,1)/odpref/taudust(ig)
     88                             
     89     ! Upper limit for dust_rad_adjust
     90     ! this is done to avoid skyrocketing taus given to the radiative transfer
     91     ! when there is a huge diurnal variation of tau (especially in polar nights)
     92      dust_rad_adjust_next(ig)=min(dust_rad_adjust_next(ig), 5.) !3.
     93     
     94   endif
     95  enddo
    10096 
    10197  ! 4. Compute dust_rad_adjust using linear interpolation
     
    116112                        (dust_rad_adjust_next(ig)-dust_rad_adjust_prev(ig))
    117113  enddo! of do=ig=1,ngrid
    118  
    119   ! update zday_prev_call
    120   zday_prev_call=zday
    121114 
    122115  end subroutine compute_dust_rad_adjust
  • trunk/LMDZ.MARS/libf/phymars/dust_scaling_mod.F90

    r2417 r2634  
    2525                       ! opacity column at odpref reference pressure
    2626    real,intent(out) :: tauscaling(ngrid) ! dust scaling factor
    27     real,intent(out) :: dust_rad_adjust(ngrid) ! Radiative adjustment
     27    real,intent(inout) :: dust_rad_adjust(ngrid) ! Radiative adjustment
    2828                          ! factor for dust
    2929    real,intent(inout) :: aerosol(ngrid,nlayer,naerkind) ! opacities
     
    3131    integer :: ig, l , iaer
    3232    real :: taudust(ngrid)
     33    real,save :: zday_prev_call=-666. ! stored value of zday from previous call
     34!$OMP THREADPRIVATE(zday_prev_call)
    3335
    3436  ! 1. compute/set tauscaling
     
    6870    elseif (dustscaling_mode==2) then
    6971      ! GCM v6 style, compute dust_rad_adjust
    70       call compute_dust_rad_adjust(ngrid,nlayer,zday,pplev, &
    71                                    taudust,dust_rad_adjust)
     72      ! only when this routine is called for the first time
     73      ! during this time step
     74      if (zday/=zday_prev_call) then
     75        call compute_dust_rad_adjust(ngrid,nlayer,zday,pplev, &
     76                                     taudust,dust_rad_adjust)
     77       endif ! of if (zday/=zday_prev_call)
     78                                     
     79      ! update zday_prev_call
     80      zday_prev_call=zday
    7281    endif
    7382
  • trunk/LMDZ.MARS/libf/phymars/rocketduststorm_mod.F90

    r2628 r2634  
    7979      REAL, INTENT(IN) :: totstormfract(ngrid)
    8080      REAL, INTENT(INOUT) :: tauscaling(ngrid)
    81       REAL,INTENT(OUT) :: dust_rad_adjust(ngrid)
     81      REAL,INTENT(INOUT) :: dust_rad_adjust(ngrid)
    8282     
    8383!     subgrid scale water ice clouds
  • trunk/LMDZ.MARS/libf/phymars/topmons_mod.F90

    r2628 r2634  
    7474      INTEGER, INTENT(IN) :: igout
    7575      REAL, INTENT(INOUT) :: tauscaling(ngrid)
    76       REAL,INTENT(OUT) :: dust_rad_adjust(ngrid)
     76      REAL,INTENT(INOUT) :: dust_rad_adjust(ngrid)
    7777!     input sub-grid scale rocket dust storm
    7878      LOGICAL, INTENT(IN) :: clearatm
Note: See TracChangeset for help on using the changeset viewer.