Ignore:
Timestamp:
Aug 22, 2013, 4:02:07 PM (11 years ago)
Author:
emillour
Message:

Common dynamics: (and collateral adaptations in Venus physics)
Improved cpdet routines in and additional sponge mode:

  • Additionnal sponge mode (trigered with "callsponge" flag), in line with the one used in the Generic and Martian GCM. This sponge is called whenever there is a dissipation step.
  • Improvement of the cpdet routines : created routines tpot2t_glo_p and t2tpot_glo_p which handle fields on the whole dynamics (scaler) grid, which are more efficient than calling tpot2t_p or t2tpot_p with slabs of data (generated use of intermediate copies of these chunks of data at every call)
  • Turned cpdet.F into a module cpdet_mod.F90 (and correspondingly adapted all routines in the Venus physics).

EM

Location:
trunk/LMDZ.COMMON/libf/dyn3d
Files:
1 added
6 edited
1 moved

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.COMMON/libf/dyn3d/bilan_dyn.F

    r97 r1017  
    2020
    2121      USE control_mod, ONLY: planet_type
     22      USE cpdet_mod, only: tpot2t
    2223
    2324      IMPLICIT NONE
  • trunk/LMDZ.COMMON/libf/dyn3d/calfis.F

    r849 r1017  
    3131      USE infotrac
    3232      USE control_mod
    33  
     33      USE cpdet_mod, only: t2tpot,tpot2t
    3434
    3535      IMPLICIT NONE
  • trunk/LMDZ.COMMON/libf/dyn3d/conf_gcm.F

    r1012 r1017  
    1515      USE infotrac, ONLY : type_trac
    1616      use assert_m, only: assert
     17      use sponge_mod, only: callsponge,mode_sponge,nsponge,tetasponge
    1718
    1819      IMPLICIT NONE
     
    361362       tau_top_bound=1.e-5
    362363       CALL getin('tau_top_bound',tau_top_bound)
     364
     365! the other possible sponge layer (sponge_mod)
     366       callsponge=.false. ! default value; don't use the sponge
     367       call getin("callsponge",callsponge)
     368       ! check that user is not trying to use both sponge models
     369       if ((iflag_top_bound.ge.1).and.callsponge) then
     370         write(lunout,*)'Bad choice of options:'
     371         write(lunout,*)' iflag_top_bound=',iflag_top_bound
     372         write(lunout,*)' and callsponge=.true.'
     373         write(lunout,*)'But both sponge models should not be',
     374     &                  ' used simultaneously!'
     375         stop
     376       endif
     377       
     378! nsponge: number of atmospheric layers over which the sponge extends
     379       nsponge=3 ! default value
     380       call getin("nsponge",nsponge)
     381
     382! mode_sponge: (quenching is towards ... over the upper nsponge layers)
     383!      0: (h=hmean,u=v=0)
     384!      1: (h=hmean,u=umean,v=0)
     385!      2: (h=hmean,u=umean,v=vmean)"
     386       mode_sponge=2 ! default value
     387       call getin("mode_sponge",mode_sponge)
     388
     389! tetasponge: characteristic time scale (seconds) at topmost layer
     390!            (time scale then doubles with decreasing layer index)."
     391       tetasponge=50000.0
     392       call getin("tetasponge",tetasponge)
    363393
    364394! FOR TITAN: tidal forces
  • trunk/LMDZ.COMMON/libf/dyn3d/cpdet_mod.F90

    r1016 r1017  
    1 ! ADAPTATION GCM POUR CP(T)
    2 c======================================================================
    3 c S. Lebonnois, 10/2010
    4 c
    5 c Cp doit être calculé par cpdet(t) pour être valable partout
    6 c
    7 c La fonction d'Exner reste pk = RCPD*(play/pref)**RKAPPA
    8 c (RCPD=cpp, RKAPPA=kappa)
    9 c
    10 c On passe de T a teta (temperature potentielle) par t2tpot(t,teta,pk)
    11 c On passe de teta a T par tpot2t(teta,t,pk)
    12 c
    13 c======================================================================
     1      module cpdet_mod
     2
     3      implicit none
     4
     5! ADAPTATION OF GCM TO CP(T)
     6!======================================================================
     7! S. Lebonnois, 10/2010
     8!
     9! Cp must be computed using cpdet(t) to be valid
     10!
     11! The Exner function is still pk = RCPD*(play/pref)**RKAPPA
     12! (RCPD=cpp, RKAPPA=kappa)
     13!
     14! One goes from T to teta (potential temperature) using t2tpot(t,teta,pk)
     15! One goes from teta to T using tpot2t(teta,t,pk)
     16!
     17!======================================================================
     18
     19      contains
    1420
    1521      SUBROUTINE ini_cpdet
     
    1723      USE control_mod, ONLY: planet_type
    1824      IMPLICIT none
    19 c======================================================================
    20 c Initialisation de nu_venus et t0_venus
    21 c======================================================================
     25!======================================================================
     26! Initialization of nu_venus and t0_venus
     27!======================================================================
    2228
    2329! for cpp, nu_venus and t0_venus:
     
    3339
    3440      return
    35       end
     41      end subroutine ini_cpdet
    3642
    37 c======================================================================
    38 c======================================================================
     43!======================================================================
     44!======================================================================
    3945
    4046      FUNCTION cpdet(t)
     
    4652#include "comconst.h"
    4753
    48       real cpdet,t
     54      real,intent(in) :: t
     55      real cpdet
    4956
    5057      if (planet_type.eq."venus") then
     
    5562
    5663      return
    57       end
     64      end function cpdet
    5865     
    59 c======================================================================
    60 c======================================================================
     66!======================================================================
     67!======================================================================
    6168
    6269      SUBROUTINE t2tpot(npoints, yt, yteta, ypk)
    63 c======================================================================
    64 c Arguments:
    65 c
    66 c yt   --------input-R- Temperature
    67 c yteta-------output-R- Temperature potentielle
    68 c ypk  --------input-R- Fonction d'Exner: RCPD*(pplay/pref)**RKAPPA
    69 c
    70 c======================================================================
     70!======================================================================
     71! Arguments:
     72!
     73! yt   --------input-R- Temperature
     74! yteta-------output-R- Temperature potentielle
     75! ypk  --------input-R- Fonction d'Exner: RCPD*(pplay/pref)**RKAPPA
     76!
     77!======================================================================
    7178
    7279      USE control_mod, ONLY: planet_type
     
    7683#include "comconst.h"
    7784
    78       integer npoints
    79       REAL    yt(npoints), yteta(npoints), ypk(npoints)
     85      integer,intent(in) :: npoints
     86      REAL,intent(in) :: yt(npoints), ypk(npoints)
     87      REAL,intent(out) :: yteta(npoints)
    8088     
    8189      if (planet_type.eq."venus") then
     
    8896
    8997      return
    90       end
     98      end subroutine t2tpot
    9199
    92 c======================================================================
    93 c======================================================================
     100!======================================================================
     101!======================================================================
    94102
    95103      SUBROUTINE tpot2t(npoints,yteta, yt, ypk)
    96 c======================================================================
    97 c Arguments:
    98 c
    99 c yteta--------input-R- Temperature potentielle
    100 c yt   -------output-R- Temperature
    101 c ypk  --------input-R- Fonction d'Exner: RCPD*(pplay/pref)**RKAPPA
    102 c
    103 c======================================================================
     104!======================================================================
     105! Arguments:
     106!
     107! yteta--------input-R- Temperature potentielle
     108! yt   -------output-R- Temperature
     109! ypk  --------input-R- Fonction d'Exner: RCPD*(pplay/pref)**RKAPPA
     110!
     111!======================================================================
    104112
    105113      USE control_mod, ONLY: planet_type
     
    109117#include "comconst.h"
    110118
    111       integer npoints
    112       REAL    yt(npoints), yteta(npoints), ypk(npoints)
     119      integer,intent(in) :: npoints
     120      REAL,intent(in) :: yteta(npoints), ypk(npoints)
     121      REAL,intent(out) :: yt(npoints)
    113122     
    114123      if (planet_type.eq."venus") then
     
    121130 
    122131      return
    123       end
     132      end subroutine tpot2t
    124133
    125 c======================================================================
    126 c======================================================================
    127 c
    128 c ATTENTION
    129 c
    130 c Si un jour on a besoin, il faudra coder les routines
    131 c    dt2dtpot / dtpto2dt
    132 c
    133 c======================================================================
    134 c======================================================================
     134!======================================================================
     135!======================================================================
     136!
     137! ATTENTION
     138!
     139! Si un jour on a besoin, il faudra coder les routines
     140!    dt2dtpot / dtpto2dt
     141!
     142!======================================================================
     143!======================================================================
     144      end module cpdet_mod
  • trunk/LMDZ.COMMON/libf/dyn3d/gcm.F

    r979 r1017  
    1616      USE infotrac
    1717      USE control_mod
     18      use cpdet_mod, only: ini_cpdet
    1819
    1920!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  • trunk/LMDZ.COMMON/libf/dyn3d/leapfrog.F

    r1012 r1017  
    1616      USE write_field
    1717      USE control_mod
     18      use cpdet_mod, only: cpdet,tpot2t,t2tpot
     19      use sponge_mod, only: callsponge,mode_sponge,sponge
    1820      IMPLICIT NONE
    1921
     
    197199      ! for CP(T)
    198200      real :: dtec
    199       real,external :: cpdet
    200201      real :: ztetaec(ip1jmp1,llm)
    201202
    202203c dummy: sinon cette routine n'est jamais compilee...
    203204      if(1.eq.0) then
     205#ifdef CPP_PHYS
    204206        CALL init_phys_lmdz(iim,jjp1,llm,1,(/(jjm-1)*iim+2/))
     207#endif
    205208      endif
    206209
     
    562565      IF(apdiss) THEN
    563566
     567        ! sponge layer
     568        if (callsponge) then
     569          CALL sponge(ucov,vcov,teta,ps,dtdiss,mode_sponge)
     570        endif
    564571
    565572c   calcul de l'energie cinetique avant dissipation
  • trunk/LMDZ.COMMON/libf/dyn3d/vlspltqs.F

    r109 r1017  
    2222c     pk exner au milieu des couches necessaire pour calculer Qsat
    2323c   --------------------------------------------------------------------
     24      use cpdet_mod, only: tpot2t
    2425      IMPLICIT NONE
    2526c
Note: See TracChangeset for help on using the changeset viewer.