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 503 for trunk/NEMO/OPA_SRC/TRA/traadv.F90 – NEMO

Ignore:
Timestamp:
2006-09-27T10:52:29+02:00 (18 years ago)
Author:
opalod
Message:

nemo_v1_update_064 : CT : general trends update including the addition of mean windows analysis possibility in the mixed layer

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/NEMO/OPA_SRC/TRA/traadv.F90

    r474 r503  
    44   !! Ocean active tracers:  advection trend  
    55   !!============================================================================== 
    6    !! History : 
    7    !!   9.0  !  05-11  (G. Madec)  Original code 
     6   !! History :  9.0  !  05-11  (G. Madec)  Original code 
     7   !!---------------------------------------------------------------------- 
     8 
    89   !!---------------------------------------------------------------------- 
    910   !!   tra_adv      : compute ocean tracer advection trend 
    1011   !!   tra_adv_ctl  : control the different options of advection scheme 
    1112   !!---------------------------------------------------------------------- 
    12    !! * Modules used 
    1313   USE oce             ! ocean dynamics and active tracers 
    1414   USE dom_oce         ! ocean space and time domain 
    15    USE traadv_cen2     ! 2nd order centered scheme   (tra_adv_cen2 routine) 
    16    USE traadv_cen2_jki ! 2nd order centered scheme   (tra_adv_cen2 routine) 
    17    USE traadv_tvd      ! TVD scheme                (tra_adv_tvd    routine) 
    18    USE traadv_muscl    ! MUSCL scheme              (tra_adv_muscl  routine) 
     15   USE traadv_cen2     ! 2nd order centered scheme (tra_adv_cen2  routine) 
     16   USE traadv_cen2_jki ! 2nd order centered scheme (tra_adv_cen2  routine) 
     17   USE traadv_tvd      ! TVD    scheme             (tra_adv_tvd    routine) 
     18   USE traadv_muscl    ! MUSCL  scheme             (tra_adv_muscl  routine) 
    1919   USE traadv_muscl2   ! MUSCL2 scheme             (tra_adv_muscl2 routine) 
     20   USE traadv_ubs      ! UBS    scheme             (tra_adv_ubs    routine) 
    2021   USE traadv_eiv      ! eddy induced velocity     (tra_adv_eiv    routine) 
    21    USE trabbl          ! ??? 
    22    USE ldftra_oce      ! ??? 
     22   USE trabbl          ! tracers: bottom boundary layer 
     23   USE ldftra_oce      ! lateral diffusion coefficient on tracers 
    2324   USE in_out_manager  ! I/O manager 
    2425   USE prtctl          ! Print control 
     
    2728   PRIVATE 
    2829 
    29    !! * Accessibility 
    30    PUBLIC tra_adv         ! routine called by step module 
     30   PUBLIC   tra_adv    ! routine called by step module 
    3131  
    32    !! * Share module variables 
    33    LOGICAL, PUBLIC ::   & 
    34       ln_traadv_cen2   = .TRUE.  ,   &  ! 2nd order centered scheme flag 
    35       ln_traadv_tvd    = .FALSE. ,   &  ! TVD scheme flag 
    36       ln_traadv_muscl  = .FALSE. ,   &  ! MUSCL scheme flag 
    37       ln_traadv_muscl2 = .FALSE.        ! MUSCL2 scheme flag 
     32   !!* Namelist nam_traadv 
     33   LOGICAL, PUBLIC ::   ln_traadv_cen2   = .TRUE.       ! 2nd order centered scheme flag 
     34   LOGICAL, PUBLIC ::   ln_traadv_tvd    = .FALSE.      ! TVD scheme flag 
     35   LOGICAL, PUBLIC ::   ln_traadv_muscl  = .FALSE.      ! MUSCL scheme flag 
     36   LOGICAL, PUBLIC ::   ln_traadv_muscl2 = .FALSE.      ! MUSCL2 scheme flag 
     37   LOGICAL, PUBLIC ::   ln_traadv_ubs    = .FALSE.      ! UBS scheme flag 
     38   NAMELIST/nam_traadv/ ln_traadv_cen2 , ln_traadv_tvd,   & 
     39      &                 ln_traadv_muscl, ln_traadv_muscl2, ln_traadv_ubs 
    3840 
    39    !! * Module variables 
    40    INTEGER ::   & 
    41       nadv                              ! choice of the type of advection scheme 
     41   INTEGER ::   nadv   ! choice of the type of advection scheme 
    4242 
    4343   !! * Substitutions 
    44  include "domzgr_substitute.h90" 
    45  include "vectopt_loop_substitute.h90" 
     44include "domzgr_substitute.h90" 
     45include "vectopt_loop_substitute.h90" 
    4646   !!---------------------------------------------------------------------- 
    47    !!   OPA 9.0 , LOCEAN-IPSL (2005)  
     47   !!   OPA 9.0 , LOCEAN-IPSL (2006)  
     48   !! $Header$  
     49   !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt) 
    4850   !!---------------------------------------------------------------------- 
    4951 
     
    5658      !! ** Purpose :   compute the ocean tracer advection trend. 
    5759      !! 
     60      !! ** Method  : - Update (ua,va) with the advection term following nadv 
    5861      !!---------------------------------------------------------------------- 
    5962#if ( defined key_trabbl_adv || defined key_traldf_eiv ) 
    60       REAL(wp), DIMENSION(jpi,jpj,jpk) ::  &  ! temporary arrays 
    61          &         zun, zvn, zwn 
     63      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zun, zvn, zwn   ! effective velocity 
    6264#else 
    63       USE oce                , zun => un,  &  ! When no advective bbl, zun == un 
    64          &                     zvn => vn,  &  !      "        "      , zvn == vn 
    65          &                     zwn => wn      !      "        "      , zwn == wn 
     65      USE oce, ONLY :                       zun => un       ! the effective velocity is the 
     66      USE oce, ONLY :                       zvn => vn       ! Eulerian velocity 
     67      USE oce, ONLY :                       zwn => wn       !  
    6668#endif 
    67  
    68       !! * Arguments 
    69       INTEGER, INTENT( in ) ::   kt     ! ocean time-step index 
     69      !! 
     70      INTEGER, INTENT( in ) ::   kt   ! ocean time-step index 
    7071      !!---------------------------------------------------------------------- 
    7172 
    72       IF( kt == nit000 )   CALL tra_adv_ctl                  ! initialisation & control of options 
     73      IF( kt == nit000 )   CALL tra_adv_ctl          ! initialisation & control of options 
    7374 
    7475#if defined key_trabbl_adv 
    75       ! Advective bottom boundary layer                      ! add the bbl velocity 
    76       ! ------------------------------- 
    77       zun(:,:,:) = un(:,:,:) - u_bbl(:,:,:) 
     76      zun(:,:,:) = un(:,:,:) - u_bbl(:,:,:)          ! add the bbl velocity 
    7877      zvn(:,:,:) = vn(:,:,:) - v_bbl(:,:,:) 
    7978      zwn(:,:,:) = wn(:,:,:) + w_bbl(:,:,:) 
    8079#endif 
    81       IF( lk_traldf_eiv ) THEN                               ! add the eiv velocity 
     80      IF( lk_traldf_eiv ) THEN                       ! commpute and add the eiv velocity 
    8281         IF( .NOT. lk_trabbl_adv ) THEN  
    8382            zun(:,:,:) = un(:,:,:) 
     
    8584            zwn(:,:,:) = wn(:,:,:) 
    8685         ENDIF 
    87          CALL tra_adv_eiv( kt, zun, zvn, zwn )                  ! compute and add the eiv velocity  
     86         CALL tra_adv_eiv( kt, zun, zvn, zwn )  
    8887      ENDIF 
    8988 
    90       SELECT CASE ( nadv )                       ! compute advection trend and add it to general trend 
    91       CASE ( -1 )                                       ! esopa: test all possibility with control print 
    92          CALL tra_adv_cen2    ( kt, zun, zvn, zwn ) 
    93          CALL prt_ctl( tab3d_1=ta, clinfo1=' ldf0 - Ta: ', mask1=tmask,               & 
    94             &          tab3d_2=sa, clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' ) 
    95          CALL tra_adv_cen2_jki( kt, zun, zvn, zwn ) 
    96          CALL prt_ctl( tab3d_1=ta, clinfo1=' ldf1 - Ta: ', mask1=tmask,               & 
    97             &          tab3d_2=sa, clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' ) 
    98          CALL tra_adv_tvd     ( kt, zun, zvn, zwn ) 
    99          CALL prt_ctl( tab3d_1=ta, clinfo1=' ldf2 - Ta: ', mask1=tmask,               & 
    100             &          tab3d_2=sa, clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' ) 
    101          CALL tra_adv_muscl   ( kt, zun, zvn, zwn ) 
    102          CALL prt_ctl( tab3d_1=ta, clinfo1=' ldf3 - Ta: ', mask1=tmask,               & 
    103             &          tab3d_2=sa, clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' ) 
    104          CALL tra_adv_muscl2  ( kt, zun, zvn, zwn ) 
    105          CALL prt_ctl( tab3d_1=ta, clinfo1=' ldf4 - Ta: ', mask1=tmask,               & 
    106             &          tab3d_2=sa, clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' ) 
    107  
    108       CASE ( 0 )                                       ! 2nd order centered scheme k-j-i loops 
    109          CALL tra_adv_cen2    ( kt, zun, zvn, zwn ) 
    110       CASE ( 1 )                                       ! 2nd order centered scheme 
    111          CALL tra_adv_cen2_jki( kt, zun, zvn, zwn ) 
    112       CASE ( 2 )                                       ! TVD scheme 
    113          CALL tra_adv_tvd    ( kt, zun, zvn, zwn ) 
    114       CASE ( 3 )                                       ! MUSCL  scheme 
    115          CALL tra_adv_muscl  ( kt, zun, zvn, zwn ) 
    116       CASE ( 4 )                                       ! MUSCL2 scheme 
    117          CALL tra_adv_muscl2 ( kt, zun, zvn, zwn ) 
     89      SELECT CASE ( nadv )                           ! compute advection trend and add it to general trend 
     90      CASE ( 0 )   ;   CALL tra_adv_cen2    ( kt, zun, zvn, zwn )    ! 2nd order centered scheme k-j-i loops 
     91      CASE ( 1 )   ;   CALL tra_adv_cen2_jki( kt, zun, zvn, zwn )    ! 2nd order centered scheme 
     92      CASE ( 2 )   ;   CALL tra_adv_tvd     ( kt, zun, zvn, zwn )    ! TVD    scheme 
     93      CASE ( 3 )   ;   CALL tra_adv_muscl   ( kt, zun, zvn, zwn )    ! MUSCL  scheme 
     94      CASE ( 4 )   ;   CALL tra_adv_muscl2  ( kt, zun, zvn, zwn )    ! MUSCL2 scheme 
     95      CASE ( 5 )   ;   CALL tra_adv_ubs     ( kt, zun, zvn, zwn )    ! UBS    scheme 
     96      ! 
     97      CASE (-1 )                                                     ! esopa: test all possibility with control print 
     98         ;             CALL tra_adv_cen2    ( kt, zun, zvn, zwn ) 
     99         ;             CALL prt_ctl( tab3d_1=ta, clinfo1=' ldf0 - Ta: ', mask1=tmask,               & 
     100            &                        tab3d_2=sa, clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' ) 
     101         ;             CALL tra_adv_cen2_jki( kt, zun, zvn, zwn ) 
     102         ;             CALL prt_ctl( tab3d_1=ta, clinfo1=' ldf1 - Ta: ', mask1=tmask,               & 
     103            &                        tab3d_2=sa, clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' ) 
     104         ;             CALL tra_adv_tvd     ( kt, zun, zvn, zwn ) 
     105         ;             CALL prt_ctl( tab3d_1=ta, clinfo1=' ldf2 - Ta: ', mask1=tmask,               & 
     106            &                        tab3d_2=sa, clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' ) 
     107         ;             CALL tra_adv_muscl   ( kt, zun, zvn, zwn ) 
     108         ;             CALL prt_ctl( tab3d_1=ta, clinfo1=' ldf3 - Ta: ', mask1=tmask,               & 
     109            &                        tab3d_2=sa, clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' ) 
     110         ;             CALL tra_adv_muscl2  ( kt, zun, zvn, zwn ) 
     111         ;             CALL prt_ctl( tab3d_1=ta, clinfo1=' ldf4 - Ta: ', mask1=tmask,               & 
     112            &                        tab3d_2=sa, clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' ) 
     113         ;             CALL tra_adv_ubs     ( kt, zun, zvn, zwn ) 
     114         ;             CALL prt_ctl( tab3d_1=ta, clinfo1=' adv5 - Ta: ', mask1=tmask,               & 
     115            &                        tab3d_2=sa, clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' ) 
    118116      END SELECT 
    119  
    120       !                                                ! print mean trends (used for debugging) 
    121 !      IF(ln_ctl)   CALL prt_ctl( tab3d_1=ta, clinfo1=' adv  - Ta: ', mask1=tmask,               & 
    122 !         &                       tab3d_2=sa, clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' ) 
    123  
     117      !                                              ! print mean trends (used for debugging) 
     118      IF(ln_ctl)   CALL prt_ctl( tab3d_1=ta, clinfo1=' adv  - Ta: ', mask1=tmask,               & 
     119         &                       tab3d_2=sa, clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' ) 
     120      ! 
    124121   END SUBROUTINE tra_adv 
    125122 
     
    129126      !!                  ***  ROUTINE tra_adv_ctl  *** 
    130127      !!                 
    131       !! ** Purpose :   Control the consistency between cpp options for  
    132       !!      tracer advection schemes 
    133       !! 
     128      !! ** Purpose :   Control the consistency between namelist options for  
     129      !!              tracer advection schemes and set nadv 
    134130      !!---------------------------------------------------------------------- 
    135131      INTEGER ::   ioptio 
    136       NAMELIST/nam_traadv/ ln_traadv_cen2 , ln_traadv_tvd,   & 
    137          &                 ln_traadv_muscl, ln_traadv_muscl2 
    138132      !!---------------------------------------------------------------------- 
    139133 
    140       ! Read Namelist nam_traadv : tracer advection scheme 
    141       ! ------------------------- 
    142       REWIND ( numnam ) 
     134      REWIND ( numnam )               ! Read Namelist nam_traadv : tracer advection scheme 
    143135      READ   ( numnam, nam_traadv ) 
    144136 
    145       ! Parameter control and print 
    146       ! --------------------------- 
    147       ! Control print 
    148       IF(lwp) THEN 
     137      IF(lwp) THEN                    ! Namelist print 
    149138         WRITE(numout,*) 
    150139         WRITE(numout,*) 'tra_adv_ctl : choice/control of the tracer advection scheme' 
    151140         WRITE(numout,*) '~~~~~~~~~~~' 
    152          WRITE(numout,*) '          Namelist nam_tra_adv : chose a advection scheme for tracers' 
    153          WRITE(numout,*) 
    154          WRITE(numout,*) '             2nd order advection scheme     ln_traadv_cen2   = ', ln_traadv_cen2 
    155          WRITE(numout,*) '             TVD advection scheme           ln_traadv_tvd    = ', ln_traadv_tvd 
    156          WRITE(numout,*) '             MUSCL  advection scheme        ln_traadv_muscl  = ', ln_traadv_muscl 
    157          WRITE(numout,*) '             MUSCL2 advection scheme        ln_traadv_muscl2 = ', ln_traadv_muscl2 
     141         WRITE(numout,*) '       Namelist nam_traadv : chose a advection scheme for tracers' 
     142         WRITE(numout,*) '          2nd order advection scheme     ln_traadv_cen2   = ', ln_traadv_cen2 
     143         WRITE(numout,*) '          TVD advection scheme           ln_traadv_tvd    = ', ln_traadv_tvd 
     144         WRITE(numout,*) '          MUSCL  advection scheme        ln_traadv_muscl  = ', ln_traadv_muscl 
     145         WRITE(numout,*) '          MUSCL2 advection scheme        ln_traadv_muscl2 = ', ln_traadv_muscl2 
     146         WRITE(numout,*) '          UBS    advection scheme        ln_traadv_ubs    = ', ln_traadv_ubs 
    158147      ENDIF 
    159148 
    160       ! Control of Advection scheme options 
    161       ! ----------------------------------- 
    162       ioptio = 0 
     149      ioptio = 0                      ! Parameter control 
    163150      IF( ln_traadv_cen2   )   ioptio = ioptio + 1 
    164151      IF( ln_traadv_tvd    )   ioptio = ioptio + 1 
    165152      IF( ln_traadv_muscl  )   ioptio = ioptio + 1 
    166153      IF( ln_traadv_muscl2 )   ioptio = ioptio + 1 
     154      IF( ln_traadv_ubs    )   ioptio = ioptio + 1 
     155      IF( lk_esopa         )   ioptio =          1 
    167156 
    168       IF( .NOT.lk_esopa .AND. ( ioptio > 1 .OR. ioptio == 0 ) ) & 
    169            &   CALL ctl_stop( ' Choose ONE advection scheme in namelist nam_traadv' ) 
     157      IF( ioptio /= 1 )   CALL ctl_stop( 'Choose ONE advection scheme in namelist nam_traadv' ) 
    170158 
    171       IF( n_cla == 1 .AND. .NOT. ln_traadv_cen2 ) & 
    172            &   CALL ctl_stop( '     cross-land advection only with 2nd order advection scheme' ) 
     159      IF( n_cla == 1 .AND. .NOT. ln_traadv_cen2 )   & 
     160         &                CALL ctl_stop( 'cross-land advection only with 2nd order advection scheme' ) 
    173161 
    174       ! Set nadv 
    175       ! -------- 
    176       IF( ln_traadv_cen2   )   nadv = 0 
     162      !                              ! Set nadv 
     163      IF( ln_traadv_cen2   )   nadv =  0 
    177164#if defined key_mpp_omp 
    178       IF( ln_traadv_cen2   )   nadv = 1 
     165      IF( ln_traadv_cen2   )   nadv =  1 
    179166#endif 
    180       IF( ln_traadv_tvd    )   nadv = 2 
    181       IF( ln_traadv_muscl  )   nadv = 3 
    182       IF( ln_traadv_muscl2 )   nadv = 4 
     167      IF( ln_traadv_tvd    )   nadv =  2 
     168      IF( ln_traadv_muscl  )   nadv =  3 
     169      IF( ln_traadv_muscl2 )   nadv =  4 
     170      IF( ln_traadv_ubs    )   nadv =  5 
     171      IF( lk_esopa         )   nadv = -1 
    183172 
    184       IF( lk_esopa ) THEN 
    185          IF(lwp) WRITE(numout,*) '          esopa control: use all lateral physics options' 
    186          nadv = -1 
     173      IF(lwp) THEN                   ! Print the choice 
     174         WRITE(numout,*) 
     175         IF( nadv ==  0 )   WRITE(numout,*) '         2nd order scheme is used' 
     176         IF( nadv ==  1 )   WRITE(numout,*) '         2nd order scheme is usedi, k-j-i case' 
     177         IF( nadv ==  2 )   WRITE(numout,*) '         TVD       scheme is used' 
     178         IF( nadv ==  3 )   WRITE(numout,*) '         MUSCL     scheme is used' 
     179         IF( nadv ==  4 )   WRITE(numout,*) '         MUSCL2    scheme is used' 
     180         IF( nadv ==  5 )   WRITE(numout,*) '         UBS       scheme is used' 
     181         IF( nadv == -1 )   WRITE(numout,*) '         esopa test: use all advection scheme' 
    187182      ENDIF 
    188       IF(lwp) WRITE(numout,*) '             choice of tra_adv_...                nadv   = ', nadv 
    189  
     183      ! 
    190184   END SUBROUTINE tra_adv_ctl 
    191185 
Note: See TracChangeset for help on using the changeset viewer.