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.
traadv.F90 in branches/2014/dev_r4743_NOC2_ZTS/NEMOGCM/NEMO/OPA_SRC/TRA – NEMO

source: branches/2014/dev_r4743_NOC2_ZTS/NEMOGCM/NEMO/OPA_SRC/TRA/traadv.F90 @ 4736

Last change on this file since 4736 was 4736, checked in by acc, 10 years ago

Branch 2014/dev_r4743_NOC2_ZTS, #1367. Code changes to introduce optional sub-timestepping for vertical advection. Changes in dynadv.F90, dynzad.F90, traadv.F90, traadv_tvd.F90 and namelist_ref only.

  • Property svn:keywords set to Id
File size: 13.6 KB
Line 
1MODULE traadv
2   !!==============================================================================
3   !!                       ***  MODULE  traadv  ***
4   !! Ocean active tracers:  advection trend
5   !!==============================================================================
6   !! History :  2.0  !  2005-11  (G. Madec)  Original code
7   !!            3.3  !  2010-09  (C. Ethe, G. Madec)  merge TRC-TRA + switch from velocity to transport
8   !!            4.0  !  2011-06  (G. Madec)  Addition of Mixed Layer Eddy parameterisation
9   !!----------------------------------------------------------------------
10
11   !!----------------------------------------------------------------------
12   !!   tra_adv      : compute ocean tracer advection trend
13   !!   tra_adv_ctl  : control the different options of advection scheme
14   !!----------------------------------------------------------------------
15   USE oce             ! ocean dynamics and active tracers
16   USE dom_oce         ! ocean space and time domain
17   USE domvvl          ! variable vertical scale factors
18   USE traadv_cen2     ! 2nd order centered scheme (tra_adv_cen2   routine)
19   USE traadv_tvd      ! TVD      scheme           (tra_adv_tvd    routine)
20   USE traadv_muscl    ! MUSCL    scheme           (tra_adv_muscl  routine)
21   USE traadv_muscl2   ! MUSCL2   scheme           (tra_adv_muscl2 routine)
22   USE traadv_ubs      ! UBS      scheme           (tra_adv_ubs    routine)
23   USE traadv_qck      ! QUICKEST scheme           (tra_adv_qck    routine)
24   USE traadv_eiv      ! eddy induced velocity     (tra_adv_eiv    routine)
25   USE traadv_mle      ! ML eddy induced velocity  (tra_adv_mle    routine)
26   USE cla             ! cross land advection      (cla_traadv     routine)
27   USE ldftra_oce      ! lateral diffusion coefficient on tracers
28   USE in_out_manager  ! I/O manager
29   USE iom             ! I/O module
30   USE prtctl          ! Print control
31   USE lib_mpp         ! MPP library
32   USE wrk_nemo        ! Memory Allocation
33   USE timing          ! Timing
34
35
36   IMPLICIT NONE
37   PRIVATE
38
39   PUBLIC   tra_adv        ! routine called by step module
40   PUBLIC   tra_adv_init   ! routine called by opa module
41
42   !                              !!* Namelist namtra_adv *
43   LOGICAL ::   ln_traadv_cen2     ! 2nd order centered scheme flag
44   LOGICAL ::   ln_traadv_tvd      ! TVD scheme flag
45   LOGICAL ::   ln_traadv_tvd_zts  ! TVD scheme flag with vertical sub time-stepping
46   LOGICAL ::   ln_traadv_muscl    ! MUSCL scheme flag
47   LOGICAL ::   ln_traadv_muscl2   ! MUSCL2 scheme flag
48   LOGICAL ::   ln_traadv_ubs      ! UBS scheme flag
49   LOGICAL ::   ln_traadv_qck      ! QUICKEST scheme flag
50   LOGICAL ::   ln_traadv_msc_ups  ! use upstream scheme within muscl
51
52
53   INTEGER ::   nadv   ! choice of the type of advection scheme
54
55   !! * Substitutions
56#  include "domzgr_substitute.h90"
57#  include "vectopt_loop_substitute.h90"
58   !!----------------------------------------------------------------------
59   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
60   !! $Id$
61   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
62   !!----------------------------------------------------------------------
63CONTAINS
64
65   SUBROUTINE tra_adv( kt )
66      !!----------------------------------------------------------------------
67      !!                  ***  ROUTINE tra_adv  ***
68      !!
69      !! ** Purpose :   compute the ocean tracer advection trend.
70      !!
71      !! ** Method  : - Update (ua,va) with the advection term following nadv
72      !!----------------------------------------------------------------------
73      !
74      INTEGER, INTENT( in ) ::   kt   ! ocean time-step index
75      !
76      INTEGER ::   jk   ! dummy loop index
77      REAL(wp), POINTER, DIMENSION(:,:,:) :: zun, zvn, zwn
78      !!----------------------------------------------------------------------
79      !
80      IF( nn_timing == 1 )  CALL timing_start('tra_adv')
81      !
82      CALL wrk_alloc( jpi, jpj, jpk, zun, zvn, zwn )
83      !                                          ! set time step
84      IF( neuler == 0 .AND. kt == nit000 ) THEN     ! at nit000
85         r2dtra(:) =  rdttra(:)                          ! = rdtra (restarting with Euler time stepping)
86      ELSEIF( kt <= nit000 + 1) THEN                ! at nit000 or nit000+1
87         r2dtra(:) = 2._wp * rdttra(:)                   ! = 2 rdttra (leapfrog)
88      ENDIF
89      !
90      IF( nn_cla == 1 .AND. cp_cfg == 'orca' .AND. jp_cfg == 2 )   CALL cla_traadv( kt )       !==  Cross Land Advection  ==! (hor. advection)
91      !
92      !                                               !==  effective transport  ==!
93      DO jk = 1, jpkm1
94         zun(:,:,jk) = e2u(:,:) * fse3u(:,:,jk) * un(:,:,jk)                  ! eulerian transport only
95         zvn(:,:,jk) = e1v(:,:) * fse3v(:,:,jk) * vn(:,:,jk)
96         zwn(:,:,jk) = e1t(:,:) * e2t(:,:)      * wn(:,:,jk)
97      END DO
98      !
99      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN
100         zun(:,:,:) = zun(:,:,:) + un_td(:,:,:)
101         zvn(:,:,:) = zvn(:,:,:) + vn_td(:,:,:)
102      ENDIF
103      !
104      zun(:,:,jpk) = 0._wp                                                     ! no transport trough the bottom
105      zvn(:,:,jpk) = 0._wp                                                     ! no transport trough the bottom
106      zwn(:,:,jpk) = 0._wp                                                     ! no transport trough the bottom
107      !
108      IF( lk_traldf_eiv .AND. .NOT. ln_traldf_grif )   &
109         &              CALL tra_adv_eiv( kt, nit000, zun, zvn, zwn, 'TRA' )    ! add the eiv transport (if necessary)
110      !
111      IF( ln_mle    )   CALL tra_adv_mle( kt, nit000, zun, zvn, zwn, 'TRA' )    ! add the mle transport (if necessary)
112      CALL iom_put( "uocetr_eff", zun )                                         ! output effective transport     
113      CALL iom_put( "vocetr_eff", zvn )
114      CALL iom_put( "wocetr_eff", zwn )
115
116      SELECT CASE ( nadv )                            !==  compute advection trend and add it to general trend  ==!
117      CASE ( 1 )   ;    CALL tra_adv_cen2  ( kt, nit000, 'TRA',         zun, zvn, zwn, tsb, tsn, tsa, jpts )   !  2nd order centered
118      CASE ( 2 )   ;    CALL tra_adv_tvd   ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts )   !  TVD
119      CASE ( 3 )   ;    CALL tra_adv_muscl ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb,      tsa, jpts, ln_traadv_msc_ups )   !  MUSCL
120      CASE ( 4 )   ;    CALL tra_adv_muscl2( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts )   !  MUSCL2
121      CASE ( 5 )   ;    CALL tra_adv_ubs   ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts )   !  UBS
122      CASE ( 6 )   ;    CALL tra_adv_qck   ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts )   !  QUICKEST
123      CASE ( 7 )   ;   CALL tra_adv_tvd_zts( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts )   !  TVD ZTS
124      !
125      CASE (-1 )                                      !==  esopa: test all possibility with control print  ==!
126         CALL tra_adv_cen2  ( kt, nit000, 'TRA',         zun, zvn, zwn, tsb, tsn, tsa, jpts )         
127         CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' adv0 - Ta: ', mask1=tmask,               &
128            &          tab3d_2=tsa(:,:,:,jp_sal), clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' )
129         CALL tra_adv_tvd   ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts )         
130         CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' adv1 - Ta: ', mask1=tmask,               &
131            &          tab3d_2=tsa(:,:,:,jp_sal), clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' )
132         CALL tra_adv_muscl ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb,      tsa, jpts, ln_traadv_msc_ups )         
133         CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' adv3 - Ta: ', mask1=tmask,               &
134            &          tab3d_2=tsa(:,:,:,jp_sal), clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' )
135         CALL tra_adv_muscl2( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts )         
136         CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' adv4 - Ta: ', mask1=tmask,               &
137            &          tab3d_2=tsa(:,:,:,jp_sal), clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' )
138         CALL tra_adv_ubs   ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts )         
139         CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' adv5 - Ta: ', mask1=tmask,               &
140            &          tab3d_2=tsa(:,:,:,jp_sal), clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' )
141         CALL tra_adv_qck   ( kt, nit000, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts )         
142         CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' adv6 - Ta: ', mask1=tmask,               &
143            &          tab3d_2=tsa(:,:,:,jp_sal), clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' )
144      END SELECT
145      !
146      !                                              ! print mean trends (used for debugging)
147      IF(ln_ctl)   CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' adv  - Ta: ', mask1=tmask,               &
148         &                       tab3d_2=tsa(:,:,:,jp_sal), clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' )
149      !
150      IF( nn_timing == 1 )  CALL timing_stop( 'tra_adv' )
151      !
152      CALL wrk_dealloc( jpi, jpj, jpk, zun, zvn, zwn )
153      !                                         
154   END SUBROUTINE tra_adv
155
156
157   SUBROUTINE tra_adv_init
158      !!---------------------------------------------------------------------
159      !!                  ***  ROUTINE tra_adv_init  ***
160      !!               
161      !! ** Purpose :   Control the consistency between namelist options for
162      !!              tracer advection schemes and set nadv
163      !!----------------------------------------------------------------------
164      INTEGER ::   ioptio
165      INTEGER ::   ios                 ! Local integer output status for namelist read
166      !!
167      NAMELIST/namtra_adv/ ln_traadv_cen2 , ln_traadv_tvd,     &
168         &                 ln_traadv_muscl, ln_traadv_muscl2,  &
169         &                 ln_traadv_ubs  , ln_traadv_qck,     &
170         &                 ln_traadv_msc_ups, ln_traadv_tvd_zts
171      !!----------------------------------------------------------------------
172
173      REWIND( numnam_ref )              ! Namelist namtra_adv in reference namelist : Tracer advection scheme
174      READ  ( numnam_ref, namtra_adv, IOSTAT = ios, ERR = 901)
175901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namtra_adv in reference namelist', lwp )
176
177      REWIND( numnam_cfg )              ! Namelist namtra_adv in configuration namelist : Tracer advection scheme
178      READ  ( numnam_cfg, namtra_adv, IOSTAT = ios, ERR = 902 )
179902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namtra_adv in configuration namelist', lwp )
180      IF(lwm) WRITE ( numond, namtra_adv )
181
182      IF(lwp) THEN                    ! Namelist print
183         WRITE(numout,*)
184         WRITE(numout,*) 'tra_adv_init : choice/control of the tracer advection scheme'
185         WRITE(numout,*) '~~~~~~~~~~~'
186         WRITE(numout,*) '   Namelist namtra_adv : chose a advection scheme for tracers'
187         WRITE(numout,*) '      2nd order advection scheme     ln_traadv_cen2    = ', ln_traadv_cen2
188         WRITE(numout,*) '      TVD advection scheme           ln_traadv_tvd     = ', ln_traadv_tvd
189         WRITE(numout,*) '      MUSCL  advection scheme        ln_traadv_muscl   = ', ln_traadv_muscl
190         WRITE(numout,*) '      MUSCL2 advection scheme        ln_traadv_muscl2  = ', ln_traadv_muscl2
191         WRITE(numout,*) '      UBS    advection scheme        ln_traadv_ubs     = ', ln_traadv_ubs
192         WRITE(numout,*) '      QUICKEST advection scheme      ln_traadv_qck     = ', ln_traadv_qck
193         WRITE(numout,*) '      upstream scheme within muscl   ln_traadv_msc_ups = ', ln_traadv_msc_ups
194         WRITE(numout,*) '      TVD advection scheme with zts  ln_traadv_tvd_zts = ', ln_traadv_tvd_zts
195      ENDIF
196
197      ioptio = 0                      ! Parameter control
198      IF( ln_traadv_cen2   )   ioptio = ioptio + 1
199      IF( ln_traadv_tvd    )   ioptio = ioptio + 1
200      IF( ln_traadv_muscl  )   ioptio = ioptio + 1
201      IF( ln_traadv_muscl2 )   ioptio = ioptio + 1
202      IF( ln_traadv_ubs    )   ioptio = ioptio + 1
203      IF( ln_traadv_qck    )   ioptio = ioptio + 1
204      IF( ln_traadv_tvd_zts)   ioptio = ioptio + 1
205      IF( lk_esopa         )   ioptio =          1
206
207      IF( ioptio /= 1 )   CALL ctl_stop( 'Choose ONE advection scheme in namelist namtra_adv' )
208
209      !                              ! Set nadv
210      IF( ln_traadv_cen2   )   nadv =  1
211      IF( ln_traadv_tvd    )   nadv =  2
212      IF( ln_traadv_muscl  )   nadv =  3
213      IF( ln_traadv_muscl2 )   nadv =  4
214      IF( ln_traadv_ubs    )   nadv =  5
215      IF( ln_traadv_qck    )   nadv =  6
216      IF( ln_traadv_tvd_zts)   nadv =  7
217      IF( lk_esopa         )   nadv = -1
218
219      IF(lwp) THEN                   ! Print the choice
220         WRITE(numout,*)
221         IF( nadv ==  1 )   WRITE(numout,*) '         2nd order scheme is used'
222         IF( nadv ==  2 )   WRITE(numout,*) '         TVD       scheme is used'
223         IF( nadv ==  3 )   WRITE(numout,*) '         MUSCL     scheme is used'
224         IF( nadv ==  4 )   WRITE(numout,*) '         MUSCL2    scheme is used'
225         IF( nadv ==  5 )   WRITE(numout,*) '         UBS       scheme is used'
226         IF( nadv ==  6 )   WRITE(numout,*) '         QUICKEST  scheme is used'
227         IF( nadv ==  7 )   WRITE(numout,*) '         TVD ZTS   scheme is used'
228         IF( nadv == -1 )   WRITE(numout,*) '         esopa test: use all advection scheme'
229      ENDIF
230      !
231      CALL tra_adv_mle_init          ! initialisation of the Mixed Layer Eddy parametrisation (MLE)
232      !
233   END SUBROUTINE tra_adv_init
234
235  !!======================================================================
236END MODULE traadv
Note: See TracBrowser for help on using the repository browser.