Changeset 2706


Ignore:
Timestamp:
Jun 21, 2022, 11:05:38 AM (3 years ago)
Author:
aslmd
Message:

is_condensable has replaced is_generic, some changes in getting back tracers names

Location:
trunk/LMDZ.GENERIC
Files:
7 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.GENERIC/deftank/traceur.def.modern

    r2693 r2706  
    5252#               is_recomb_qotf ! 1 if tracer recombination is done on-the-fly, else 0
    5353#                                (if 1, must have is_recomb_qset=0)
    54 #               is_generic     ! 1 if tracer is generic, else 0
     54#               is_condensable  ! 1 if tracer is_condensable, else 0
    5555#               radius          ! radius value in meter
    5656#               rho             ! density in kg/m3
  • trunk/LMDZ.GENERIC/libf/phystd/condensation_generic_mod.F90

    r2704 r2706  
    88        use ioipsl_getin_p_mod, only: getin_p !-> to get the metallicity
    99        use generic_cloud_common_h, only : RLVTT, cpp, &
    10         Psat_generic,Lcpdqsat_generic,specie_parameters
     10        Psat_generic,Lcpdqsat_generic,specie_parameters,specie_parameters_table
    1111        USE tracer_h
    1212        IMPLICIT none
     
    8686        ! Let's loop on tracers
    8787        do iq=1,nq
    88                 if((is_generic(iq)==1) .and. (index(noms(iq),"vap") .ne. 0)) then
    89                         ! Hyp : vap tracer index comes right before ice tracer index in traceur.def
     88                if((is_condensable(iq)==1) .and. (index(noms(iq),"vap") .ne. 0)) then
     89                        write(*,*) "the specie ", noms(iq)," is condensable, for generic condensation"
     90!                       Let's get the index of our tracers (we look for igcm _generic_vap and igcm_generic_ice)
     91                        ! tname_ice = trim(noms(iq)(1:len(tname_ice)-3))//"ice"
     92                        ! print*,trim(adjustl(trim(noms(iq))(9:len(trim(noms(iq)))-4))) !testing here, should go away
     93                        print*,noms(iq)(9:len(trim(noms(iq)))-4)
     94                        ! stop
    9095                        igcm_generic_vap=iq
    91                         igcm_generic_ice = iq+1
     96
     97                        ! We look for the corresponding ice traceur (before or after in the list of traceurs, maybe could be generalised to the whole list)
     98                        if ((noms(iq)(9:len(trim(noms(iq)))-4) .eq. noms(iq+1)(9:len(trim(noms(iq+1)))-4)) .and. (index(noms(iq+1),"ice") .ne. 0)) then
     99                                igcm_generic_ice = iq+1
     100                        else if ((iq .gt. 1) .and. (noms(iq)(9:len(trim(noms(iq)))-4) .eq. noms(iq-1)(9:len(trim(noms(iq-1)))-4)) .and. (index(noms(iq-1),"ice") .ne. 0)) then
     101                                igcm_generic_ice = iq-1
     102                        else
     103                                write(*,*) "ERROR : You set a vap traceur but you forgot to set the corresponding ice traceur"
     104                        endif
     105
    92106                        ! Need to call specie_parameters of the considered specie
    93107                        call specie_parameters(noms(iq)(9:len(trim(noms(iq)))-4))
     108                        !call specie_parameters_table(noms(iq)(9:len(trim(noms(iq)))-4))
     109                       
    94110                        Lcp=RLVTT/cpp !need to be init here, otherwise we don't know RLVTT yet
    95111
     
    130146
    131147                        Enddo ! k= nlayer, 1, -1
    132                 endif !(is_generic(iq)==1) .and. (index(noms(iq),"vap") .ne. 0))
     148                endif !(is_condensable(iq)==1) .and. (index(noms(iq),"vap") .ne. 0))
    133149        enddo ! iq=1,nq
    134150
  • trunk/LMDZ.GENERIC/libf/phystd/dyn1d/rcm1d.F

    r2699 r2706  
    66      use mod_grid_phy_lmdz, only : regular_lonlat
    77      use infotrac, only: nqtot, tname
    8       use tracer_h, only: noms
     8      use tracer_h, only: noms, is_condensable
    99      use surfdat_h, only: albedodat, phisfi, dryness, watercaptag,
    1010     &                     zmea, zstd, zsig, zgam, zthe,
  • trunk/LMDZ.GENERIC/libf/phystd/generic_cloud_common_h.F90

    r2705 r2706  
    9494        epsi = m/mugaz
    9595    end subroutine specie_parameters
     96
     97
     98    subroutine specie_parameters_table(specname)
     99       
     100        implicit none
     101    !============================================================================
     102    !   Load the adequate set of parameters for specname
     103    !   From a table of traceurs
     104    !============================================================================
     105
     106        character(*), intent(in) :: specname
     107        integer k
     108        character(len=500):: table_traceurs_line ! table_traceurs_line lines with parameters
     109
     110        open(117,file='table_traceurs',form='formatted',status='old')
     111
     112        read(117,'(A)') table_traceurs_line
     113
     114        do
     115            read(117,'(A)') table_traceurs_line
     116           
     117            if (index(table_traceurs_line,specname) /= 0) then
     118
     119                write(*,*) table_traceurs_line
     120               
     121                if (index(table_traceurs_line,'deltavapH='   ) /= 0) then
     122                    read(table_traceurs_line(index(table_traceurs_line,'deltavapH=')+len('deltavapH='):),*) delta_vapH
     123                    write(*,*) 'delta_vapH ', delta_vapH
     124                end if
     125                if (index(table_traceurs_line,'Tref='   ) /= 0) then
     126                    read(table_traceurs_line(index(table_traceurs_line,'Tref=')+len('Tref='):),*) Tref
     127                end if
     128                if (index(table_traceurs_line,'Pref='   ) /= 0) then
     129                    read(table_traceurs_line(index(table_traceurs_line,'Pref=')+len('Pref='):),*) Pref
     130                end if
     131                if (index(table_traceurs_line,'mass='   ) /= 0) then
     132                    read(table_traceurs_line(index(table_traceurs_line,'mass=')+len('mass='):),*) m
     133                end if
     134                if (index(table_traceurs_line,'metallicity_coeff='   ) /= 0) then
     135                   read(table_traceurs_line(index(table_traceurs_line,'metallicity_coeff=')+len('metallicity_coeff='):),*) metallicity_coeff
     136                end if
     137            end if
     138       
     139        end do
     140        ! RLVTT
     141
     142        !if (is_master)
     143        close(117)
     144
     145    end subroutine specie_parameters_table
    96146
    97147    subroutine Psat_generic(T,p,metallicity,psat,qsat)
  • trunk/LMDZ.GENERIC/libf/phystd/initracer.F90

    r2703 r2706  
    116116         ALLOCATE(is_recomb_qotf(nqtot))
    117117       ENDIF
    118        IF (.NOT. allocated(is_generic)) allocate(is_generic(nq)) !LT
     118       IF (.NOT. allocated(is_condensable)) allocate(is_condensable(nq)) !LT
    119119       !! initialization
    120120       alpha_lift(:)     = 0.
     
    133133       radius(:)=0.
    134134       qext(:)=0.
    135        is_generic(:)= 0 !LT
     135       is_condensable(:)= 0 !LT
    136136       rho_q(:) = 0. !need to be init here if we want to read it from modern traceur with get_tracdat
    137137
     
    400400
    401401      enddo ! of do iq=1,nq
     402
     403      ! 3. find condensable traceurs different from h2o and co2
     404      do iq=1,nq
     405        if ((index(noms(iq),"vap") .ne. 0) .and. (index(noms(iq),"h2o") .eq. 0) .and. (index(noms(iq),"co2") .eq. 0)) then
     406          count=count+1
     407        endif
     408        if ((index(noms(iq),"ice") .ne. 0) .and. (index(noms(iq),"h2o") .eq. 0) .and. (index(noms(iq),"co2") .eq. 0)) then
     409          count=count+1
     410        endif
     411
     412      enddo ! of do iq=1,nq
    402413     
    403414      ! check that we identified all tracers:
     
    430441
    431442      ! Calculate number of generic tracers
    432       ngt = sum(is_generic)
     443      ngt = sum(is_condensable)
    433444      write(*,*) 'Number of generic tracer is  ngt = ',ngt
    434445!     Processing modern traceur options
     
    608619              is_recomb_qotf(iq)
    609620          end if
    610           !option is_generic (LT)
    611           if (index(tracline,'is_generic=') /=0) then
    612             read(tracline(index(tracline,'is_generic=') &
    613               +len('is_generic='):),*) is_generic(iq)
     621          !option is_condensable (LT)
     622          if (index(tracline,'is_condensable=') /=0) then
     623            read(tracline(index(tracline,'is_condensable=') &
     624              +len('is_condensable='):),*) is_condensable(iq)
    614625            write(*,*) ' Parameter value (traceur.def) :'// &
    615               ' is_generic=', is_generic(iq)       
     626              ' is_condensable=', is_condensable(iq)       
    616627          else
    617628              write(*,*) ' Parameter value (default)     :'// &
    618               ' is_generic=', is_generic(iq)       
     629              ' is_condensable=', is_condensable(iq)       
    619630          endif
    620631          !option radius
  • trunk/LMDZ.GENERIC/libf/phystd/physiq_mod.F90

    r2701 r2706  
    3232                          alpha_lift, alpha_devil, qextrhor, &
    3333                          igcm_h2o_ice, igcm_h2o_vap, igcm_dustbin, &
    34                           igcm_co2_ice, nesp, is_chim, is_generic
     34                          igcm_co2_ice, nesp, is_chim, is_condensable
    3535      use time_phylmdz_mod, only: ecritphy, iphysiq, nday
    3636      use phyetat0_mod, only: phyetat0
  • trunk/LMDZ.GENERIC/libf/phystd/tracer_h.F90

    r2690 r2706  
    4040       integer, save, allocatable :: is_recomb_qotf(:) ! 1 if tracer recombination is done on-the-fly, else 0 (if 1, must have is_recomb_qset=0)
    4141!$OMP THREADPRIVATE(is_recomb,is_recomb_qset,is_recomb_qotf)
    42        integer, save, allocatable :: is_generic(:)      ! 1 if tracer is generic, else 0 (added LT)
    43 !$OMP THREADPRIVATE(is_generic)   !also added by LT
     42       integer, save, allocatable :: is_condensable(:)      ! 1 if tracer is generic, else 0 (added LT)
     43!$OMP THREADPRIVATE(is_condensable)   !also added by LT
    4444! tracer indexes: these are initialized in initracer and should be 0 if the
    4545!                 corresponding tracer does not exist
Note: See TracChangeset for help on using the changeset viewer.