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

major changes in largescale_generic and how it's called by physiq_mod

File:
1 edited

Legend:

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

    r2696 r2700  
    55   
    66    subroutine largescale_generic(ngrid,nlayer,nq,ptimestep, pplev, pplay,    &
    7                 pt, pq, pdt, pdq,tname_vap, pdtlsc, pdqvaplsc, pdqliqlsc)
     7                pt, pq, pdt, pdq, pdtlsc, pdqvaplsc, pdqliqlsc)
    88        use ioipsl_getin_p_mod, only: getin_p !-> to get the metallicity
    99        use genericcommon_h, only : RLVTT, cpp, &
     
    3838        REAL, intent(in) :: pdt(ngrid,nlayer)     ! physical temperature tendency (K/s)
    3939        REAL, intent(in) :: pdq(ngrid,nlayer,nq)  ! physical tracer tendency (K/s)
    40         CHARACTER(*), intent(in) :: tname_vap     ! name of the tracer we consider. BY convention, it ends with _vap !!!
     40        ! CHARACTER(*), intent(in) :: tname_vap     ! name of the tracer we consider. BY convention, it ends with _vap !!!
    4141        REAL, intent(out) :: pdtlsc(ngrid,nlayer)  ! incrementation de la temperature (K)
    42         REAL, intent(out) :: pdqvaplsc(ngrid,nlayer) ! incrementation de la vapeur du traceur
    43         REAL, intent(out) :: pdqliqlsc(ngrid,nlayer) ! incrementation du traceur liquide
     42        REAL, intent(out) :: pdqvaplsc(ngrid,nlayer,nq) ! incrementation de la vapeur du traceur
     43        REAL, intent(out) :: pdqliqlsc(ngrid,nlayer,nq) ! incrementation du traceur liquide
    4444
    4545!       Options :
     
    5858        real zt(ngrid),local_p,psat_tmp,dlnpsat_tmp,Lcp,zqs_temp,zdqs
    5959        integer igcm_generic_vap, igcm_generic_ice! index of the vap and ice of generic_tracer
    60         CHARACTER(len=len(trim(tname_vap))) :: tname_ice
     60        ! CHARACTER(len=*) :: tname_ice
    6161        ! evaporation calculations
    6262        REAL dqevap(ngrid,nlayer),dtevap(ngrid,nlayer)     
     
    7474                firstcall = .false.
    7575        ENDIF
    76 !       Initialisation of outputs and Lcp
     76!       Initialisation of outputs and local variables
    7777        pdtlsc(1:ngrid,1:nlayer)  = 0.0
    78         pdqvaplsc(1:ngrid,1:nlayer)  = 0.0
    79         pdqliqlsc(1:ngrid,1:nlayer) = 0.0
     78        pdqvaplsc(1:ngrid,1:nlayer,1:nq)  = 0.0
     79        pdqliqlsc(1:ngrid,1:nlayer,1:nq) = 0.0
     80        dqevap(1:ngrid,1:nlayer)=0.0
     81        dtevap(1:ngrid,1:nlayer)=0.0
     82        qevap(1:ngrid,1:nlayer,1:nq)=0.0
     83        tevap(1:ngrid,1:nlayer)=0.0
    8084        igcm_generic_vap = -1
    8185        igcm_generic_ice = -1
     86        ! Let's loop on tracers
     87        do iq=1,nq
     88                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
     94                        ! Hyp : vap tracer index comes right before ice tracer index in traceur.def
     95                        igcm_generic_vap=iq
     96                        igcm_generic_ice = iq+1
     97                        ! Need to call specie_parameters of the considered specie
     98                        call specie_parameters(noms(iq)(9:len(trim(noms(iq)))-4))
     99                        ! Evaporate generic clouds (all of them)
     100                        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
    82106
    83 !       Let's get the index of our tracers (we look for igcm _generic_vap and igcm_generic_ice)
    84         tname_ice = trim(tname_vap(1:len(tname_ice)-3))//"ice"
    85         print*,tname_ice !testing here, should go away
    86         do iq=1,nq
    87             if (noms(iq) .eq. tname_vap) igcm_generic_vap = iq
    88             if (noms(iq) .eq. tname_ice) igcm_generic_ice = iq
    89             if ((igcm_generic_vap .ne. -1) .and. (igcm_generic_ice .ne. -1)) exit
    90         Enddo
    91         ! Need to call specie_parameters of the considered specie
    92         call specie_parameters(trim(adjustl(tname_vap(9:len(tname_ice)-4))))
    93         ! Evaporate generic clouds (all of them)
    94         Lcp=RLVTT/cpp !need to be init here, otherwise we don't know RLVTT yet
    95         call evap_generic(ngrid,nlayer,nq,ptimestep,pt,pq,pdq,pdt,&
    96                 igcm_generic_vap,igcm_generic_ice, &
    97                 dqevap,dtevap,qevap,tevap)
    98         ! note: we use qevap but not tevap in largescale/moistadj
    99         ! otherwise is a big mess
     107                        !  Vertical loop (from top to bottom)
     108                        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)
    100111
    101         !  Vertical loop (from top to bottom)
    102         DO k = nlayer, 1, -1
    103             zt(1:ngrid)=pt(1:ngrid,k)+(pdt(1:ngrid,k)+dtevap(1:ngrid,k))*ptimestep
    104             zq(1:ngrid)=qevap(1:ngrid,k,igcm_generic_vap)
     112                                ! Computes Psat and the partial condensation
     113                                DO i = 1, ngrid
     114                                        local_p=pplay(i,k)
     115                                        if(zt(i).le.15.) then
     116                                                print*,'in lsc',i,k,zt(i)
     117                                        !           zt(i)=15.   ! check too low temperatures
     118                                        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
     122                                        zx_q(i) = pq(i,k,igcm_generic_vap)+pdq(i,k,igcm_generic_vap)*ptimestep
     123                                        ! dqevap(i,k)=0.
     124                                        ! iterative process to stabilize the scheme when large water amounts JL12
     125                                        zcond(i) = 0.0d0
     126                                        Do nn=1,nitermax 
     127                                                call Psat_generic(zt(i),local_p,metallicity,psat_tmp,zqs_temp)
     128                                                zqs(i)=zqs_temp
     129                                                call Lcpdqsat_generic(zt(i),local_p,psat_tmp,zqs_temp,zdqs,dlnpsat_tmp)
     130                                                zcond_iter = alpha*(zx_q(i)-zqs(i))/(1.d0+zdqs)
     131                                                !zcond can be negative here
     132                                                zx_q(i) = zx_q(i) - zcond_iter
     133                                                zcond(i) = zcond(i) + zcond_iter
     134                                                zt(i) = zt(i) + zcond_iter*Lcp
     135                                                if (ABS(zcond_iter/alpha/zqs(i)).lt.qthreshold) exit
     136                                                if (nn.eq.nitermax) print*,'itermax in largescale'
     137                                        End do ! niter
     138                                        zcond(i)=MAX(zcond(i),-(pq(i,k,igcm_generic_ice)+pdq(i,k,igcm_generic_ice)*ptimestep))
     139                                        zcond(i) = zcond(i)/ptimestep
     140                                ENDDO ! i=1,ngrid
    105141
    106 !           Computes Psat and the partial condensation
    107             DO i = 1, ngrid
    108                 local_p=pplay(i,k)
    109                 if(zt(i).le.15.) then
    110                     print*,'in lsc',i,k,zt(i)
    111                      !      zt(i)=15.   ! check too low temperatures
    112                 endif
    113                 call Psat_generic(zt(i),local_p,metallicity,psat_tmp,zqs_temp) !useless ?
    114                 zqs(i)=zqs_temp ! useless ?
    115                 zt(i)=pt(i,k)+pdt(i,k)*ptimestep
    116                 zx_q(i) = pq(i,k,igcm_generic_vap)+pdq(i,k,igcm_generic_vap)*ptimestep
    117                 ! dqevap(i,k)=0.
    118         !           iterative process to stabilize the scheme when large water amounts JL12
    119                 zcond(i) = 0.0d0
    120                 Do nn=1,nitermax 
    121                     call Psat_generic(zt(i),local_p,metallicity,psat_tmp,zqs_temp)
    122                     zqs(i)=zqs_temp
    123                     call Lcpdqsat_generic(zt(i),local_p,psat_tmp,zqs_temp,zdqs,dlnpsat_tmp)
    124                     zcond_iter = alpha*(zx_q(i)-zqs(i))/(1.d0+zdqs)
    125                           !zcond can be negative here
    126                     zx_q(i) = zx_q(i) - zcond_iter
    127                     zcond(i) = zcond(i) + zcond_iter
    128                     zt(i) = zt(i) + zcond_iter*Lcp
    129                     if (ABS(zcond_iter/alpha/zqs(i)).lt.qthreshold) exit
    130                     if (nn.eq.nitermax) print*,'itermax in largescale'
    131                 End do ! niter
    132                 zcond(i)=MAX(zcond(i),-(pq(i,k,igcm_generic_ice)+pdq(i,k,igcm_generic_ice)*ptimestep))
    133                 zcond(i) = zcond(i)/ptimestep
    134             ENDDO ! i=1,ngrid
     142                                !Tendances de t et q
     143                                pdqvaplsc(1:ngrid,k,igcm_generic_vap)  = dqevap(1:ngrid,k) - zcond(1:ngrid)
     144                                pdqliqlsc(1:ngrid,k,igcm_generic_ice) = - pdqvaplsc(1:ngrid,k,igcm_generic_vap)
     145                                pdtlsc(1:ngrid,k)  = pdtlsc(1:ngrid,k) + pdqliqlsc(1:ngrid,k,igcm_generic_ice)*Lcp
    135146
    136 !           Tendances de t et q
    137             pdqvaplsc(1:ngrid,k)  = dqevap(1:ngrid,k) - zcond(1:ngrid)
    138             pdqliqlsc(1:ngrid,k) = - pdqvaplsc(1:ngrid,k)
    139             pdtlsc(1:ngrid,k)  = pdqliqlsc(1:ngrid,k)*Lcp
    140 
    141         Enddo ! k= nlayer, 1, -1
    142         print*, "pdqvaplsc = ",pdqvaplsc(:,:)
    143         print*,"pdtlsc = ",pdtlsc(:,:)
     147                        Enddo ! k= nlayer, 1, -1
     148                endif !(is_generic(iq)==1) .and. (index(noms(iq),"vap") .ne. 0))
     149        enddo ! iq=1,nq
     150        ! print*, "pdqvaplsc = ",pdqvaplsc(:,:)
     151        ! print*,"pdtlsc = ",pdtlsc(:,:)
    144152    end subroutine largescale_generic
    145153end module largescale_generic_mod
Note: See TracChangeset for help on using the changeset viewer.