Ignore:
Timestamp:
Mar 29, 2016, 3:39:55 PM (9 years ago)
Author:
fhourdin
Message:

Introduction d'une dependance epmax=f(Cape) sur proposition de Camille Risi

File:
1 edited

Legend:

Unmodified
Added
Removed
  • LMDZ5/trunk/libf/phylmd/cv30_routines.F90

    r2311 r2481  
    839839    q, qs, gz, p, h, tv, lv, pbase, buoybase, plcl, inb, tp, tvp, clw, hp, &
    840840    ep, sigp, buoy)
     841    ! epmax_cape: ajout arguments
    841842  IMPLICIT NONE
    842843
     
    12421243  REAL dtmin(nloc, nd), sigold(nloc, nd)
    12431244
    1244 
    12451245  ! -------------------------------------------------------
    12461246  ! -- Initialization
     
    13481348
    13491349  ! the interval on which cape is computed starts at pbase :
    1350 
    13511350  DO k = 1, nl
    13521351    DO i = 1, ncum
     
    31463145    vprecip, evap, ep, sig, w0, ft, fq, fu, fv, ftra, inb, ma, upwd, dnwd, &
    31473146    dnwd0, qcondc, wd, cape, da, phi, mp, phi2, d1a, dam, sij, elij, clw, &
    3148     epmlmmm, eplamm, wdtraina, wdtrainm, iflag1, precip1, vprecip1, evap1, &
     3147    epmlmmm, eplamm, wdtraina, wdtrainm,epmax_diag, iflag1, precip1, vprecip1, evap1, &
    31493148    ep1, sig1, w01, ft1, fq1, fu1, fv1, ftra1, inb1, ma1, upwd1, dnwd1, &
    31503149    dnwd01, qcondc1, wd1, cape1, da1, phi1, mp1, phi21, d1a1, dam1, sij1, &
    3151     elij1, clw1, epmlmmm1, eplamm1, wdtraina1, wdtrainm1)
     3150    elij1, clw1, epmlmmm1, eplamm1, wdtraina1, wdtrainm1,epmax_diag1) ! epmax_cape
    31523151  IMPLICIT NONE
    31533152
     
    31703169  REAL wd(nloc), cape(nloc)
    31713170  REAL da(nloc, nd), phi(nloc, nd, nd), mp(nloc, nd)
     3171  REAL epmax_diag(nloc) ! epmax_cape
    31723172  ! RomP >>>
    31733173  REAL phi2(nloc, nd, nd)
     
    31933193  REAL wd1(nloc), cape1(nloc)
    31943194  REAL da1(nloc, nd), phi1(nloc, nd, nd), mp1(nloc, nd)
     3195  REAL epmax_diag1(len) ! epmax_cape
    31953196  ! RomP >>>
    31963197  REAL phi21(len, nd, nd)
     
    32113212    inb1(idcum(i)) = inb(i)
    32123213    cape1(idcum(i)) = cape(i)
     3214    epmax_diag1(idcum(i))=epmax_diag(i) ! epmax_cape
    32133215  END DO
    32143216
     
    32693271END SUBROUTINE cv30_uncompress
    32703272
     3273        subroutine cv30_epmax_fn_cape(nloc,ncum,nd &
     3274                ,cape,ep,hp,icb,inb,clw,nk,t,h,lv &
     3275                ,epmax_diag)
     3276        implicit none
     3277
     3278        ! On fait varier epmax en fn de la cape
     3279        ! Il faut donc recalculer ep, et hp qui a déjà été calculé et
     3280        ! qui en dépend
     3281        ! Toutes les autres variables fn de ep sont calculées plus bas.
     3282
     3283#include "cvthermo.h"
     3284#include "cv30param.h"
     3285#include "conema3.h"
     3286
     3287! inputs:
     3288      integer ncum, nd, nloc
     3289      integer icb(nloc), inb(nloc)
     3290      real cape(nloc)
     3291      real clw(nloc,nd),lv(nloc,nd),t(nloc,nd),h(nloc,nd)
     3292      integer nk(nloc)
     3293! inouts:
     3294      real ep(nloc,nd)
     3295      real hp(nloc,nd)
     3296! outputs ou local
     3297      real epmax_diag(nloc)
     3298! locals
     3299      integer i,k   
     3300      real hp_bak(nloc,nd)
     3301
     3302        ! on recalcule ep et hp
     3303       
     3304        if (coef_epmax_cape.gt.1e-12) then
     3305        do i=1,ncum
     3306           epmax_diag(i)=epmax-coef_epmax_cape*sqrt(cape(i))
     3307           do k=1,nl
     3308                ep(i,k)=ep(i,k)/epmax*epmax_diag(i)
     3309                ep(i,k)=amax1(ep(i,k),0.0)
     3310                ep(i,k)=amin1(ep(i,k),epmax_diag(i))
     3311           enddo
     3312        enddo
     3313
     3314! On recalcule hp:
     3315      do k=1,nl
     3316        do i=1,ncum
     3317          hp_bak(i,k)=hp(i,k)
     3318        enddo
     3319      enddo
     3320      do k=1,nlp
     3321        do i=1,ncum
     3322          hp(i,k)=h(i,k)
     3323        enddo
     3324      enddo
     3325      do k=minorig+1,nl
     3326       do i=1,ncum
     3327        if((k.ge.icb(i)).and.(k.le.inb(i)))then
     3328          hp(i,k)=h(i,nk(i))+(lv(i,k)+(cpd-cpv)*t(i,k))*ep(i,k)*clw(i,k)
     3329        endif
     3330       enddo
     3331      enddo !do k=minorig+1,n
     3332!     write(*,*) 'cv30_routines 6218: hp(1,20)=',hp(1,20)
     3333      do i=1,ncum 
     3334       do k=1,nl
     3335        if (abs(hp_bak(i,k)-hp(i,k)).gt.0.01) then
     3336           write(*,*) 'i,k=',i,k
     3337           write(*,*) 'coef_epmax_cape=',coef_epmax_cape
     3338           write(*,*) 'epmax_diag(i)=',epmax_diag(i)
     3339           write(*,*) 'ep(i,k)=',ep(i,k)
     3340           write(*,*) 'hp(i,k)=',hp(i,k)
     3341           write(*,*) 'hp_bak(i,k)=',hp_bak(i,k)
     3342           write(*,*) 'h(i,k)=',h(i,k)
     3343           write(*,*) 'nk(i)=',nk(i)
     3344           write(*,*) 'h(i,nk(i))=',h(i,nk(i))
     3345           write(*,*) 'lv(i,k)=',lv(i,k)
     3346           write(*,*) 't(i,k)=',t(i,k)
     3347           write(*,*) 'clw(i,k)=',clw(i,k)
     3348           write(*,*) 'cpd,cpv=',cpd,cpv
     3349           stop
     3350        endif
     3351       enddo !do k=1,nl
     3352      enddo !do i=1,ncum 
     3353      endif !if (coef_epmax_cape.gt.1e-12) then
     3354
     3355      return
     3356      end subroutine cv30_epmax_fn_cape
     3357
     3358
Note: See TracChangeset for help on using the changeset viewer.