Changeset 2803


Ignore:
Timestamp:
Oct 19, 2022, 11:24:43 AM (2 years ago)
Author:
aslmd
Message:

Initialisation of Radiative Generic Condensable Aerosols

We can activate the scheme by putting aerogeneric = # of aerosols in callphys.def. This is the only needed thing for activating the radiative effects. They must be tracer in modern tracer.def

Commented out the abort if we use more than 4 aerosols

Added reading of optprop files for Radiative Generic Condensable Aerosols

We use the same file for IR and VI channel. For now, only MnS, Cr,Fe and Mg2SiO4 can be read. If you want to add another specie, check the code, it is explained how to quickly do that (right above the Initializations)

Added radii calculation for Radiative Generic Condensable aerosols

Changed the hardcoded size of the totalemit array

The hardcoded size is now 1900 instead of 100 so we don't exceed the array size when working at high spectral resolution (very rare case)

Added opacity computation for Radiative Generic Condensable aerosols

We do this computation in the same fashion as what's performed on water and dust.

switch iniaerosol and initracer order, to prepare for the RGCS scheme

Needed to switch the order of initialization so we can use the RGCS scheme without the assumption that ice and vap tracer of the same specie are following each other in the traceur.def file

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

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.GENERIC/libf/phystd/aerave_new.F

    r2625 r2803  
    5656c
    5757      INTEGER iir
    58       INTEGER,PARAMETER :: nirmx=100
     58      INTEGER,PARAMETER :: nirmx=1900
    5959      INTEGER idata,ndata
    6060c
  • trunk/LMDZ.GENERIC/libf/phystd/aeropacity.F90

    r2297 r2803  
    44       use radinc_h, only : L_TAUMAX,naerkind
    55       use aerosol_mod
    6        USE tracer_h, only: noms,rho_co2,rho_ice
    7        use comcstfi_mod, only: g, pi
     6       USE tracer_h, only: noms,rho_co2,rho_ice,rho_q,mmol
     7       use comcstfi_mod, only: g, pi, mugaz, avocado
    88       use geometry_mod, only: latitude
    99       use callkeys_mod, only: aerofixco2,aerofixh2o,kastprof,cloudlvl, &
     
    1313                tau_nh3_cloud, pres_nh3_cloud,                          &
    1414                nlayaero, aeronlay_tauref, aeronlay_choice,             &
    15                 aeronlay_pbot, aeronlay_ptop, aeronlay_sclhght
    16                  
     15                aeronlay_pbot, aeronlay_ptop, aeronlay_sclhght,         &
     16                aerogeneric
     17        use generic_tracer_index_mod, only: generic_tracer_index
    1718       implicit none
    1819
     
    3031!     dust removal, simplification by Robin Wordsworth (2009)
    3132!     Generic n-layer aerosol - J. Vatant d'Ollone (2020)
     33!     Radiative Generic Condensable aerosols - Lucas Teinturier (2022)
    3234!
    3335!     Input
     
    9092
    9193      real CLFtot
    92 
     94      integer igen_ice,igen_vap ! to store the index of generic tracer
     95      logical dummy_bool ! dummy boolean just in case we need one
    9396      ! identify tracers
    9497      IF (firstcall) THEN
     
    141144          print*,'iaero_aurora= ',iaero_aurora
    142145        endif
    143 
     146        if (iaero_generic(1) .ne. 0) then
     147          print*,"iaero_generic= ",iaero_generic(:)
     148        endif
    144149        firstcall=.false.
    145150      ENDIF ! of IF (firstcall)
     
    362367            ENDDO
    363368         ENDDO
    364 
     369         
    365370! 1/700. is assuming a "sulfurtau" of 1
    366371! Sulfur aerosol routine to be improved.
     
    645650
    646651      end if ! if Auroral aerosols 
    647 
     652!===========================================================================
     653!    Radiative Generic Condensable aerosols scheme
     654!    Only used when we give aerogeneric != 0 in callphys.def
     655!    Computes the generic aerosols' opacity in the same fashion as water of
     656!    dust, using the QREFvis3d of the concerned specie
     657!    Lucas Teinturier (2022)
     658!===========================================================================
     659      if (iaero_generic(1) .ne. 0) then ! we enter the scheme
     660        do ia=1,aerogeneric
     661          iaer = iaero_generic(ia)
     662          ! Initialization
     663          aerosol(1:ngrid,1:nlayer,iaer) = 0.D0
     664          igen_ice = -1
     665          igen_vap = -1
     666          ! Get index of the ice tracer (for pq, mol, rho_q)
     667          call generic_tracer_index(nq,2*iaer,igen_vap,igen_ice,dummy_bool)
     668          if (igen_ice .eq. -1) then
     669            !if igen_ice not changed by generic_tracer_index, it means it's an ice tracer => igen_ice=2*iaer
     670            igen_ice = 2*iaer
     671          endif
     672          ! Let's loop on the horizontal and vertical grid
     673          do ig=1,ngrid
     674            do l=1,nlayer
     675              aerosol(ig,l,iaer) = ( 0.75*QREFvis3d(ig,l,iaer)  / &
     676                                  (rho_q(igen_ice) * reffrad(ig,l,iaer)) ) * &
     677                                  (pq(ig,l,igen_ice)+1E-9 ) *                &
     678                                  (pplev(ig,l) - pplev(ig,l+1)) /g
     679            enddo !l=1,nlayer
     680          enddo !ig=1,ngrid
     681        enddo !ia=1,aerogeneric
     682      endif !iaergo_generic(1) .ne. 0
    648683
    649684! --------------------------------------------------------------------------
  • trunk/LMDZ.GENERIC/libf/phystd/aerosol_mod.F90

    r2297 r2803  
    2222! Auroral aerosols
    2323      integer :: iaero_aurora = 0
    24 !$OMP THREADPRIVATE(iaero_co2,iaero_h2o,iaero_dust,iaero_h2so4,noaero,iaero_back2lay,iaero_nh3,iaero_nlay,iaero_aurora)
     24! Generic aerosols
     25      integer, dimension(:),allocatable,save :: iaero_generic
     26!$OMP THREADPRIVATE(iaero_co2,iaero_h2o,iaero_dust,iaero_h2so4,noaero,iaero_back2lay,iaero_nh3,iaero_nlay,iaero_aurora,iaero_generic)
    2527     
    2628!==================================================================
  • trunk/LMDZ.GENERIC/libf/phystd/callcorrk.F90

    r2736 r2803  
    297297#endif
    298298
    299          if(naerkind.gt.4)then
    300             message='Code not general enough to deal with naerkind > 4 yet.'
    301             call abort_physic(subname,message,1)
    302          endif
     299         ! if(naerkind.gt.4)then
     300         !    message='Code not general enough to deal with naerkind > 4 yet.'
     301         !    call abort_physic(subname,message,1)
     302         ! endif
    303303         call su_aer_radii(ngrid,nlayer,reffrad,nueffrad)
    304304         
  • trunk/LMDZ.GENERIC/libf/phystd/callkeys_mod.F90

    r2721 r2803  
    5050      logical,save :: aerofixco2,aerofixh2o
    5151!$OMP THREADPRIVATE(aerofixco2,aerofixh2o)
     52      integer,save :: aerogeneric
     53!$OMP THREADPRIVATE(aerogeneric)
    5254      logical,save :: hydrology
    5355      logical,save :: sourceevol
  • trunk/LMDZ.GENERIC/libf/phystd/iniaerosol.F

    r2297 r2803  
    55      use aerosol_mod
    66      use callkeys_mod, only: aeroco2,aeroh2o,dusttau,aeroh2so4,
    7      &          aeroback2lay,aeronh3, nlayaero, aeronlay, aeroaurora
     7     &          aeroback2lay,aeronh3, nlayaero, aeronlay, aeroaurora,
     8     &      aerogeneric
    89
    910      IMPLICIT NONE
     
    2425      IF(.NOT.ALLOCATED(iaero_nlay)) ALLOCATE(iaero_nlay(nlayaero))
    2526      iaero_nlay(:) = 0
    26 
     27      ! Do the same for iaero_generic
     28      IF (.not. allocated(iaero_generic))
     29     &  allocate(iaero_generic(aerogeneric))
     30      iaero_generic(:)=0
    2731      ia=0
    2832      if (aeroco2) then
     
    7680      write(*,*) '--- Auroral aerosols = ', iaero_aurora
    7781
     82      if (aerogeneric .ne. 0) then
     83         do i=1,aerogeneric
     84            ia = ia+1
     85            iaero_generic(i) = ia
     86         enddo
     87      endif
     88      write(*,*)'--- Generic Condensable Aerosols',iaero_generic
    7889      write(*,*) '=== Number of aerosols= ', ia
    79      
     90
    8091! For the zero aerosol case, we currently make a dummy co2 aerosol which is zero everywhere.
    8192! (See aeropacity.F90 for how this works). A better solution would be to turn off the
  • trunk/LMDZ.GENERIC/libf/phystd/inifis_mod.F90

    r2721 r2803  
    655655     if (is_master) write(*,*)trim(rname)//": aeroaurora = ",aeroaurora
    656656
     657     if (is_master) write(*,*)trim(rname)//&
     658     ": Radiatively active Generic Condensable aerosols?"
     659     aerogeneric=0     ! default value
     660     call getin_p("aerogeneric",aerogeneric)
     661     if (is_master) write(*,*)trim(rname)//":aerogeneric",aerogeneric
     662
     663
    657664!=================================
    658665! TWOLAY scheme and NH3 cloudare left for retrocompatibility only,
  • trunk/LMDZ.GENERIC/libf/phystd/physiq_mod.F90

    r2802 r2803  
    495495         zdtlw(:,:) = 0.0
    496496
    497 
    498 !        Initialize aerosol indexes.
    499 !        ~~~~~~~~~~~~~~~~~~~~~~~~~~~
    500          call iniaerosol()
    501 
    502      
    503497!        Initialize tracer names, indexes and properties.
    504498!        ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
     
    512506            endif
    513507         endif
     508!        Initialize aerosol indexes.
     509!        ~~~~~~~~~~~~~~~~~~~~~~~~~~~
     510         call iniaerosol()
    514511
    515512#ifdef CPP_XIOS
     
    25372534            endif
    25382535
    2539             do ig=1,ngrid
    2540                if(tau_col(ig).gt.1.e3)then
    2541                   print*,'WARNING: tau_col=',tau_col(ig)
    2542                   print*,'at ig=',ig,'in PHYSIQ'
    2543                endif         
    2544             end do
     2536            ! do ig=1,ngrid
     2537            !    if(tau_col(ig).gt.1.e3)then
     2538            !       print*,'WARNING: tau_col=',tau_col(ig)
     2539            !       print*,'at ig=',ig,'in PHYSIQ'
     2540            !    endif         
     2541            ! end do
    25452542
    25462543            call writediagfi(ngrid,"tau_col","Total aerosol optical depth","[]",2,tau_col)
  • trunk/LMDZ.GENERIC/libf/phystd/radii_mod.F90

    r2340 r2803  
    2828!     -------
    2929!     Compute the effective radii of liquid and icy water particles
     30!     Jeremy Leconte (2012)
     31!     Extended to dust, CO2, NH3, 2-lay,Nlay,auroral aerosols by ??
     32!     Added Radiative Generic Condensable Aerosols effective radii
     33!     calculations  (Lucas Teinturier 2022)
    3034!
    3135!     Authors
     
    3741      use radinc_h, only: naerkind
    3842      use aerosol_mod, only: iaero_back2lay, iaero_co2, iaero_dust, &
    39                              iaero_h2o, iaero_h2so4, iaero_nh3, iaero_nlay, iaero_aurora
    40       use callkeys_mod, only: size_nh3_cloud, nlayaero, aeronlay_size, aeronlay_nueff
     43                             iaero_h2o, iaero_h2so4, iaero_nh3, iaero_nlay, iaero_aurora, &
     44                             iaero_generic
     45      use callkeys_mod, only: size_nh3_cloud, nlayaero, aeronlay_size, aeronlay_nueff,aerogeneric
    4146
    4247      Implicit none
     
    5156!$OMP THREADPRIVATE(firstcall)
    5257      integer :: iaer, ia   
    53      
    54       print*,'enter su_aer_radii'
    55           do iaer=1,naerkind
     58      real :: temprad !temporary variable for radius when we have Radiative Generic Condensable aerosols
     59      do iaer=1,naerkind
    5660!     these values will change once the microphysics gets to work
    5761!     UNLESS tracer=.false., in which case we should be working with
     
    6064!                |-> Done in th n-layer aerosol case (JVO 20)
    6165
    62             if(iaer.eq.iaero_co2)then ! CO2 ice
    63                reffrad(1:ngrid,1:nlayer,iaer) = 1.e-4
    64                nueffrad(1:ngrid,1:nlayer,iaer) = 0.1
     66         if(iaer.eq.iaero_co2)then ! CO2 ice
     67            reffrad(1:ngrid,1:nlayer,iaer) = 1.e-4
     68            nueffrad(1:ngrid,1:nlayer,iaer) = 0.1
     69         endif
     70
     71         if(iaer.eq.iaero_h2o)then ! H2O ice
     72            reffrad(1:ngrid,1:nlayer,iaer) = 1.e-5
     73            nueffrad(1:ngrid,1:nlayer,iaer) = 0.1
     74         endif
     75
     76         if(iaer.eq.iaero_dust)then ! dust
     77            reffrad(1:ngrid,1:nlayer,iaer) = 1.e-5
     78            nueffrad(1:ngrid,1:nlayer,iaer) = 0.1
     79         endif
     80 
     81         if(iaer.eq.iaero_h2so4)then ! H2O ice
     82            reffrad(1:ngrid,1:nlayer,iaer) = 1.e-6
     83            nueffrad(1:ngrid,1:nlayer,iaer) = 0.1
     84         endif
     85           
     86         if(iaer.eq.iaero_back2lay)then ! Two-layer aerosols
     87            reffrad(1:ngrid,1:nlayer,iaer) = 2.e-6
     88            nueffrad(1:ngrid,1:nlayer,iaer) = 0.1
     89         endif
     90
     91
     92              if(iaer.eq.iaero_nh3)then ! Nh3 cloud
     93            reffrad(1:ngrid,1:nlayer,iaer) = size_nh3_cloud
     94            nueffrad(1:ngrid,1:nlayer,iaer) = 0.1
     95         endif
     96
     97         do ia=1,nlayaero
     98            if(iaer.eq.iaero_nlay(ia))then ! N-layer aerosols
     99               reffrad(1:ngrid,1:nlayer,iaer) = aeronlay_size(ia)
     100               nueffrad(1:ngrid,1:nlayer,iaer) = aeronlay_nueff(ia)
    65101            endif
    66 
    67             if(iaer.eq.iaero_h2o)then ! H2O ice
    68                reffrad(1:ngrid,1:nlayer,iaer) = 1.e-5
    69                nueffrad(1:ngrid,1:nlayer,iaer) = 0.1
    70             endif
    71 
    72             if(iaer.eq.iaero_dust)then ! dust
    73                reffrad(1:ngrid,1:nlayer,iaer) = 1.e-5
    74                nueffrad(1:ngrid,1:nlayer,iaer) = 0.1
    75             endif
    76  
    77             if(iaer.eq.iaero_h2so4)then ! H2O ice
    78                reffrad(1:ngrid,1:nlayer,iaer) = 1.e-6
    79                nueffrad(1:ngrid,1:nlayer,iaer) = 0.1
    80             endif
    81            
    82             if(iaer.eq.iaero_back2lay)then ! Two-layer aerosols
    83                reffrad(1:ngrid,1:nlayer,iaer) = 2.e-6
    84                nueffrad(1:ngrid,1:nlayer,iaer) = 0.1
    85             endif
    86 
    87 
    88             if(iaer.eq.iaero_nh3)then ! Nh3 cloud
    89                reffrad(1:ngrid,1:nlayer,iaer) = size_nh3_cloud
    90                nueffrad(1:ngrid,1:nlayer,iaer) = 0.1
    91             endif
    92 
    93             do ia=1,nlayaero
    94               if(iaer.eq.iaero_nlay(ia))then ! N-layer aerosols
    95                  reffrad(1:ngrid,1:nlayer,iaer) = aeronlay_size(ia)
    96                  nueffrad(1:ngrid,1:nlayer,iaer) = aeronlay_nueff(ia)
    97               endif
    98             enddo
    99 
    100             if(iaer.eq.iaero_aurora)then ! Auroral aerosols
    101                reffrad(1:ngrid,1:nlayer,iaer) = 3.e-7
    102                nueffrad(1:ngrid,1:nlayer,iaer) = 0.1
    103             endif
    104 
    105 
    106             if(iaer.gt.5)then
    107                print*,'Error in callcorrk, naerkind is too high (>5).'
    108                print*,'The code still needs generalisation to arbitrary'
    109                print*,'aerosol kinds and number.'
    110                call abort
    111             endif
    112 
    113102         enddo
    114103
    115 
    116          if (radfixed) then
    117 
    118             write(*,*)"radius of H2O water particles:"
    119             rad_h2o=13. ! default value
    120             call getin_p("rad_h2o",rad_h2o)
    121             write(*,*)" rad_h2o = ",rad_h2o
    122 
    123             write(*,*)"radius of H2O ice particles:"
    124             rad_h2o_ice=35. ! default value
    125             call getin_p("rad_h2o_ice",rad_h2o_ice)
    126             write(*,*)" rad_h2o_ice = ",rad_h2o_ice
    127 
    128          else
    129 
    130             write(*,*)"Number mixing ratio of H2O water particles:"
    131             Nmix_h2o=1.e6 ! default value
    132             call getin_p("Nmix_h2o",Nmix_h2o)
    133             write(*,*)" Nmix_h2o = ",Nmix_h2o
    134 
    135             write(*,*)"Number mixing ratio of H2O ice particles:"
    136             Nmix_h2o_ice=Nmix_h2o ! default value
    137             call getin_p("Nmix_h2o_ice",Nmix_h2o_ice)
    138             write(*,*)" Nmix_h2o_ice = ",Nmix_h2o_ice
    139          endif
    140 
    141       print*,'exit su_aer_radii'
     104              if(iaer.eq.iaero_aurora)then ! Auroral aerosols
     105            reffrad(1:ngrid,1:nlayer,iaer) = 3.e-7
     106            nueffrad(1:ngrid,1:nlayer,iaer) = 0.1
     107         endif
     108
     109         do ia=1,aerogeneric     ! Radiative Generic Condensable aerosols
     110            if (iaer .eq. iaero_generic(ia)) then
     111               call generic_get_traceur_radius(iaer,temprad)
     112               reffrad(1:ngrid,1:nlayer,iaer)=temprad
     113               nueffrad(1:ngrid,1:nlayer,iaer) = 0.1
     114            endif
     115         enddo  ! generic radiative condensable aerosols
     116
     117         if(iaer.gt.5)then ! I think this could be commented out by this needs testing. It works with 5 aerosols (LT 2022)
     118            print*,'Error in callcorrk, naerkind is too high (>5).'
     119            print*,'The code still needs generalisation to arbitrary'
     120            print*,'aerosol kinds and number.'
     121            call abort
     122         endif
     123      enddo ! iaer=1,naerkind
     124
     125
     126      if (radfixed) then
     127
     128         write(*,*)"radius of H2O water particles:"
     129         rad_h2o=13. ! default value
     130         call getin_p("rad_h2o",rad_h2o)
     131         write(*,*)" rad_h2o = ",rad_h2o
     132
     133         write(*,*)"radius of H2O ice particles:"
     134         rad_h2o_ice=35. ! default value
     135         call getin_p("rad_h2o_ice",rad_h2o_ice)
     136         write(*,*)" rad_h2o_ice = ",rad_h2o_ice
     137
     138      else
     139
     140         write(*,*)"Number mixing ratio of H2O water particles:"
     141         Nmix_h2o=1.e6 ! default value
     142         call getin_p("Nmix_h2o",Nmix_h2o)
     143         write(*,*)" Nmix_h2o = ",Nmix_h2o
     144
     145         write(*,*)"Number mixing ratio of H2O ice particles:"
     146         Nmix_h2o_ice=Nmix_h2o ! default value
     147         call getin_p("Nmix_h2o_ice",Nmix_h2o_ice)
     148         write(*,*)" Nmix_h2o_ice = ",Nmix_h2o_ice
     149      endif
     150
    142151
    143152   end subroutine su_aer_radii
     
    381390!==================================================================
    382391
     392!==================================================================
     393   subroutine generic_get_traceur_radius(iaer,radii)
     394!==================================================================
     395!     Purpose
     396!     -------
     397!     Get the generic traceur radius inputted in modern_traceur.def
     398!      for the specie "specie"
     399!
     400!     Authors
     401!     -------
     402!     Lucas Teinturier (2022)
     403!
     404!==================================================================
     405      Use tracer_h, only: radius, noms, nqtot
     406      Implicit none
     407     
     408      integer, intent(in) :: iaer !index of aerosol which radius we're looking for
     409      real, intent(out) :: radii !radii(ngrid,nlayer,naerkind) ! Particle radii (m)
     410      integer :: iq !just a loop index
     411
     412      if (index(noms(1),'ice') .ne. 0) then !traceur.def is organized with _ice then _vap
     413         radii=radius(2*iaer-1)
     414      elseif (index(noms(1),'vap') .ne. 0) then !traceur.def is organized with _vap then _ice
     415         radii = radius(2*iaer)
     416      else
     417         print*,"radii_mod: tracer is neither ice or vap :",noms(1)
     418         print*,"Aborting"
     419         call abort
     420      endif
     421     
     422   end subroutine generic_get_traceur_radius
    383423
    384424end module radii_mod
  • trunk/LMDZ.GENERIC/libf/phystd/suaer_corrk.F90

    r2297 r2803  
    1313      use callkeys_mod, only: tplanet, optprop_back2lay_vis, optprop_back2lay_ir, &
    1414                              optprop_aeronlay_vis, optprop_aeronlay_ir,          &
    15                               aeronlay_lamref, nlayaero
    16 
     15                              aeronlay_lamref, nlayaero,aerogeneric
     16      use tracer_h, only: noms
    1717      implicit none
    1818
     
    119119
    120120!--------------------------------------------------------------
    121 
    122121      if (noaero) then
    123122        print*, 'naerkind= 0'
     
    233232! added by SG
    234233       endif
    235  
    236        
     234! the following was added by LT
     235       do ia=1,aerogeneric ! Read Radiative Generic Condensable Aerosols data
     236         if (iaer .eq. iaero_generic(ia)) then
     237            if (index(noms(2*iaer),'Fe') .ne. 0) then
     238               print*,"Reading Fe file"
     239               file_id(iaer,1)='Fe.ocst.txt'
     240               file_id(iaer,2)='Fe.ocst.txt'
     241               lamrefvis(iaer) = 1.0E-6 !random pick
     242               lamrefir(iaer) = 1.0E-6 !dummy but random pick
     243            else if (index(noms(2*iaer),'Mn') .ne. 0) then
     244               print*,"Reading MnS file"
     245               file_id(iaer,1)='MnS.ocst.txt'
     246               file_id(iaer,2)='MnS.ocst.txt'
     247               lamrefvis(iaer) = 1.0E-6 !random pick
     248               lamrefir(iaer) = 1.0E-6 !dummy but random pick   
     249            else if (index(noms(2*iaer),'Mg') .ne. 0) then 
     250               print*,"Reading Mg2SiO4 file"
     251               file_id(iaer,1)='Mg2SiO4.ocst.txt'
     252               file_id(iaer,2)='Mg2SiO4.ocst.txt'
     253               lamrefvis(iaer) = 1.0E-6 !random pick
     254               lamrefir(iaer) = 1.0E-6 !dummy but random pick 
     255            else if (index(noms(2*iaer),'Cr') .ne. 0) then
     256               print*,"Reading Cr file"
     257               file_id(iaer,1)='Cr.ocst.txt'
     258               file_id(iaer,2)='Cr.ocst.txt'
     259               lamrefvis(iaer) = 1.0E-6 !random pick
     260               lamrefir(iaer) = 1.0E-6 !dummy but random pick
     261            else
     262! If you want to add another specie, copy,paste & adapt the elseif block up here to your new specie (LT 2022)
     263               call abort_physic("suaer_corrk", "Unknown specie in generic condensable aerosols",1)
     264            endif
     265         endif
     266       enddo ! ia=1,aerogeneric
    237267      enddo
    238268
     
    415445
    416446
    417 
    418447!==================================================================
    419448!     2. AVERAGED PROPERTIES AND VARIABLE ASSIGNMENTS
     
    453482
    454483         ENDDO
    455 
    456484!==================================================================
    457485      CASE(2) domain            !      INFRARED DOMAIN (idomain=2)
    458486!==================================================================
    459 
    460487
    461488         DO isize=1,nsize(iaer,idomain) ! ----------------------------------
     
    517544      END DO                    ! Loop on idomain
    518545!========================================================================
    519 
    520546      RETURN
    521547
Note: See TracChangeset for help on using the changeset viewer.