Changeset 2670 in lmdz_wrf


Ignore:
Timestamp:
Jul 12, 2019, 6:14:01 PM (5 years ago)
Author:
lfita
Message:

Adding pressure !

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/tools/module_ForDiagnosticsVars.f90

    r2669 r2670  
    20472047  END SUBROUTINE var_front_R04
    20482048
    2049   SUBROUTINE var_Frontogenesis(dx, dy, dz, dt, theta, ua, va, wa, ddx, ddy, ddz, xdiab, ydiab, zdiab, \
    2050     xdef, ydef, zdef, xtilt, ytilt, zdiv, f)
     2049  SUBROUTINE var_Frontogenesis(dx, dy, dz, dt, theta, ua, va, wa, press, ddx, ddy, ddz, ddt, xdiab,  \
     2050    ydiab, zdiab, xdef, ydef, zdef, xtilt, ytilt, zdiv, f)
    20512051  ! Subroutine to compute the terms of the equation of frontogenesis
    20522052  ! AFTER: https://en.wikipedia.org/wiki/Frontogenesis, http://glossary.ametsoc.org/wiki/Frontogenetical_function
     
    20552055
    20562056    INTEGER, INTENT(in)                                  :: dx, dy, dz, dt
     2057    REAL(r_k), INTENT(in)                                :: ddt
    20572058    REAL(r_k), DIMENSION(dx,dy), INTENT(in)              :: ddx,ddy
    20582059    REAL(r_k), DIMENSION(dz), INTENT(in)                 :: ddz
    2059     REAL(r_k), DIMENSION(dx,dy,dz,dt), INTENT(in)        :: theta, ua, va, wa
     2060    REAL(r_k), DIMENSION(dx,dy,dz,dt), INTENT(in)        :: theta, ua, va, wa, press
    20602061    REAL(r_k), DIMENSION(dx,dy,dz,dt), INTENT(out)       :: xdiab, ydiab, zdiab
    20612062    REAL(r_k), DIMENSION(dx,dy,dz,dt), INTENT(out)       :: xdef, ydef, zdef
     
    20652066    ! Local
    20662067    INTEGER                                              :: i,j,k,l,it
    2067     REAL(r_k), DIMENSION(dx,dy,dz)                       :: modthetagrad
    2068     REAL(r_k), DIMENSION(dx,dy,dz,3)                     :: thetagrad, uagrad, vagrad, wagrad,        &
    2069       thetadef, thetatilt
     2068    REAL(r_k), DIMENSION(dx,dy,dz)                       :: modthetagrad, thetadt
     2069    REAL(r_k), DIMENSION(dx,dy,dz,3)                     :: thetagrad, thetadtgrad, uagrad, vagrad,   &
     2070      wagrad, pressgrad, thetadef, thetatilt
    20702071
    20712072!!!!!!! Variables
     
    20732074! theta: potential temperature [K]
    20742075! ua, va, wa: x/y/z wind direction [ms-1]
     2076! press: pressure [Pa]
     2077! ddt: temporal delta [s]
    20752078! ddx, ddy, ddz: x/y/z distances among grid points [m]
    20762079! xdiab, ydiab, zdiab: x/y/z diabatic term [Ks-1m-1]
     
    20822085    fname = 'var_Frontogenesis'
    20832086
     2087    CALL deltat3D(dx, dy, dz, dt, theta, ddt, 'forward', thetadt)
     2088
    20842089    ! Computing components separately by time-step
    20852090    DO it=1, dt
     2091      CALL gradient3D_1o(dx, dy, dz, thetadt(:,:,:,it), ddx, ddy, ddz, thetadtgrad)
    20862092      CALL gradient3D_1o(dx, dy, dz, theta(:,:,:,it), ddx, ddy, ddz, thetagrad)
    20872093      CALL gradient3D_1o(dx, dy, dz, ua(:,:,:,it), ddx, ddy, ddz, uagrad)
    20882094      CALL gradient3D_1o(dx, dy, dz, va(:,:,:,it), ddx, ddy, ddz, vagrad)
    20892095      CALL gradient3D_1o(dx, dy, dz, wa(:,:,:,it), ddx, ddy, ddz, wagrad)
     2096      CALL gradient3D_1o(dx, dy, dz, press(:,:,:,it), ddx, ddy, ddz, pressgrad)
    20902097      CALL deformation3D(dx, dy, dz, thetagrad, uagrad, vagrad, thetadef)
    20912098      CALL tilting3D(dx, dy, dz, thetagrad, wagrad, thetatilt)
     2099      xdiab(:,:,:,it) = thetadtgrad(:,:,:,1)
     2100
    20922101      xdef(:,:,:,it) = thetadef(:,:,:,1)
    20932102      ydef(:,:,:,it) = thetadef(:,:,:,2)
Note: See TracChangeset for help on using the changeset viewer.