Changeset 2704


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

fixed condensation_generic to conserve tracers mass. evap_generic is now useless

File:
1 edited

Legend:

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

    r2701 r2704  
    8787        do iq=1,nq
    8888                if((is_generic(iq)==1) .and. (index(noms(iq),"vap") .ne. 0)) then
    89 !                       Let's get the index of our tracers (we look for igcm _generic_vap and igcm_generic_ice)
    90                         ! tname_ice = trim(noms(iq)(1:len(tname_ice)-3))//"ice"
    91                         ! print*,trim(adjustl(trim(noms(iq))(9:len(trim(noms(iq)))-4))) !testing here, should go away
    92                         print*,noms(iq)(9:len(trim(noms(iq)))-4)
    93                         ! stop
    9489                        ! Hyp : vap tracer index comes right before ice tracer index in traceur.def
    9590                        igcm_generic_vap=iq
     
    9792                        ! Need to call specie_parameters of the considered specie
    9893                        call specie_parameters(noms(iq)(9:len(trim(noms(iq)))-4))
    99                         ! Evaporate generic clouds (all of them)
    10094                        Lcp=RLVTT/cpp !need to be init here, otherwise we don't know RLVTT yet
    101                         call evap_generic(ngrid,nlayer,nq,ptimestep,pt,pq,pdq,pdt,&
    102                                         igcm_generic_vap,igcm_generic_ice, &
    103                                         dqevap,dtevap,qevap,tevap)
    104                                 ! note: we use qevap but not tevap in largescale/moistadj
    105                                 ! otherwise is a big mess
    10695
    10796                        !  Vertical loop (from top to bottom)
    10897                        DO k = nlayer, 1, -1
    109                                 zt(1:ngrid)=pt(1:ngrid,k)+(pdt(1:ngrid,k)+dtevap(1:ngrid,k))*ptimestep
    110                                 zq(1:ngrid)=qevap(1:ngrid,k,igcm_generic_vap)
     98                                zt(1:ngrid)=pt(1:ngrid,k)+pdt(1:ngrid,k)*ptimestep
    11199
    112100                                ! Computes Psat and the partial condensation
     
    117105                                        !           zt(i)=15.   ! check too low temperatures
    118106                                        endif
    119                                         call Psat_generic(zt(i),local_p,metallicity,psat_tmp,zqs_temp) !useless ?
    120                                         zqs(i)=zqs_temp ! useless ?
    121                                         zt(i)=pt(i,k)+pdt(i,k)*ptimestep
    122107                                        zx_q(i) = pq(i,k,igcm_generic_vap)+pdq(i,k,igcm_generic_vap)*ptimestep
    123                                         ! dqevap(i,k)=0.
    124108                                        ! iterative process to stabilize the scheme when large water amounts JL12
    125109                                        zcond(i) = 0.0d0
     
    141125
    142126                                !Tendances de t et q
    143                                 pdqvaplsc(1:ngrid,k,igcm_generic_vap)  = dqevap(1:ngrid,k) - zcond(1:ngrid)
     127                                pdqvaplsc(1:ngrid,k,igcm_generic_vap)  = - zcond(1:ngrid)
    144128                                pdqliqlsc(1:ngrid,k,igcm_generic_ice) = - pdqvaplsc(1:ngrid,k,igcm_generic_vap)
    145129                                pdtlsc(1:ngrid,k)  = pdtlsc(1:ngrid,k) + pdqliqlsc(1:ngrid,k,igcm_generic_ice)*Lcp
Note: See TracChangeset for help on using the changeset viewer.