module cpdet_mod implicit none ! ADAPTATION OF GCM TO CP(T) !====================================================================== ! S. Lebonnois, 10/2010 ! ! Cp must be computed using cpdet(t) to be valid ! ! The Exner function is still pk = RCPD*(play/pref)**RKAPPA ! (RCPD=cpp, RKAPPA=kappa) ! ! One goes from T to teta (potential temperature) using t2tpot(t,teta,pk) ! One goes from teta to T using tpot2t(teta,t,pk) ! !====================================================================== contains SUBROUTINE ini_cpdet USE control_mod, ONLY: planet_type IMPLICIT none !====================================================================== ! Initialization of nu_venus and t0_venus !====================================================================== ! for cpp, nu_venus and t0_venus: #include "comconst.h" if (planet_type.eq."venus") then nu_venus=0.35 t0_venus=460. else nu_venus=0. t0_venus=0. endif return end subroutine ini_cpdet !====================================================================== !====================================================================== FUNCTION cpdet(t) USE control_mod, ONLY: planet_type IMPLICIT none ! for cpp, nu_venus and t0_venus: #include "comconst.h" real,intent(in) :: t real cpdet if (planet_type.eq."venus") then cpdet = cpp*(t/t0_venus)**nu_venus else cpdet = cpp endif return end function cpdet !====================================================================== !====================================================================== SUBROUTINE t2tpot(npoints, yt, yteta, ypk) !====================================================================== ! Arguments: ! ! yt --------input-R- Temperature ! yteta-------output-R- Temperature potentielle ! ypk --------input-R- Fonction d'Exner: RCPD*(pplay/pref)**RKAPPA ! !====================================================================== USE control_mod, ONLY: planet_type IMPLICIT NONE ! for cpp, nu_venus and t0_venus: #include "comconst.h" integer,intent(in) :: npoints REAL,intent(in) :: yt(npoints), ypk(npoints) REAL,intent(out) :: yteta(npoints) if (planet_type.eq."venus") then yteta = yt**nu_venus & & - nu_venus * t0_venus**nu_venus * log(ypk/cpp) yteta = yteta**(1./nu_venus) else yteta = yt * cpp/ypk endif return end subroutine t2tpot !====================================================================== !====================================================================== SUBROUTINE tpot2t(npoints,yteta, yt, ypk) !====================================================================== ! Arguments: ! ! yteta--------input-R- Temperature potentielle ! yt -------output-R- Temperature ! ypk --------input-R- Fonction d'Exner: RCPD*(pplay/pref)**RKAPPA ! !====================================================================== USE control_mod, ONLY: planet_type IMPLICIT NONE ! for cpp, nu_venus and t0_venus: #include "comconst.h" integer,intent(in) :: npoints REAL,intent(in) :: yteta(npoints), ypk(npoints) REAL,intent(out) :: yt(npoints) if (planet_type.eq."venus") then yt = yteta**nu_venus & & + nu_venus * t0_venus**nu_venus * log(ypk/cpp) yt = yt**(1./nu_venus) else yt = yteta * ypk/cpp endif return end subroutine tpot2t !====================================================================== !====================================================================== ! ! ATTENTION ! ! Si un jour on a besoin, il faudra coder les routines ! dt2dtpot / dtpto2dt ! !====================================================================== !====================================================================== end module cpdet_mod