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 NEMO/branches/2020/dev_r12377_KERNEL-06_techene_e3/src/OCE/TRA – NEMO

source: NEMO/branches/2020/dev_r12377_KERNEL-06_techene_e3/src/OCE/TRA/traadv.F90 @ 12624

Last change on this file since 12624 was 12624, checked in by techene, 4 years ago

all: add e3 substitute and limit precompiled files lines to about 130 character, change key_LF into key_QCO, change module name (dynatfQCO, traatfQCO, stepLF)

  • Property svn:keywords set to Id
File size: 15.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   !!            3.6  !  2011-06  (G. Madec)  Addition of Mixed Layer Eddy parameterisation
9   !!            3.7  !  2014-05  (G. Madec)  Add 2nd/4th order cases for CEN and FCT schemes
10   !!             -   !  2014-12  (G. Madec) suppression of cross land advection option
11   !!            3.6  !  2015-06  (E. Clementi) Addition of Stokes drift in case of wave coupling
12   !!----------------------------------------------------------------------
13
14   !!----------------------------------------------------------------------
15   !!   tra_adv       : compute ocean tracer advection trend
16   !!   tra_adv_init  : control the different options of advection scheme
17   !!----------------------------------------------------------------------
18   USE oce            ! ocean dynamics and active tracers
19   USE dom_oce        ! ocean space and time domain
20   USE domvvl         ! variable vertical scale factors
21   USE sbcwave        ! wave module
22   USE sbc_oce        ! surface boundary condition: ocean
23   USE traadv_cen     ! centered scheme            (tra_adv_cen  routine)
24   USE traadv_fct     ! FCT      scheme            (tra_adv_fct  routine)
25   USE traadv_mus     ! MUSCL    scheme            (tra_adv_mus  routine)
26   USE traadv_ubs     ! UBS      scheme            (tra_adv_ubs  routine)
27   USE traadv_qck     ! QUICKEST scheme            (tra_adv_qck  routine)
28   USE tramle         ! Mixed Layer Eddy transport (tra_mle_trp  routine)
29   USE ldftra         ! Eddy Induced transport     (ldf_eiv_trp  routine)
30   USE ldfslp         ! Lateral diffusion: slopes of neutral surfaces
31   USE trd_oce        ! trends: ocean variables
32   USE trdtra         ! trends manager: tracers
33   USE diaptr         ! Poleward heat transport
34   !
35   USE in_out_manager ! I/O manager
36   USE iom            ! I/O module
37   USE prtctl         ! Print control
38   USE lib_mpp        ! MPP library
39   USE timing         ! Timing
40
41   IMPLICIT NONE
42   PRIVATE
43
44   PUBLIC   tra_adv        ! called by step.F90
45   PUBLIC   tra_adv_init   ! called by nemogcm.F90
46
47   !                            !!* Namelist namtra_adv *
48   LOGICAL ::   ln_traadv_OFF    ! no advection on T and S
49   LOGICAL ::   ln_traadv_cen    ! centered scheme flag
50   INTEGER ::      nn_cen_h, nn_cen_v   ! =2/4 : horizontal and vertical choices of the order of CEN scheme
51   LOGICAL ::   ln_traadv_fct    ! FCT scheme flag
52   INTEGER ::      nn_fct_h, nn_fct_v   ! =2/4 : horizontal and vertical choices of the order of FCT scheme
53   LOGICAL ::   ln_traadv_mus    ! MUSCL scheme flag
54   LOGICAL ::      ln_mus_ups           ! use upstream scheme in vivcinity of river mouths
55   LOGICAL ::   ln_traadv_ubs    ! UBS scheme flag
56   INTEGER ::      nn_ubs_v             ! =2/4 : vertical choice of the order of UBS scheme
57   LOGICAL ::   ln_traadv_qck    ! QUICKEST scheme flag
58
59   INTEGER ::   nadv             ! choice of the type of advection scheme
60   !                             ! associated indices:
61   INTEGER, PARAMETER ::   np_NO_adv  = 0   ! no T-S advection
62   INTEGER, PARAMETER ::   np_CEN     = 1   ! 2nd/4th order centered scheme
63   INTEGER, PARAMETER ::   np_FCT     = 2   ! 2nd/4th order Flux Corrected Transport scheme
64   INTEGER, PARAMETER ::   np_MUS     = 3   ! MUSCL scheme
65   INTEGER, PARAMETER ::   np_UBS     = 4   ! 3rd order Upstream Biased Scheme
66   INTEGER, PARAMETER ::   np_QCK     = 5   ! QUICK scheme
67
68#  include "domzgr_substitute.h90"
69   !!----------------------------------------------------------------------
70   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
71   !! $Id$
72   !! Software governed by the CeCILL license (see ./LICENSE)
73   !!----------------------------------------------------------------------
74CONTAINS
75
76   SUBROUTINE tra_adv( kt, Kbb, Kmm, pts, Krhs )
77      !!----------------------------------------------------------------------
78      !!                  ***  ROUTINE tra_adv  ***
79      !!
80      !! ** Purpose :   compute the ocean tracer advection trend.
81      !!
82      !! ** Method  : - Update (uu(:,:,:,Krhs),vv(:,:,:,Krhs)) with the advection term following nadv
83      !!----------------------------------------------------------------------
84      INTEGER                                  , INTENT(in)    :: kt             ! ocean time-step index
85      INTEGER                                  , INTENT(in)    :: Kbb, Kmm, Krhs ! time level indices
86      REAL(wp), DIMENSION(jpi,jpj,jpk,jpts,jpt), INTENT(inout) :: pts            ! active tracers and RHS of tracer equation
87      !
88      INTEGER ::   jk   ! dummy loop index
89      REAL(wp), DIMENSION(jpi,jpj,jpk)        :: zuu, zvv, zww   ! 3D workspace
90      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::   ztrdt, ztrds
91      !!----------------------------------------------------------------------
92      !
93      IF( ln_timing )   CALL timing_start('tra_adv')
94      !
95      !                                          ! set time step
96      IF( neuler == 0 .AND. kt == nit000 ) THEN   ;   r2dt =         rdt   ! at nit000             (Euler)
97      ELSEIF( kt <= nit000 + 1 )           THEN   ;   r2dt = 2._wp * rdt   ! at nit000 or nit000+1 (Leapfrog)
98      ENDIF
99      !
100      !                                         !==  effective transport  ==!
101      zuu(:,:,jpk) = 0._wp
102      zvv(:,:,jpk) = 0._wp
103      zww(:,:,jpk) = 0._wp
104      IF( ln_wave .AND. ln_sdw )  THEN
105         DO jk = 1, jpkm1                                                       ! eulerian transport + Stokes Drift
106            zuu(:,:,jk) =   &
107               &  e2u  (:,:) * e3u(:,:,jk,Kmm) * ( uu(:,:,jk,Kmm) + usd(:,:,jk) )
108            zvv(:,:,jk) =   & 
109               &  e1v  (:,:) * e3v(:,:,jk,Kmm) * ( vv(:,:,jk,Kmm) + vsd(:,:,jk) )
110            zww(:,:,jk) =   & 
111               &  e1e2t(:,:)                 * ( ww(:,:,jk) + wsd(:,:,jk) )
112         END DO
113      ELSE
114         DO jk = 1, jpkm1
115            zuu(:,:,jk) = e2u  (:,:) * e3u(:,:,jk,Kmm) * uu(:,:,jk,Kmm)               ! eulerian transport only
116            zvv(:,:,jk) = e1v  (:,:) * e3v(:,:,jk,Kmm) * vv(:,:,jk,Kmm)
117            zww(:,:,jk) = e1e2t(:,:)                 * ww(:,:,jk)
118         END DO
119      ENDIF
120      !
121      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN                                ! add z-tilde and/or vvl corrections
122         zuu(:,:,:) = zuu(:,:,:) + un_td(:,:,:)
123         zvv(:,:,:) = zvv(:,:,:) + vn_td(:,:,:)
124      ENDIF
125      !
126      zuu(:,:,jpk) = 0._wp                                                      ! no transport trough the bottom
127      zvv(:,:,jpk) = 0._wp
128      zww(:,:,jpk) = 0._wp
129      !
130      IF( ln_ldfeiv .AND. .NOT. ln_traldf_triad )   &
131         &              CALL ldf_eiv_trp( kt, nit000, zuu, zvv, zww, 'TRA', Kmm, Krhs )   ! add the eiv transport (if necessary)
132      !
133      IF( ln_mle    )   CALL tra_mle_trp( kt, nit000, zuu, zvv, zww, 'TRA', Kmm       )   ! add the mle transport (if necessary)
134      !
135      CALL iom_put( "uocetr_eff", zuu )                                        ! output effective transport     
136      CALL iom_put( "vocetr_eff", zvv )
137      CALL iom_put( "wocetr_eff", zww )
138      !
139!!gm ???
140      CALL dia_ptr( kt, Kmm, zvv )                                    ! diagnose the effective MSF
141!!gm ???
142      !
143
144      IF( l_trdtra )   THEN                    !* Save ta and sa trends
145         ALLOCATE( ztrdt(jpi,jpj,jpk), ztrds(jpi,jpj,jpk) )
146         ztrdt(:,:,:) = pts(:,:,:,jp_tem,Krhs)
147         ztrds(:,:,:) = pts(:,:,:,jp_sal,Krhs)
148      ENDIF
149      !
150      SELECT CASE ( nadv )                      !==  compute advection trend and add it to general trend  ==!
151      !
152      CASE ( np_CEN )                                 ! Centered scheme : 2nd / 4th order
153         CALL tra_adv_cen    ( kt, nit000, 'TRA',         zuu, zvv, zww, Kmm, pts, jpts, Krhs, nn_cen_h, nn_cen_v )
154      CASE ( np_FCT )                                 ! FCT scheme      : 2nd / 4th order
155         CALL tra_adv_fct    ( kt, nit000, 'TRA', r2dt, zuu, zvv, zww, Kbb, Kmm, pts, jpts, Krhs, nn_fct_h, nn_fct_v )
156      CASE ( np_MUS )                                 ! MUSCL
157         CALL tra_adv_mus    ( kt, nit000, 'TRA', r2dt, zuu, zvv, zww, Kbb, Kmm, pts, jpts, Krhs, ln_mus_ups ) 
158      CASE ( np_UBS )                                 ! UBS
159         CALL tra_adv_ubs    ( kt, nit000, 'TRA', r2dt, zuu, zvv, zww, Kbb, Kmm, pts, jpts, Krhs, nn_ubs_v   )
160      CASE ( np_QCK )                                 ! QUICKEST
161         CALL tra_adv_qck    ( kt, nit000, 'TRA', r2dt, zuu, zvv, zww, Kbb, Kmm, pts, jpts, Krhs )
162      !
163      END SELECT
164      !
165      IF( l_trdtra )   THEN                      ! save the advective trends for further diagnostics
166         DO jk = 1, jpkm1
167            ztrdt(:,:,jk) = pts(:,:,jk,jp_tem,Krhs) - ztrdt(:,:,jk)
168            ztrds(:,:,jk) = pts(:,:,jk,jp_sal,Krhs) - ztrds(:,:,jk)
169         END DO
170         CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_tem, jptra_totad, ztrdt )
171         CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_sal, jptra_totad, ztrds )
172         DEALLOCATE( ztrdt, ztrds )
173      ENDIF
174      !                                              ! print mean trends (used for debugging)
175      IF(sn_cfctl%l_prtctl)   CALL prt_ctl( tab3d_1=pts(:,:,:,jp_tem,Krhs), clinfo1=' adv  - Ta: ', mask1=tmask,               &
176         &                                  tab3d_2=pts(:,:,:,jp_sal,Krhs), clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' )
177      !
178      IF( ln_timing )   CALL timing_stop( 'tra_adv' )
179      !
180   END SUBROUTINE tra_adv
181
182
183   SUBROUTINE tra_adv_init
184      !!---------------------------------------------------------------------
185      !!                  ***  ROUTINE tra_adv_init  ***
186      !!               
187      !! ** Purpose :   Control the consistency between namelist options for
188      !!              tracer advection schemes and set nadv
189      !!----------------------------------------------------------------------
190      INTEGER ::   ioptio, ios   ! Local integers
191      !
192      NAMELIST/namtra_adv/ ln_traadv_OFF,                        &   ! No advection
193         &                 ln_traadv_cen , nn_cen_h, nn_cen_v,   &   ! CEN
194         &                 ln_traadv_fct , nn_fct_h, nn_fct_v,   &   ! FCT
195         &                 ln_traadv_mus , ln_mus_ups,           &   ! MUSCL
196         &                 ln_traadv_ubs ,           nn_ubs_v,   &   ! UBS
197         &                 ln_traadv_qck                             ! QCK
198      !!----------------------------------------------------------------------
199      !
200      !                                !==  Namelist  ==!
201      READ  ( numnam_ref, namtra_adv, IOSTAT = ios, ERR = 901)
202901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namtra_adv in reference namelist' )
203      !
204      READ  ( numnam_cfg, namtra_adv, IOSTAT = ios, ERR = 902 )
205902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namtra_adv in configuration namelist' )
206      IF(lwm) WRITE( numond, namtra_adv )
207      !
208      IF(lwp) THEN                           ! Namelist print
209         WRITE(numout,*)
210         WRITE(numout,*) 'tra_adv_init : choice/control of the tracer advection scheme'
211         WRITE(numout,*) '~~~~~~~~~~~~'
212         WRITE(numout,*) '   Namelist namtra_adv : chose a advection scheme for tracers'
213         WRITE(numout,*) '      No advection on T & S                     ln_traadv_OFF = ', ln_traadv_OFF
214         WRITE(numout,*) '      centered scheme                           ln_traadv_cen = ', ln_traadv_cen
215         WRITE(numout,*) '            horizontal 2nd/4th order               nn_cen_h   = ', nn_fct_h
216         WRITE(numout,*) '            vertical   2nd/4th order               nn_cen_v   = ', nn_fct_v
217         WRITE(numout,*) '      Flux Corrected Transport scheme           ln_traadv_fct = ', ln_traadv_fct
218         WRITE(numout,*) '            horizontal 2nd/4th order               nn_fct_h   = ', nn_fct_h
219         WRITE(numout,*) '            vertical   2nd/4th order               nn_fct_v   = ', nn_fct_v
220         WRITE(numout,*) '      MUSCL scheme                              ln_traadv_mus = ', ln_traadv_mus
221         WRITE(numout,*) '            + upstream scheme near river mouths    ln_mus_ups = ', ln_mus_ups
222         WRITE(numout,*) '      UBS scheme                                ln_traadv_ubs = ', ln_traadv_ubs
223         WRITE(numout,*) '            vertical   2nd/4th order               nn_ubs_v   = ', nn_ubs_v
224         WRITE(numout,*) '      QUICKEST scheme                           ln_traadv_qck = ', ln_traadv_qck
225      ENDIF
226      !
227      !                                !==  Parameter control & set nadv ==!
228      ioptio = 0                       
229      IF( ln_traadv_OFF ) THEN   ;   ioptio = ioptio + 1   ;   nadv = np_NO_adv   ;   ENDIF
230      IF( ln_traadv_cen ) THEN   ;   ioptio = ioptio + 1   ;   nadv = np_CEN      ;   ENDIF
231      IF( ln_traadv_fct ) THEN   ;   ioptio = ioptio + 1   ;   nadv = np_FCT      ;   ENDIF
232      IF( ln_traadv_mus ) THEN   ;   ioptio = ioptio + 1   ;   nadv = np_MUS      ;   ENDIF
233      IF( ln_traadv_ubs ) THEN   ;   ioptio = ioptio + 1   ;   nadv = np_UBS      ;   ENDIF
234      IF( ln_traadv_qck ) THEN   ;   ioptio = ioptio + 1   ;   nadv = np_QCK      ;   ENDIF
235      !
236      IF( ioptio /= 1 )   CALL ctl_stop( 'tra_adv_init: Choose ONE advection option in namelist namtra_adv' )
237      !
238      IF( ln_traadv_cen .AND. ( nn_cen_h /= 2 .AND. nn_cen_h /= 4 )   &          ! Centered
239                        .AND. ( nn_cen_v /= 2 .AND. nn_cen_v /= 4 )   ) THEN
240        CALL ctl_stop( 'tra_adv_init: CEN scheme, choose 2nd or 4th order' )
241      ENDIF
242      IF( ln_traadv_fct .AND. ( nn_fct_h /= 2 .AND. nn_fct_h /= 4 )   &          ! FCT
243                        .AND. ( nn_fct_v /= 2 .AND. nn_fct_v /= 4 )   ) THEN
244        CALL ctl_stop( 'tra_adv_init: FCT scheme, choose 2nd or 4th order' )
245      ENDIF
246      IF( ln_traadv_ubs .AND. ( nn_ubs_v /= 2 .AND. nn_ubs_v /= 4 )   ) THEN     ! UBS
247        CALL ctl_stop( 'tra_adv_init: UBS scheme, choose 2nd or 4th order' )
248      ENDIF
249      IF( ln_traadv_ubs .AND. nn_ubs_v == 4 ) THEN
250         CALL ctl_warn( 'tra_adv_init: UBS scheme, only 2nd FCT scheme available on the vertical. It will be used' )
251      ENDIF
252      IF( ln_isfcav ) THEN                                                       ! ice-shelf cavities
253         IF(  ln_traadv_cen .AND. nn_cen_v == 4    .OR.   &                            ! NO 4th order with ISF
254            & ln_traadv_fct .AND. nn_fct_v == 4   )   CALL ctl_stop( 'tra_adv_init: 4th order COMPACT scheme not allowed with ISF' )
255      ENDIF
256      !
257      !                                !==  Print the choice  ==! 
258      IF(lwp) THEN
259         WRITE(numout,*)
260         SELECT CASE ( nadv )
261         CASE( np_NO_adv  )   ;   WRITE(numout,*) '   ==>>>   NO T-S advection'
262         CASE( np_CEN     )   ;   WRITE(numout,*) '   ==>>>   CEN      scheme is used. Horizontal order: ', nn_cen_h,   &
263            &                                                                        ' Vertical   order: ', nn_cen_v
264         CASE( np_FCT     )   ;   WRITE(numout,*) '   ==>>>   FCT      scheme is used. Horizontal order: ', nn_fct_h,   &
265            &                                                                        ' Vertical   order: ', nn_fct_v
266         CASE( np_MUS     )   ;   WRITE(numout,*) '   ==>>>   MUSCL    scheme is used'
267         CASE( np_UBS     )   ;   WRITE(numout,*) '   ==>>>   UBS      scheme is used'
268         CASE( np_QCK     )   ;   WRITE(numout,*) '   ==>>>   QUICKEST scheme is used'
269         END SELECT
270      ENDIF
271      !
272      CALL tra_mle_init            !== initialisation of the Mixed Layer Eddy parametrisation (MLE)  ==!
273      !
274   END SUBROUTINE tra_adv_init
275
276  !!======================================================================
277END MODULE traadv
Note: See TracBrowser for help on using the repository browser.