Changeset 1677


Ignore:
Timestamp:
Mar 22, 2017, 4:18:23 PM (8 years ago)
Author:
sglmd
Message:

Two aerosol kinds added for giant planets: one is a compact (NH3) cloud where the opacity, particle size and bottom pressure level are taken as inputs (a scale height of 0.2 is hard-coded to simulate a compact cloud). Corresponds to option aeronh3. The other one is not generic at all and corresponds to the auroral, stratospheric aerosols observed on Jupiter...(option aeroaurora=.false. by default)

Location:
trunk/LMDZ.GENERIC/libf/phystd
Files:
8 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.GENERIC/libf/phystd/aeropacity.F90

    r1542 r1677  
    55       use aerosol_mod
    66       USE tracer_h, only: noms,rho_co2,rho_ice
    7        use comcstfi_mod, only: g
     7       use comcstfi_mod, only: g, pi
     8       use geometry_mod, only: latitude
    89       use callkeys_mod, only: aerofixco2,aerofixh2o,kastprof,cloudlvl, &
    910                CLFvarying,CLFfixval,dusttau,                           &
    1011                pres_bottom_tropo,pres_top_tropo,obs_tau_col_tropo,     &
    11                 pres_bottom_strato,pres_top_strato,obs_tau_col_strato
     12                pres_bottom_strato,pres_top_strato,obs_tau_col_strato,  &
     13                tau_nh3_cloud, pres_nh3_cloud
    1214                 
    1315       implicit none
     
    6062      logical,intent(in) :: clearsky
    6163
    62       real aerosol0
     64      real aerosol0, obs_tau_col_aurora, pm, pente_cloud
    6365
    6466      INTEGER l,ig,iq,iaer
     
    121123        if (iaero_back2lay.ne.0) then
    122124          print*,'iaero_back2lay= ',iaero_back2lay
     125        endif
     126        if (iaero_nh3.ne.0) then
     127          print*,'iaero_nh3= ',iaero_nh3
     128        endif
     129        if (iaero_aurora.ne.0) then
     130          print*,'iaero_aurora= ',iaero_aurora
    123131        endif
    124132
     
    409417      end if ! if Two-layer aerosols 
    410418
     419!==================================================================
     420!    Saturn/Jupiter ammonia cloud = thin cloud (scale height 0.2 hard coded...)
     421!    S. Guerlet (2013)
     422!==================================================================
     423
     424      if (iaero_nh3 .ne.0) then
     425           iaer=iaero_nh3
     426!       1. Initialization
     427            aerosol(1:ngrid,1:nlayer,iaer)=0.D0
     428!       2. Opacity calculation
     429          DO ig=1,ngrid
     430
     431           DO l=1,nlayer-1
     432            !! 1. below cloud layer: no opacity
     433           
     434            IF (pplev(ig,l) .gt. pres_nh3_cloud ) THEN
     435            aerosol(ig,l,iaer) = 0.D0           
     436
     437             ELSEIF (pplev(ig,l) .le. pres_nh3_cloud ) THEN
     438             pente_cloud=5. !!(hard-coded, correspond to scale height 0.2)
     439             aerosol(ig,l,iaer) = ((pplev(ig,l)/pres_nh3_cloud)**(pente_cloud))*tau_nh3_cloud
     440
     441             ENDIF
     442            ENDDO
     443
     444          END DO
     445         
     446!       3. Re-normalize to observed total column
     447         tau_col(:)=0.0
     448         DO l=1,nlayer
     449          DO ig=1,ngrid
     450               tau_col(ig) = tau_col(ig) &
     451                     + aerosol(ig,l,iaer)/tau_nh3_cloud
     452            ENDDO
     453         ENDDO
     454
     455         DO ig=1,ngrid
     456           DO l=1,nlayer-1
     457                aerosol(ig,l,iaer)=aerosol(ig,l,iaer)/tau_col(ig)
     458           ENDDO
     459         ENDDO
     460
     461     end if ! if NH3 cloud 
     462!==================================================================
     463!    Jovian auroral aerosols (unknown composition) NON-GENERIC: vertical and meridional profile tuned to observations
     464!    S. Guerlet (2015)
     465!==================================================================
     466
     467
     468      if (iaero_aurora .ne.0) then
     469           iaer=iaero_aurora
     470!       1. Initialization
     471            aerosol(1:ngrid,1:nlayer,iaer)=0.D0
     472         pm = 2000. !!case study: maxi aerosols at 20 hPa
     473!       2. Opacity calculation
     474          DO ig=1,ngrid
     475
     476          !! Test Jupiter (based on Zhang et al 2013 observations, but a bit different), decembre 2015
     477              DO l=1,nlayer
     478                aerosol(ig,l,iaer) = (pplev(ig,l)/pm)**2 * exp(-(pplev(ig,l)/pm)**2)
     479              ENDDO
     480          ENDDO
     481         
     482 !       3. Meridional distribution, and re-normalize to observed total column
     483         tau_col(:)=0.D0
     484         DO ig=1,ngrid
     485          !!Jupiter
     486          !!Hem sud:
     487          IF (latitude(ig)*180.D0/pi .lt. -45.D0 .and. latitude(ig)*180.D0/pi .gt. -70.) THEN
     488          obs_tau_col_aurora= 10.D0**(-0.06D0*latitude(ig)*180.D0/pi-3.4D0)
     489          ELSEIF (latitude(ig)*180.D0/pi .lt. -37.D0 .and. latitude(ig)*180.D0/pi .ge. -45.) THEN
     490          obs_tau_col_aurora= 10.D0**(-0.3D0*latitude(ig)*180.D0/pi-14.3D0)
     491           ELSEIF (latitude(ig)*180./pi .le. -70. ) THEN
     492          obs_tau_col_aurora= 10**(0.06*70.-3.4D0)
     493          !!Hem Nord: 
     494          ELSEIF (latitude(ig)*180.D0/pi .gt. 30.D0 .and. latitude(ig)*180.D0/pi .lt. 70.) THEN
     495          obs_tau_col_aurora= 10.D0**(0.03D0*latitude(ig)*180.D0/pi-1.17D0) 
     496          ELSEIF (latitude(ig)*180.D0/pi .gt. 22.D0 .and. latitude(ig)*180.D0/pi .le. 30.) THEN
     497          obs_tau_col_aurora= 10.D0**(0.3D0*latitude(ig)*180.D0/pi-9.4D0) 
     498          ELSEIF (latitude(ig)*180.D0/pi .ge. 70.) THEN
     499          obs_tau_col_aurora= 10**(0.03*70.-1.17D0) 
     500          ELSEIF (latitude(ig)*180.D0/pi .ge. -37. .and. latitude(ig)*180.D0/pi .le. 22.) THEN
     501         obs_tau_col_aurora = 0.001D0    !!Jupiter: mini pas a zero
     502          ENDIF
     503
     504          DO l=1,nlayer 
     505               tau_col(ig) = tau_col(ig) + aerosol(ig,l,iaer)/obs_tau_col_aurora
     506          ENDDO
     507         ENDDO
     508
     509         DO ig=1,ngrid
     510           DO l=1,nlayer-1
     511                aerosol(ig,l,iaer)=aerosol(ig,l,iaer)/tau_col(ig)
     512           ENDDO
     513         ENDDO
     514
     515
     516      end if ! if Auroral aerosols 
     517
    411518
    412519! --------------------------------------------------------------------------
  • trunk/LMDZ.GENERIC/libf/phystd/aerosol_mod.F90

    r1315 r1677  
    1616! two-layer simple aerosol model
    1717      integer :: iaero_back2lay = 0
    18 !$OMP THREADPRIVATE(iaero_co2,iaero_h2o,iaero_dust,iaero_h2so4,noaero,iaero_back2lay)
     18 ! NH3 cloud
     19      integer :: iaero_nh3 = 0
     20 ! Auroral aerosols
     21      integer :: iaero_aurora = 0
     22!$OMP THREADPRIVATE(iaero_co2,iaero_h2o,iaero_dust,iaero_h2so4,noaero,iaero_back2lay,iaero_nh3,iaero_aurora)
    1923     
    2024!==================================================================
  • trunk/LMDZ.GENERIC/libf/phystd/callcorrk.F90

    r1529 r1677  
    1616      use gases_h
    1717      use radii_mod, only : su_aer_radii,co2_reffrad,h2o_reffrad,dust_reffrad,h2so4_reffrad,back2lay_reffrad
    18       use aerosol_mod, only : iaero_co2,iaero_h2o,iaero_dust,iaero_h2so4, iaero_back2lay
     18      use aerosol_mod, only : iaero_co2,iaero_h2o,iaero_dust,iaero_h2so4, iaero_back2lay, iaero_nh3, iaero_aurora
    1919      USE tracer_h
    2020      use comcstfi_mod, only: pi, mugaz, cpp
     
    305305            call back2lay_reffrad(ngrid,reffrad(1,1,iaero_back2lay),nlayer,pplev)
    306306         endif
    307          
     307!         if(iaer.eq.iaero_nh3)then
     308!           call nh3_reffrad(ngrid,nlayer,reffrad(1,1,iaero_nh3))
     309!         endif
     310!         if(iaer.eq.iaero_aurora)then
     311!           call aurora_reffrad(ngrid,nlayer,reffrad(1,1,iaero_aurora))
     312!         endif
     313       
    308314     end do !iaer=1,naerkind.
    309315
  • trunk/LMDZ.GENERIC/libf/phystd/callkeys_mod.F90

    r1669 r1677  
    4242      logical,save :: aeroco2,aeroh2o,aeroh2so4,aeroback2lay
    4343!$OMP THREADPRIVATE(aeroco2,aeroh2o,aeroh2so4,aeroback2lay)
     44      logical,save :: aeronh3, aeroaurora
     45!$OMP THREADPRIVATE(aeronh3,aeroaurora)
    4446      logical,save :: aerofixco2,aerofixh2o
    4547!$OMP THREADPRIVATE(aerofixco2,aerofixh2o)
     
    8486      real,save :: n2mixratio
    8587!$OMP THREADPRIVATE(size_tropo,size_strato,satval,CLFfixval,n2mixratio)
     88      real,save :: size_nh3_cloud
     89      real,save :: pres_nh3_cloud
     90      real,save :: tau_nh3_cloud
     91!$OMP THREADPRIVATE(size_nh3_cloud, pres_nh3_cloud, tau_nh3_cloud)
    8692      real,save :: co2supsat
    8793      real,save :: pceil
  • trunk/LMDZ.GENERIC/libf/phystd/iniaerosol.F

    r1397 r1677  
    55      use aerosol_mod
    66      use callkeys_mod, only: aeroco2,aeroh2o,dusttau,aeroh2so4,
    7      &          aeroback2lay
     7     &          aeroback2lay,aeronh3, aeroaurora
    88
    99      IMPLICIT NONE
     
    5252      write(*,*) '--- Two-layer aerosol = ', iaero_back2lay
    5353
     54      if (aeronh3) then
     55         ia=ia+1
     56         iaero_nh3=ia
     57      endif
     58      write(*,*) '--- NH3 Cloud = ', iaero_nh3
     59
     60      if (aeroaurora) then
     61         ia=ia+1
     62         iaero_aurora=ia
     63      endif
     64      write(*,*) '--- Auroral aerosols = ', iaero_aurora
     65
    5466      write(*,*) '=== Number of aerosols= ', ia
    5567     
  • trunk/LMDZ.GENERIC/libf/phystd/inifis_mod.F90

    r1669 r1677  
    498498     write(*,*)" aerofixh2o = ",aerofixh2o
    499499
    500      write(*,*)"Radiatively active sulfuric acid aersols?"
     500     write(*,*)"Radiatively active sulfuric acid aerosols?"
    501501     aeroh2so4=.false.     ! default value
    502502     call getin_p("aeroh2so4",aeroh2so4)
     
    505505!=================================
    506506
    507      write(*,*)"Radiatively active two-layer aersols?"
     507     write(*,*)"Radiatively active two-layer aerosols?"
    508508     aeroback2lay=.false.     ! default value
    509509     call getin_p("aeroback2lay",aeroback2lay)
    510510     write(*,*)" aeroback2lay = ",aeroback2lay
     511
     512     write(*,*)"Radiatively active ammonia cloud?"
     513     aeronh3=.false.     ! default value
     514     call getin_p("aeronh3",aeronh3)
     515     write(*,*)" aeronh3 = ",aeronh3
     516
     517     write(*,*)"Radiatively active auroral aerosols?"
     518     aeroaurora=.false.     ! default value
     519     call getin_p("aeroaurora",aeroaurora)
     520     write(*,*)" aeroaurora = ",aeroaurora
    511521
    512522     write(*,*)"TWOLAY AEROSOL: total optical depth ", &
     
    553563     call getin_p("size_strato",size_strato)
    554564     write(*,*)" size_strato = ",size_strato
     565
     566     write(*,*)"NH3 (thin) cloud: total optical depth"
     567     tau_nh3_cloud=7.D0
     568     call getin_p("tau_nh3_cloud",tau_nh3_cloud)
     569     write(*,*)" tau_nh3_cloud = ",tau_nh3_cloud
     570
     571     write(*,*)"NH3 (thin) cloud pressure level"
     572     pres_nh3_cloud=7.D0
     573     call getin_p("pres_nh3_cloud",pres_nh3_cloud)
     574     write(*,*)" pres_nh3_cloud = ",pres_nh3_cloud
     575
     576     write(*,*)"NH3 (thin) cloud: particle sizes"
     577     size_nh3_cloud=3.e-6
     578     call getin_p("size_nh3_cloud",size_nh3_cloud)
     579     write(*,*)" size_nh3_cloud = ",size_nh3_cloud
    555580
    556581!=================================
  • trunk/LMDZ.GENERIC/libf/phystd/radii_mod.F90

    r1529 r1677  
    3939      use radinc_h, only: naerkind
    4040      use aerosol_mod, only: iaero_back2lay, iaero_co2, iaero_dust, &
    41                              iaero_h2o, iaero_h2so4
     41                             iaero_h2o, iaero_h2so4,iaero_nh3,iaero_aurora
    4242      Implicit none
    4343
     
    8484            endif
    8585
     86
     87            if(iaer.eq.iaero_nh3)then ! Nh3 cloud
     88               reffrad(1:ngrid,1:nlayer,iaer) = 1e-6
     89               nueffrad(1:ngrid,1:nlayer,iaer) = 0.1
     90            endif
     91
     92            if(iaer.eq.iaero_aurora)then ! Auroral aerosols
     93               reffrad(1:ngrid,1:nlayer,iaer) = 3e-7
     94               nueffrad(1:ngrid,1:nlayer,iaer) = 0.1
     95            endif
    8696
    8797
     
    363373
    364374
    365 
    366375end module radii_mod
    367376!==================================================================
  • trunk/LMDZ.GENERIC/libf/phystd/suaer_corrk.F90

    r1529 r1677  
    180180! added by SG
    181181       endif
    182        
     182     
     183       if (iaer.eq.iaero_nh3) then
     184         print*, 'naerkind= nh3', iaer
     185
     186!     visible
     187         file_id(iaer,1) = 'optprop_nh3ice_vis.dat'
     188         lamrefvis(iaer)=0.756E-6  !
     189!     infrared
     190         file_id(iaer,2) = 'optprop_nh3ice_ir.dat'
     191         lamrefir(iaer)=6.E-6  !
     192! added by SG
     193       endif
     194
     195
     196       if (iaer.eq.iaero_aurora) then
     197         print*, 'naerkind= aurora', iaer
     198
     199!     visible
     200         file_id(iaer,1) = 'optprop_aurora_vis.dat'
     201         lamrefvis(iaer)=0.3E-6  !
     202!     infrared
     203         file_id(iaer,2) = 'optprop_aurora_ir.dat'
     204         lamrefir(iaer)=6.E-6  !
     205! added by SG
     206       endif
     207 
    183208       
    184209      enddo
Note: See TracChangeset for help on using the changeset viewer.