New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
Changeset 5770 for branches/2015/dev_r5721_CNRS9_NOC3_LDF/NEMOGCM/NEMO/OPA_SRC/TRA/traadv.F90 – NEMO

Ignore:
Timestamp:
2015-10-01T14:48:08+02:00 (9 years ago)
Author:
gm
Message:

#1593: LDF-ADV, step II.2: phasing the improvements/simplifications of advective tracer trend (see wiki page)

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/2015/dev_r5721_CNRS9_NOC3_LDF/NEMOGCM/NEMO/OPA_SRC/TRA/traadv.F90

    r5758 r5770  
    77   !!            3.3  !  2010-09  (C. Ethe, G. Madec)  merge TRC-TRA + switch from velocity to transport 
    88   !!            3.6  !  2011-06  (G. Madec)  Addition of Mixed Layer Eddy parameterisation 
    9    !!            4.0  !  2014-05  (G. Madec)  Add 2nd/4th order cases for CEN and FCT schemes  
     9   !!            3.7  !  2014-05  (G. Madec)  Add 2nd/4th order cases for CEN and FCT schemes  
    1010   !!             -   !  2014-12  (G. Madec) suppression of cross land advection option 
    1111   !!---------------------------------------------------------------------- 
    1212 
    1313   !!---------------------------------------------------------------------- 
    14    !!   tra_adv      : compute ocean tracer advection trend 
    15    !!   tra_adv_ctl  : control the different options of advection scheme 
    16    !!---------------------------------------------------------------------- 
    17    USE oce             ! ocean dynamics and active tracers 
    18    USE dom_oce         ! ocean space and time domain 
    19    USE domvvl          ! variable vertical scale factors 
    20    USE traadv_cen2     ! 2nd order centered scheme (tra_adv_cen2   routine) 
    21    USE traadv_tvd      ! TVD      scheme           (tra_adv_tvd    routine) 
    22    USE traadv_muscl    ! MUSCL    scheme           (tra_adv_muscl  routine) 
    23    USE traadv_muscl2   ! MUSCL2   scheme           (tra_adv_muscl2 routine) 
    24    USE traadv_ubs      ! UBS      scheme           (tra_adv_ubs    routine) 
    25    USE traadv_qck      ! QUICKEST scheme           (tra_adv_qck    routine) 
    26    USE traadv_mle      ! ML eddy induced velocity  (tra_adv_mle    routine) 
    27    USE cla             ! cross land advection      (cla_traadv     routine) 
    28    USE ldftra          ! lateral diffusion coefficient on tracers 
    29    USE ldfslp          ! Lateral diffusion: slopes of neutral surfaces 
     14   !!   tra_adv       : compute ocean tracer advection trend 
     15   !!   tra_adv_ctl   : control the different options of advection scheme 
     16   !!---------------------------------------------------------------------- 
     17   USE oce            ! ocean dynamics and active tracers 
     18   USE dom_oce        ! ocean space and time domain 
     19   USE domvvl         ! variable vertical scale factors 
     20   USE traadv_cen     ! centered scheme           (tra_adv_cen  routine) 
     21   USE traadv_fct     ! FCT      scheme           (tra_adv_fct  routine) 
     22   USE traadv_mus     ! MUSCL    scheme           (tra_adv_mus  routine) 
     23   USE traadv_ubs     ! UBS      scheme           (tra_adv_ubs  routine) 
     24   USE traadv_qck     ! QUICKEST scheme           (tra_adv_qck  routine) 
     25   USE traadv_mle     ! ML eddy induced velocity  (tra_adv_mle  routine) 
     26   USE ldftra         ! lateral diffusion: eddy diffusivity & EIV coeff. 
     27   USE ldfslp         ! Lateral diffusion: slopes of neutral surfaces 
    3028   ! 
    31    USE in_out_manager  ! I/O manager 
    32    USE iom             ! I/O module 
    33    USE prtctl          ! Print control 
    34    USE lib_mpp         ! MPP library 
    35    USE wrk_nemo        ! Memory Allocation 
    36    USE timing          ! Timing 
    37    USE sbc_oce 
     29   USE in_out_manager ! I/O manager 
     30   USE iom            ! I/O module 
     31   USE prtctl         ! Print control 
     32   USE lib_mpp        ! MPP library 
     33   USE wrk_nemo       ! Memory Allocation 
     34   USE timing         ! Timing 
     35 
    3836   USE diaptr          ! Poleward heat transport  
    39  
    4037 
    4138   IMPLICIT NONE 
     
    4542   PUBLIC   tra_adv_init   ! routine called by opa module 
    4643 
    47    !                              !!* Namelist namtra_adv * 
    48    LOGICAL ::   ln_traadv_cen2     ! 2nd order centered scheme flag 
    49    LOGICAL ::   ln_traadv_tvd      ! TVD scheme flag 
    50    LOGICAL ::   ln_traadv_tvd_zts  ! TVD scheme flag with vertical sub time-stepping 
    51    LOGICAL ::   ln_traadv_muscl    ! MUSCL scheme flag 
    52    LOGICAL ::   ln_traadv_muscl2   ! MUSCL2 scheme flag 
    53    LOGICAL ::   ln_traadv_ubs      ! UBS scheme flag 
    54    LOGICAL ::   ln_traadv_qck      ! QUICKEST scheme flag 
    55    LOGICAL ::   ln_traadv_msc_ups  ! use upstream scheme within muscl 
    56  
    57  
    58    INTEGER ::   nadv   ! choice of the type of advection scheme 
    59  
     44   !                            !!* Namelist namtra_adv * 
     45   LOGICAL ::   ln_traadv_cen    ! centered scheme flag 
     46   INTEGER ::      nn_cen_h, nn_cen_v   ! =2/4 : horizontal and vertical choices of the order of CEN scheme 
     47   LOGICAL ::   ln_traadv_fct    ! FCT scheme flag 
     48   INTEGER ::      nn_fct_h, nn_fct_v   ! =2/4 : horizontal and vertical choices of the order of FCT scheme 
     49   INTEGER ::      nn_fct_zts           ! >=1 : 2nd order FCT with vertical sub-timestepping 
     50   LOGICAL ::   ln_traadv_mus    ! MUSCL scheme flag 
     51   LOGICAL ::      ln_mus_ups           ! use upstream scheme in vivcinity of river mouths 
     52   LOGICAL ::   ln_traadv_ubs    ! UBS scheme flag 
     53   INTEGER ::      nn_ubs_v             ! =2/4 : vertical choice of the order of UBS scheme 
     54   LOGICAL ::   ln_traadv_qck    ! QUICKEST scheme flag 
     55 
     56   INTEGER ::              nadv             ! choice of the type of advection scheme 
     57   ! 
     58   !                                        ! associated indices: 
     59   INTEGER, PARAMETER ::   np_NO_adv  = 0   ! no T-S advection 
     60   INTEGER, PARAMETER ::   np_CEN     = 1   ! 2nd/4th order centered scheme 
     61   INTEGER, PARAMETER ::   np_FCT     = 2   ! 2nd/4th order Flux Corrected Transport scheme 
     62   INTEGER, PARAMETER ::   np_FCT_zts = 3   ! 2nd order FCT scheme with vertical sub-timestepping 
     63   INTEGER, PARAMETER ::   np_MUS     = 4   ! MUSCL scheme 
     64   INTEGER, PARAMETER ::   np_UBS     = 5   ! 3rd order Upstream Biased Scheme 
     65   INTEGER, PARAMETER ::   np_QCK     = 6   ! QUICK scheme 
     66    
    6067   !! * Substitutions 
    6168#  include "domzgr_substitute.h90" 
    6269#  include "vectopt_loop_substitute.h90" 
    6370   !!---------------------------------------------------------------------- 
    64    !! NEMO/OPA 3.3 , NEMO Consortium (2010) 
     71   !! NEMO/OPA 3.7 , NEMO Consortium (2014) 
    6572   !! $Id$ 
    6673   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt) 
     
    8491      IF( nn_timing == 1 )  CALL timing_start('tra_adv') 
    8592      ! 
    86       CALL wrk_alloc( jpi, jpj, jpk, zun, zvn, zwn ) 
     93      CALL wrk_alloc( jpi,jpj,jpk,  zun, zvn, zwn ) 
    8794      ! 
    8895      !                                          ! set time step 
     
    93100      ENDIF 
    94101      ! 
    95       IF( nn_cla == 1 .AND. cp_cfg == 'orca' .AND. jp_cfg == 2 )   CALL cla_traadv( kt )       !==  Cross Land Advection  ==! (hor. advection) 
    96       ! 
    97       !                                               !==  effective transport  ==! 
     102      !                                         !==  effective transport  ==! 
    98103      DO jk = 1, jpkm1 
    99104         zun(:,:,jk) = e2u  (:,:) * fse3u(:,:,jk) * un(:,:,jk)                  ! eulerian transport only 
     
    102107      END DO 
    103108      ! 
    104       IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN 
     109      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN                                ! add z-tilde and/or vvl corrections 
    105110         zun(:,:,:) = zun(:,:,:) + un_td(:,:,:) 
    106111         zvn(:,:,:) = zvn(:,:,:) + vn_td(:,:,:) 
    107112      ENDIF 
    108113      ! 
    109       zun(:,:,jpk) = 0._wp                                                     ! no transport trough the bottom 
    110       zvn(:,:,jpk) = 0._wp                                                     ! no transport trough the bottom 
    111       zwn(:,:,jpk) = 0._wp                                                     ! no transport trough the bottom 
     114      zun(:,:,jpk) = 0._wp                                                      ! no transport trough the bottom 
     115      zvn(:,:,jpk) = 0._wp 
     116      zwn(:,:,jpk) = 0._wp 
    112117      ! 
    113118      IF( ln_ldfeiv .AND. .NOT. ln_traldf_triad )   & 
     
    120125      CALL iom_put( "wocetr_eff", zwn ) 
    121126      ! 
     127!!gm ??? 
    122128      IF( ln_diaptr )   CALL dia_ptr( zvn )                                    ! diagnose the effective MSF  
    123       ! 
    124     
    125       SELECT CASE ( nadv )                            !==  compute advection trend and add it to general trend  ==! 
    126       CASE ( 1 )   ;    CALL tra_adv_cen2   ( kt, nit000, 'TRA',         zun, zvn, zwn, tsb, tsn, tsa, jpts )   !  2nd order centered 
    127       CASE ( 2 )   ;    CALL tra_adv_tvd    ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts )   !  TVD  
    128       CASE ( 3 )   ;    CALL tra_adv_muscl  ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb,      tsa, jpts, ln_traadv_msc_ups )   !  MUSCL  
    129       CASE ( 4 )   ;    CALL tra_adv_muscl2 ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts )   !  MUSCL2  
    130       CASE ( 5 )   ;    CALL tra_adv_ubs    ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts )   !  UBS  
    131       CASE ( 6 )   ;    CALL tra_adv_qck    ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts )   !  QUICKEST  
    132       CASE ( 7 )   ;    CALL tra_adv_tvd_zts( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts )   !  TVD ZTS 
    133       ! 
    134       CASE (-1 )                                      !==  esopa: test all possibility with control print  ==! 
    135          CALL tra_adv_cen2  ( kt, nit000, 'TRA',         zun, zvn, zwn, tsb, tsn, tsa, jpts )           
    136          CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' adv0 - Ta: ', mask1=tmask,               & 
    137             &          tab3d_2=tsa(:,:,:,jp_sal), clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' ) 
    138          CALL tra_adv_tvd   ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts )           
    139          CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' adv1 - Ta: ', mask1=tmask,               & 
    140             &          tab3d_2=tsa(:,:,:,jp_sal), clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' ) 
    141          CALL tra_adv_muscl ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb,      tsa, jpts, ln_traadv_msc_ups )           
    142          CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' adv3 - Ta: ', mask1=tmask,               & 
    143             &          tab3d_2=tsa(:,:,:,jp_sal), clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' ) 
    144          CALL tra_adv_muscl2( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts )           
    145          CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' adv4 - Ta: ', mask1=tmask,               & 
    146             &          tab3d_2=tsa(:,:,:,jp_sal), clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' ) 
    147          CALL tra_adv_ubs   ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts )           
    148          CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' adv5 - Ta: ', mask1=tmask,               & 
    149             &          tab3d_2=tsa(:,:,:,jp_sal), clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' ) 
    150          CALL tra_adv_qck   ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts )           
    151          CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' adv6 - Ta: ', mask1=tmask,               & 
    152             &          tab3d_2=tsa(:,:,:,jp_sal), clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' ) 
     129!!gm ??? 
     130      ! 
     131      SELECT CASE ( nadv )                      !==  compute advection trend and add it to general trend  ==! 
     132      ! 
     133      CASE ( np_CEN )                                    ! Centered scheme : 2nd / 4th order 
     134         CALL tra_adv_cen    ( kt, nit000, 'TRA',         zun, zvn, zwn     , tsn, tsa, jpts, nn_cen_h, nn_cen_v ) 
     135      CASE ( np_FCT )                                    ! FCT scheme      : 2nd / 4th order 
     136         CALL tra_adv_fct    ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts, nn_fct_h, nn_fct_v ) 
     137      CASE ( np_FCT_zts )                                ! 2nd order FCT with vertical time-splitting 
     138         CALL tra_adv_fct_zts( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts        , nn_fct_zts ) 
     139      CASE ( np_MUS )                                    ! MUSCL 
     140         CALL tra_adv_mus    ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb,      tsa, jpts        , ln_mus_ups )  
     141      CASE ( np_UBS )                                    ! UBS 
     142         CALL tra_adv_ubs    ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts        , nn_ubs_v   ) 
     143      CASE ( np_QCK )                                    ! QUICKEST 
     144         CALL tra_adv_qck    ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts                     ) 
     145      ! 
    153146      END SELECT 
    154147      ! 
     
    159152      IF( nn_timing == 1 )  CALL timing_stop( 'tra_adv' ) 
    160153      ! 
    161       CALL wrk_dealloc( jpi, jpj, jpk, zun, zvn, zwn ) 
     154      CALL wrk_dealloc( jpi,jpj,jpk,  zun, zvn, zwn ) 
    162155      !                                           
    163156   END SUBROUTINE tra_adv 
     
    171164      !!              tracer advection schemes and set nadv 
    172165      !!---------------------------------------------------------------------- 
    173       INTEGER ::   ioptio 
    174       INTEGER ::   ios                 ! Local integer output status for namelist read 
    175       !! 
    176       NAMELIST/namtra_adv/ ln_traadv_cen2 , ln_traadv_tvd,     & 
    177          &                 ln_traadv_muscl, ln_traadv_muscl2,  & 
    178          &                 ln_traadv_ubs  , ln_traadv_qck,     & 
    179          &                 ln_traadv_msc_ups, ln_traadv_tvd_zts 
    180       !!---------------------------------------------------------------------- 
    181  
    182       REWIND( numnam_ref )              ! Namelist namtra_adv in reference namelist : Tracer advection scheme 
     166      INTEGER ::   ioptio, ios   ! Local integers 
     167      ! 
     168      NAMELIST/namtra_adv/ ln_traadv_cen, nn_cen_h, nn_cen_v,               &   ! CEN 
     169         &                 ln_traadv_fct, nn_fct_h, nn_fct_v, nn_fct_zts,   &   ! FCT 
     170         &                 ln_traadv_mus,                     ln_mus_ups,   &   ! MUSCL 
     171         &                 ln_traadv_ubs,           nn_ubs_v,               &   ! UBS 
     172         &                 ln_traadv_qck                                        ! QCK 
     173      !!---------------------------------------------------------------------- 
     174      ! 
     175      !                                !==  Namelist  ==! 
     176      ! 
     177      REWIND( numnam_ref )                   ! Namelist namtra_adv in reference namelist : Tracer advection scheme 
    183178      READ  ( numnam_ref, namtra_adv, IOSTAT = ios, ERR = 901) 
    184179901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namtra_adv in reference namelist', lwp ) 
    185  
    186       REWIND( numnam_cfg )              ! Namelist namtra_adv in configuration namelist : Tracer advection scheme 
     180      ! 
     181      REWIND( numnam_cfg )                   ! Namelist namtra_adv in configuration namelist : Tracer advection scheme 
    187182      READ  ( numnam_cfg, namtra_adv, IOSTAT = ios, ERR = 902 ) 
    188183902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namtra_adv in configuration namelist', lwp ) 
    189184      IF(lwm) WRITE ( numond, namtra_adv ) 
    190185 
    191       IF(lwp) THEN                    ! Namelist print 
     186      IF(lwp) THEN                           ! Namelist print 
    192187         WRITE(numout,*) 
    193188         WRITE(numout,*) 'tra_adv_init : choice/control of the tracer advection scheme' 
    194189         WRITE(numout,*) '~~~~~~~~~~~' 
    195190         WRITE(numout,*) '   Namelist namtra_adv : chose a advection scheme for tracers' 
    196          WRITE(numout,*) '      2nd order advection scheme     ln_traadv_cen2    = ', ln_traadv_cen2 
    197          WRITE(numout,*) '      TVD advection scheme           ln_traadv_tvd     = ', ln_traadv_tvd 
    198          WRITE(numout,*) '      MUSCL  advection scheme        ln_traadv_muscl   = ', ln_traadv_muscl 
    199          WRITE(numout,*) '      MUSCL2 advection scheme        ln_traadv_muscl2  = ', ln_traadv_muscl2 
    200          WRITE(numout,*) '      UBS    advection scheme        ln_traadv_ubs     = ', ln_traadv_ubs 
    201          WRITE(numout,*) '      QUICKEST advection scheme      ln_traadv_qck     = ', ln_traadv_qck 
    202          WRITE(numout,*) '      upstream scheme within muscl   ln_traadv_msc_ups = ', ln_traadv_msc_ups 
    203          WRITE(numout,*) '      TVD advection scheme with zts  ln_traadv_tvd_zts = ', ln_traadv_tvd_zts 
    204       ENDIF 
    205  
    206       ioptio = 0                      ! Parameter control 
    207       IF( ln_traadv_cen2   )   ioptio = ioptio + 1 
    208       IF( ln_traadv_tvd    )   ioptio = ioptio + 1 
    209       IF( ln_traadv_muscl  )   ioptio = ioptio + 1 
    210       IF( ln_traadv_muscl2 )   ioptio = ioptio + 1 
    211       IF( ln_traadv_ubs    )   ioptio = ioptio + 1 
    212       IF( ln_traadv_qck    )   ioptio = ioptio + 1 
    213       IF( ln_traadv_tvd_zts)   ioptio = ioptio + 1 
    214       IF( lk_esopa         )   ioptio =          1 
    215  
    216       IF( ( ln_traadv_muscl .OR. ln_traadv_muscl2 .OR. ln_traadv_ubs .OR. ln_traadv_qck .OR. ln_traadv_tvd_zts )   & 
    217          .AND. ln_isfcav )   CALL ctl_stop( 'Only traadv_cen2 and traadv_tvd is compatible with ice shelf cavity') 
    218  
    219       IF( ioptio /= 1 )   CALL ctl_stop( 'Choose ONE advection scheme in namelist namtra_adv' ) 
    220  
    221       !                              ! Set nadv 
    222       IF( ln_traadv_cen2   )   nadv =  1 
    223       IF( ln_traadv_tvd    )   nadv =  2 
    224       IF( ln_traadv_muscl  )   nadv =  3 
    225       IF( ln_traadv_muscl2 )   nadv =  4 
    226       IF( ln_traadv_ubs    )   nadv =  5 
    227       IF( ln_traadv_qck    )   nadv =  6 
    228       IF( ln_traadv_tvd_zts)   nadv =  7 
    229       IF( lk_esopa         )   nadv = -1 
    230  
    231       IF(lwp) THEN                   ! Print the choice 
     191         WRITE(numout,*) '      centered scheme                           ln_traadv_cen = ', ln_traadv_cen 
     192         WRITE(numout,*) '            horizontal 2nd/4th order               nn_cen_h   = ', nn_fct_h 
     193         WRITE(numout,*) '            vertical   2nd/4th order               nn_cen_v   = ', nn_fct_v 
     194         WRITE(numout,*) '      Flux Corrected Transport scheme           ln_traadv_fct = ', ln_traadv_fct 
     195         WRITE(numout,*) '            horizontal 2nd/4th order               nn_fct_h   = ', nn_fct_h 
     196         WRITE(numout,*) '            vertical   2nd/4th order               nn_fct_v   = ', nn_fct_v 
     197         WRITE(numout,*) '            2nd order + vertical sub-timestepping  nn_fct_zts = ', nn_fct_zts 
     198         WRITE(numout,*) '      MUSCL scheme                              ln_traadv_mus = ', ln_traadv_mus 
     199         WRITE(numout,*) '            + upstream scheme near river mouths    ln_mus_ups = ', ln_mus_ups 
     200         WRITE(numout,*) '      UBS scheme                                ln_traadv_ubs = ', ln_traadv_ubs 
     201         WRITE(numout,*) '            vertical   2nd/4th order               nn_ubs_v   = ', nn_ubs_v 
     202         WRITE(numout,*) '      QUICKEST scheme                           ln_traadv_qck = ', ln_traadv_qck 
     203      ENDIF 
     204 
     205      ioptio = 0                       !==  Parameter control  ==! 
     206      IF( ln_traadv_cen )   ioptio = ioptio + 1 
     207      IF( ln_traadv_fct )   ioptio = ioptio + 1 
     208      IF( ln_traadv_mus )   ioptio = ioptio + 1 
     209      IF( ln_traadv_ubs )   ioptio = ioptio + 1 
     210      IF( ln_traadv_qck )   ioptio = ioptio + 1 
     211      ! 
     212      IF( ioptio == 0 ) THEN 
     213         nadv = np_NO_adv 
     214         CALL ctl_warn( 'tra_adv_init: You are running without tracer advection.' ) 
     215      ENDIF 
     216      IF( ioptio /= 1 )   CALL ctl_stop( 'tra_adv_init: Choose ONE advection scheme in namelist namtra_adv' ) 
     217      ! 
     218      IF( ln_traadv_cen .AND. ( nn_cen_h /= 2 .AND. nn_cen_h /= 4 )   &          ! Centered 
     219                        .AND. ( nn_cen_v /= 2 .AND. nn_cen_v /= 4 )   ) THEN 
     220        CALL ctl_stop( 'tra_adv_init: CEN scheme, choose 2nd or 4th order' ) 
     221      ENDIF 
     222      IF( ln_traadv_fct .AND. ( nn_fct_h /= 2 .AND. nn_fct_h /= 4 )   &          ! FCT 
     223                        .AND. ( nn_fct_v /= 2 .AND. nn_fct_v /= 4 )   ) THEN 
     224        CALL ctl_stop( 'tra_adv_init: FCT scheme, choose 2nd or 4th order' ) 
     225      ENDIF 
     226      IF( ln_traadv_fct .AND. nn_fct_zts > 0 ) THEN 
     227         IF( nn_fct_h == 4 ) THEN 
     228            nn_fct_h = 2 
     229            CALL ctl_stop( 'tra_adv_init: force 2nd order FCT scheme, 4th order does not exist with sub-timestepping' ) 
     230         ENDIF 
     231         IF( lk_vvl ) THEN 
     232            CALL ctl_stop( 'tra_adv_init: vertical sub-timestepping not allow in non-linear free surface' ) 
     233         ENDIF 
     234         IF( nn_fct_zts == 1 )   CALL ctl_warn( 'tra_adv_init: FCT with ONE sub-timestep = FCT without sub-timestep' ) 
     235      ENDIF 
     236      IF( ln_traadv_ubs .AND. ( nn_ubs_v /= 2 .AND. nn_ubs_v /= 4 )   ) THEN     ! UBS 
     237        CALL ctl_stop( 'tra_adv_init: UBS scheme, choose 2nd or 4th order' ) 
     238      ENDIF 
     239      IF( ln_traadv_ubs .AND. nn_ubs_v == 4 ) THEN 
     240         CALL ctl_warn( 'tra_adv_init: UBS scheme, only 2nd FCT scheme available on the vertical. It will be used' ) 
     241      ENDIF 
     242      IF( ln_isfcav ) THEN                                                       ! ice-shelf cavities 
     243         IF(  ln_traadv_cen .AND. nn_cen_v /= 4    .OR.   &                            ! NO 4th order with ISF 
     244            & ln_traadv_fct .AND. nn_fct_v /= 4   )   CALL ctl_stop( 'tra_adv_init: 4th order COMPACT scheme not allowed with ISF' ) 
     245      ENDIF 
     246      ! 
     247      !                                !==  used advection scheme  ==!   
     248      !                                      ! set nadv 
     249      IF( ln_traadv_cen                      )   nadv = np_CEN 
     250      IF( ln_traadv_fct                      )   nadv = np_FCT 
     251      IF( ln_traadv_fct .AND. nn_fct_zts > 0 )   nadv = np_FCT_zts 
     252      IF( ln_traadv_mus                      )   nadv = np_MUS 
     253      IF( ln_traadv_ubs                      )   nadv = np_UBS 
     254      IF( ln_traadv_qck                      )   nadv = np_QCK 
     255 
     256      IF(lwp) THEN                           ! Print the choice 
    232257         WRITE(numout,*) 
    233          IF( nadv ==  1 )   WRITE(numout,*) '         2nd order scheme is used' 
    234          IF( nadv ==  2 )   WRITE(numout,*) '         TVD       scheme is used' 
    235          IF( nadv ==  3 )   WRITE(numout,*) '         MUSCL     scheme is used' 
    236          IF( nadv ==  4 )   WRITE(numout,*) '         MUSCL2    scheme is used' 
    237          IF( nadv ==  5 )   WRITE(numout,*) '         UBS       scheme is used' 
    238          IF( nadv ==  6 )   WRITE(numout,*) '         QUICKEST  scheme is used' 
    239          IF( nadv ==  7 )   WRITE(numout,*) '         TVD ZTS   scheme is used' 
    240          IF( nadv == -1 )   WRITE(numout,*) '         esopa test: use all advection scheme' 
    241       ENDIF 
    242       ! 
    243       CALL tra_adv_mle_init          ! initialisation of the Mixed Layer Eddy parametrisation (MLE) 
     258         IF( nadv == np_NO_adv  )   WRITE(numout,*) '         NO T-S advection' 
     259         IF( nadv == np_CEN     )   WRITE(numout,*) '         CEN      scheme is used. Horizontal order: ', nn_cen_h,   & 
     260            &                                                                        ' Vertical   order: ', nn_cen_v 
     261         IF( nadv == np_FCT     )   WRITE(numout,*) '         FCT      scheme is used. Horizontal order: ', nn_fct_h,   & 
     262            &                                                                        ' Vertical   order: ', nn_fct_v 
     263         IF( nadv == np_FCT_zts )   WRITE(numout,*) '         use 2nd order FCT with ', nn_fct_zts,'vertical sub-timestepping' 
     264         IF( nadv == np_MUS     )   WRITE(numout,*) '         MUSCL    scheme is used' 
     265         IF( nadv == np_UBS     )   WRITE(numout,*) '         UBS      scheme is used' 
     266         IF( nadv == np_QCK     )   WRITE(numout,*) '         QUICKEST scheme is used' 
     267      ENDIF 
     268      ! 
     269      CALL tra_adv_mle_init            !== initialisation of the Mixed Layer Eddy parametrisation (MLE)  ==! 
    244270      ! 
    245271   END SUBROUTINE tra_adv_init 
Note: See TracChangeset for help on using the changeset viewer.