Changeset 2694


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

fixed uninitialized variable Lcp in largescale_generic+ typo in index of tendencies after largescale_generic in physiq_mod

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

Legend:

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

    r2690 r2694  
    3636    real zlvdcp
    3737    integer l,ig ! Used for loops on the layers and grid and the tracers respectively
     38
    3839    !   Evaporate all the ice from the tracer
    3940    zlvdcp = RLVTT/cpp ! RLVTT is the latent heat of vaporization (comes from genericcommon_h attention)
     41 
    4042    DO l=1,nlayer
    4143        DO ig=1,ngrid
     
    4547                    +pdq(ig,l,igcm_generic_ice)*dtime
    4648            tevap(ig,l) = pt(ig,l)+pdt(ig,l)*dtime
    47             dqevap(ig,l) = max(0.,qevap(ig,l,igcm_generic_ice)/dtime)
     49            dqevap(ig,l) = max(0.,qevap(ig,l,igcm_generic_ice))/dtime
    4850            dtevap(ig,l) = -dqevap(ig,l)*zlvdcp
    4951            qevap(ig,l,igcm_generic_vap) = qevap(ig,l,igcm_generic_vap) &
     
    5456    ENDDO ! nlayer
    5557
    56 
    5758end subroutine evap_generic
    5859
  • trunk/LMDZ.GENERIC/libf/phystd/genericcommon_h.F90

    r2690 r2694  
    3434            Tref = 2303.0 !K
    3535            Pref =   1.0E5 !in Pa
     36            m = 24.3055
    3637            RLVTT = delta_vapH/(m/1000.)
    3738            metallicity_coeff=1.0*log(10.0)
    38             m = 24.3055
     39           
    3940        else if (trim(specname) .eq. "Na") then
    4041            print*,"Loading data for Na"
     
    4243            Tref = 1624.0
    4344            Pref =  1.0E5
     45            m = 22.9898
    4446            RLVTT = delta_vapH/(m/1000.)
    4547            metallicity_coeff=0.5*log(10.0)
    46             m = 22.9898
     48           
    4749        else if (trim(specname) .eq. "Fe") then
    4850            print*,"Loading data for Fe"
     
    5052            Tref = 2903.0
    5153            Pref  = 1.0E5
     54            m = 55.8450
    5255            RLVTT = delta_vapH/(m/1000.)
    5356            metallicity_coeff=0.0*log(10.0)
    54             m = 55.8450
     57           
    5558        else if (trim(specname) .eq. "Cr") then
    5659            print*,"Loading data for Cr"
     
    5861            Tref = 2749.0
    5962            Pref  = 1.0E5
     63            m = 51.9961
    6064            RLVTT = delta_vapH/(m/1000.)
    6165            metallicity_coeff=0.0*log(10.0)
    62             m = 51.9961
     66           
    6367        else if (trim(specname) .eq. "KCl") then
    6468            print*,"Loading data for KCl"
     
    6872            metallicity_coeff=0.0*log(10.0)
    6973            m = 74.5498
     74            RLVTT = delta_vapH/(m/1000.)
    7075        else if (trim(specname) .eq. "Mn") then
    7176            print*,"Loading data for Mn"
     
    7580            metallicity_coeff=1.0*log(10.0)
    7681            m = 54.9380
     82            RLVTT = delta_vapH/(m/1000.)
    7783        else if (trim(specname) .eq. "Zn") then
    7884            print*,"Loading data for Zn"
     
    8288            metallicity_coeff=1.0*log(10.0)
    8389            m = 65.38
     90            RLVTT = delta_vapH/(m/1000.)
    8491        else
    8592            print*,"Unknow species (not in Mg, Fe, Na, KCl, Cr, Mn or Zn)"
  • trunk/LMDZ.GENERIC/libf/phystd/largescale_generic_mod.F90

    r2690 r2694  
    7878        pdqvaplsc(1:ngrid,1:nlayer)  = 0.0
    7979        pdqliqlsc(1:ngrid,1:nlayer) = 0.0
    80         Lcp=RLVTT/cpp
    8180        igcm_generic_vap = -1
    8281        igcm_generic_ice = -1
     
    9392        call specie_parameters(trim(adjustl(tname_vap(9:len(tname_ice)-4))))
    9493        ! Evaporate generic clouds (all of them)
     94        Lcp=RLVTT/cpp !need to be init here, otherwise we don't know RLVTT yet
    9595        call evap_generic(ngrid,nlayer,nq,ptimestep,pt,pq,pdq,pdt,&
    9696                igcm_generic_vap,igcm_generic_ice, &
  • trunk/LMDZ.GENERIC/libf/phystd/physiq_mod.F90

    r2692 r2694  
    15561556                                       dqvaplscale_generic,dqcldlscale_generic)
    15571557               pdt(1:ngrid,1:nlayer) = pdt(1:ngrid,1:nlayer)+dtlscale_generic(1:ngrid,1:nlayer)
    1558                pdq(1:ngrid,1:nlayer,ig) = pdq(1:ngrid,1:nlayer,ig)+dtlscale_generic(1:ngrid,1:nlayer)
    1559                pdq(1:ngrid,1:nlayer,ig+1) = pdq(1:ngrid,1:nlayer,ig+1)+dqcldlscale_generic(1:ngrid,1:nlayer)
     1558               pdq(1:ngrid,1:nlayer,iq) = pdq(1:ngrid,1:nlayer,iq)+dtlscale_generic(1:ngrid,1:nlayer)
     1559               pdq(1:ngrid,1:nlayer,iq+1) = pdq(1:ngrid,1:nlayer,iq+1)+dqcldlscale_generic(1:ngrid,1:nlayer)
    15601560               !!!!huge assumption in the last line : ice tracer has index +1 of vap tracer
    15611561            endif ! (is_generic(iq)==1 .and. (index(noms(iq),"vap")==1))
Note: See TracChangeset for help on using the changeset viewer.