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/UKMO/dev_r5107_kara_mld/NEMOGCM/NEMO/OPA_SRC/TRA – NEMO

source: branches/UKMO/dev_r5107_kara_mld/NEMOGCM/NEMO/OPA_SRC/TRA/traadv.F90 @ 5244

Last change on this file since 5244 was 5244, checked in by davestorkey, 9 years ago

UKMO Kara MLD branch: remove svn keyword property and clear keywords.

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