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.
trcadv.F90 in NEMO/branches/2021/dev_r14318_RK3_stage1/src/TOP/TRP – NEMO

source: NEMO/branches/2021/dev_r14318_RK3_stage1/src/TOP/TRP/trcadv.F90 @ 15510

Last change on this file since 15510 was 15510, checked in by techene, 3 years ago

#2715 RK3 & TOP: musl adv scheme is only called at stage 3

  • Property svn:keywords set to Id
File size: 16.6 KB
Line 
1MODULE trcadv
2   !!==============================================================================
3   !!                       ***  MODULE  trcadv  ***
4   !! Ocean passive tracers:  advection trend
5   !!==============================================================================
6   !! History :  2.0  !  2005-11  (G. Madec)  Original code
7   !!            3.0  !  2010-06  (C. Ethe)   Adapted to passive tracers
8   !!            3.7  !  2014-05  (G. Madec, C. Ethe)  Add 2nd/4th order cases for CEN and FCT schemes
9   !!            4.0  !  2017-09  (G. Madec)  remove vertical time-splitting option
10   !!            4.5  !  2021-08  (G. Madec, S. Techene) add advective velocities as optional arguments
11   !!----------------------------------------------------------------------
12#if defined key_top
13   !!----------------------------------------------------------------------
14   !!   'key_top'                                                TOP models
15   !!----------------------------------------------------------------------
16   !!   trc_adv       : compute ocean tracer advection trend
17   !!   trc_adv_ini   : control the different options of advection scheme
18   !!----------------------------------------------------------------------
19   USE par_trc        ! need jptra, number of passive tracers
20   USE oce_trc        ! ocean dynamics and active tracers
21   USE trc            ! ocean passive tracers variables
22   USE sbcwave        ! wave module
23   USE sbc_oce        ! surface boundary condition: ocean
24   USE traadv_cen     ! centered scheme           (tra_adv_cen  routine)
25   USE traadv_fct     ! FCT      scheme           (tra_adv_fct  routine)
26   USE traadv_fct_lf  ! FCT      scheme           (tra_adv_fct  routine - loop fusion version)
27   USE traadv_mus     ! MUSCL    scheme           (tra_adv_mus  routine)
28   USE traadv_mus_lf  ! MUSCL    scheme           (tra_adv_mus  routine - loop fusion version)
29   USE traadv_ubs     ! UBS      scheme           (tra_adv_ubs  routine)
30   USE traadv_qck     ! QUICKEST scheme           (tra_adv_qck  routine)
31   USE tramle         ! ML eddy induced transport (tra_adv_mle  routine)
32   USE ldftra         ! lateral diffusion: eddy diffusivity & EIV coeff.
33   USE ldfslp         ! Lateral diffusion: slopes of neutral surfaces
34   !
35   USE prtctl         ! control print
36   USE timing         ! Timing
37
38   IMPLICIT NONE
39   PRIVATE
40
41   PUBLIC   trc_adv       ! called by trctrp.F90 and stprk3_stg.F90
42   PUBLIC   trc_adv_ini   ! called by trcini.F90
43
44   !                            !!* Namelist namtrc_adv *
45   LOGICAL ::   ln_trcadv_OFF    ! no advection on passive tracers
46   LOGICAL ::   ln_trcadv_cen    ! centered scheme flag
47   INTEGER ::      nn_cen_h, nn_cen_v   ! =2/4 : horizontal and vertical choices of the order of CEN scheme
48   LOGICAL ::   ln_trcadv_fct    ! FCT scheme flag
49   INTEGER ::      nn_fct_h, nn_fct_v   ! =2/4 : horizontal and vertical choices of the order of FCT scheme
50   LOGICAL, PUBLIC ::   ln_trcadv_mus    ! MUSCL scheme flag
51   LOGICAL ::      ln_mus_ups           ! use upstream scheme in vivcinity of river mouths
52   LOGICAL ::   ln_trcadv_ubs    ! UBS scheme flag
53   INTEGER ::      nn_ubs_v             ! =2/4 : vertical choice of the order of UBS scheme
54   LOGICAL ::   ln_trcadv_qck    ! QUICKEST scheme flag
55
56   INTEGER ::   nadv             ! choice of the type of advection scheme
57   !                             ! associated indices:
58   INTEGER, PARAMETER ::   np_NO_adv  = 0   ! no T-S advection
59   INTEGER, PARAMETER ::   np_CEN     = 1   ! 2nd/4th order centered scheme
60   INTEGER, PARAMETER ::   np_FCT     = 2   ! 2nd/4th order Flux Corrected Transport scheme
61   INTEGER, PARAMETER ::   np_MUS     = 3   ! MUSCL scheme
62   INTEGER, PARAMETER ::   np_UBS     = 4   ! 3rd order Upstream Biased Scheme
63   INTEGER, PARAMETER ::   np_QCK     = 5   ! QUICK scheme
64
65   !! * Substitutions
66#  include "do_loop_substitute.h90"
67#  include "domzgr_substitute.h90"
68   !!----------------------------------------------------------------------
69   !! NEMO/TOP 4.0 , NEMO Consortium (2018)
70   !! $Id$
71   !! Software governed by the CeCILL license (see ./LICENSE)
72   !!----------------------------------------------------------------------
73CONTAINS
74
75   SUBROUTINE trc_adv( kt, Kbb, Kmm, ptr, Krhs, pau, pav, paw )
76      !!----------------------------------------------------------------------
77      !!                  ***  ROUTINE trc_adv  ***
78      !!
79      !! ** Purpose :   compute the ocean tracer advection trend.
80      !!
81      !! ** Method  : - Update tr(Krhs) with the advective trend following nadv
82      !!----------------------------------------------------------------------
83      INTEGER                                     , INTENT(in   ) ::   kt             ! ocean time-step index
84      INTEGER                                     , INTENT(in   ) ::   Kbb, Kmm, Krhs ! time level indices
85      REAL(wp), DIMENSION(:,:,:), OPTIONAL, TARGET, INTENT(in   ) ::   pau, pav, paw  ! advective velocity
86      REAL(wp), DIMENSION(jpi,jpj,jpk,jptra,jpt)  , INTENT(inout) ::   ptr            ! passive tracers and RHS of tracer equation
87      !
88      INTEGER ::  ji, jj, jk   ! dummy loop index
89      CHARACTER (len=22) ::   charout
90      REAL(wp), DIMENSION(:,:,:), POINTER ::   zptu, zptv, zptw
91      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zuu, zvv, zww  ! effective velocity
92      !!----------------------------------------------------------------------
93      !
94      IF( ln_timing )   CALL timing_start('trc_adv')
95      !
96      !                                         !==  effective transport  ==!
97      IF( l_offline ) THEN
98         zuu(:,:,:) = uu(:,:,:,Kmm)                != already in (uu(Kmm),vv(Kmm),ww)
99         zvv(:,:,:) = vv(:,:,:,Kmm)
100         zww(:,:,:) = ww(:,:,:)
101      ELSE                                         != build the effective transport
102         zuu(:,:,jpk) = 0._wp                            ! no transport trough the bottom
103         zvv(:,:,jpk) = 0._wp
104         zww(:,:,jpk) = 0._wp
105         !
106         IF( PRESENT( pau ) ) THEN                       ! RK3: advective velocity (pau,pav,paw) /= advected velocity (uu,vv,ww)
107            zptu => pau(:,:,:)
108            zptv => pav(:,:,:)
109            zptw => paw(:,:,:)
110         ELSE                                            ! MLF: advective velocity = (uu,vv,ww)
111            zptu => uu(:,:,:,Kmm)
112            zptv => vv(:,:,:,Kmm)
113            zptw => ww(:,:,:    )
114         ENDIF
115         !
116         IF( ln_wave .AND. ln_sdw )  THEN                ! eulerian transport + Stokes Drift
117            DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, 1, jpkm1 )
118               zuu(ji,jj,jk) = e2u  (ji,jj) * e3u(ji,jj,jk,Kmm) * ( zptu(ji,jj,jk) + usd(ji,jj,jk) )
119               zvv(ji,jj,jk) = e1v  (ji,jj) * e3v(ji,jj,jk,Kmm) * ( zptv(ji,jj,jk) + vsd(ji,jj,jk) )
120               zww(ji,jj,jk) = e1e2t(ji,jj)                     * ( zptw(ji,jj,jk) + wsd(ji,jj,jk) )
121            END_3D
122         ELSE                                            ! eulerian transport only
123            DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, 1, jpkm1 )
124               zuu(ji,jj,jk) = e2u  (ji,jj) * e3u(ji,jj,jk,Kmm) * zptu(ji,jj,jk)
125               zvv(ji,jj,jk) = e1v  (ji,jj) * e3v(ji,jj,jk,Kmm) * zptv(ji,jj,jk)
126               zww(ji,jj,jk) = e1e2t(ji,jj)                     * zptw(ji,jj,jk)
127            END_3D
128         ENDIF
129         !
130         IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN      ! add z-tilde and/or vvl corrections
131            DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, 1, jpkm1 )
132               zuu(ji,jj,jk) = zuu(ji,jj,jk) + un_td(ji,jj,jk)
133               zvv(ji,jj,jk) = zvv(ji,jj,jk) + vn_td(ji,jj,jk)
134            END_3D
135         ENDIF
136         !
137         IF( ln_ldfeiv .AND. .NOT. ln_traldf_triad )   & 
138            &              CALL ldf_eiv_trp( kt, nittrc000, zuu, zvv, zww, 'TRC', Kmm, Krhs )  ! add the eiv transport
139         !
140         IF( ln_mle    )   CALL tra_mle_trp( kt, nittrc000, zuu, zvv, zww, 'TRC', Kmm       )  ! add the mle transport
141         !
142      ENDIF
143      !
144      SELECT CASE ( nadv )                      !==  compute advection trend and add it to general trend  ==!
145      !
146      CASE ( np_CEN )                                 ! Centered : 2nd / 4th order
147         IF( nn_hls == 2 ) CALL lbc_lnk( 'trcadv', ptr(:,:,:,:,Kmm), 'T', 1._wp )
148         CALL tra_adv_cen( kt, nittrc000,'TRC',          zuu, zvv, zww,      Kmm, ptr, jptra, Krhs, nn_cen_h, nn_cen_v )
149      CASE ( np_FCT )                                 ! FCT      : 2nd / 4th order
150         IF( nn_hls == 2 ) THEN
151            CALL lbc_lnk_multi( 'trcadv', ptr(:,:,:,:,Kbb), 'T', 1._wp, ptr(:,:,:,:,Kmm), 'T', 1._wp)
152            CALL lbc_lnk_multi( 'traadv', zuu(:,:,:), 'U', -1._wp, zvv(:,:,:), 'V', -1._wp, zww(:,:,:), 'W', 1._wp)
153#if defined key_loop_fusion
154            CALL tra_adv_fct_lf( kt, nittrc000,'TRC', rDt_trc, zuu, zvv, zww, Kbb, Kmm, ptr, jptra, Krhs, nn_fct_h, nn_fct_v )
155#else
156            CALL tra_adv_fct( kt, nittrc000,'TRC', rDt_trc, zuu, zvv, zww, Kbb, Kmm, ptr, jptra, Krhs, nn_fct_h, nn_fct_v )
157#endif
158         ELSE
159            CALL tra_adv_fct( kt, nittrc000,'TRC', rDt_trc, zuu, zvv, zww, Kbb, Kmm, ptr, jptra, Krhs, nn_fct_h, nn_fct_v )
160         END IF
161      CASE ( np_MUS )                                 ! MUSCL
162         IF( nn_hls == 2 ) THEN
163            CALL lbc_lnk( 'trcadv', ptr(:,:,:,:,Kbb), 'T', 1._wp)
164#if defined key_loop_fusion
165            CALL tra_adv_mus_lf( kt, nittrc000,'TRC', rDt_trc, zuu, zvv, zww, Kbb, Kmm, ptr, jptra, Krhs, ln_mus_ups ) 
166#else
167            CALL tra_adv_mus( kt, nittrc000,'TRC', rDt_trc, zuu, zvv, zww, Kbb, Kmm, ptr, jptra, Krhs, ln_mus_ups ) 
168#endif
169         ELSE
170            CALL tra_adv_mus( kt, nittrc000,'TRC', rDt_trc, zuu, zvv, zww, Kbb, Kmm, ptr, jptra, Krhs, ln_mus_ups ) 
171         END IF
172      CASE ( np_UBS )                                 ! UBS
173         IF( nn_hls == 2 )   CALL lbc_lnk( 'trcadv', ptr(:,:,:,:,Kbb), 'T', 1._wp)
174         CALL tra_adv_ubs( kt, nittrc000,'TRC', rDt_trc, zuu, zvv, zww, Kbb, Kmm, ptr, jptra, Krhs, nn_ubs_v           )
175      CASE ( np_QCK )                                 ! QUICKEST
176         IF( nn_hls == 2 ) THEN
177            CALL lbc_lnk_multi( 'trcadv', zuu(:,:,:), 'U', -1._wp, zvv(:,:,:), 'V', -1._wp)
178            CALL lbc_lnk( 'traadv', ptr(:,:,:,:,Kbb), 'T', 1._wp)
179         END IF
180         CALL tra_adv_qck( kt, nittrc000,'TRC', rDt_trc, zuu, zvv, zww, Kbb, Kmm, ptr, jptra, Krhs                     )
181      !
182      END SELECT
183      !                 
184      IF( sn_cfctl%l_prttrc ) THEN        !== print mean trends (used for debugging)
185         WRITE(charout, FMT="('adv ')")
186         CALL prt_ctl_info( charout, cdcomp = 'top' )
187         CALL prt_ctl( tab4d_1=tr(:,:,:,:,Krhs), mask1=tmask, clinfo=ctrcnm, clinfo3='trd' )
188      END IF
189      !
190      IF( ln_timing )   CALL timing_stop('trc_adv')
191      !
192   END SUBROUTINE trc_adv
193
194
195   SUBROUTINE trc_adv_ini
196      !!---------------------------------------------------------------------
197      !!                  ***  ROUTINE trc_adv_ini  ***
198      !!               
199      !! ** Purpose :   Control the consistency between namelist options for
200      !!              passive tracer advection schemes and set nadv
201      !!----------------------------------------------------------------------
202      INTEGER ::   ioptio, ios   ! Local integer
203      !!
204      NAMELIST/namtrc_adv/ ln_trcadv_OFF,                        &   ! No advection
205         &                 ln_trcadv_cen, nn_cen_h, nn_cen_v,    &   ! CEN
206         &                 ln_trcadv_fct, nn_fct_h, nn_fct_v,    &   ! FCT
207         &                 ln_trcadv_mus, ln_mus_ups,            &   ! MUSCL
208         &                 ln_trcadv_ubs,           nn_ubs_v,    &   ! UBS
209         &                 ln_trcadv_qck                             ! QCK
210      !!----------------------------------------------------------------------
211      !
212      !                                !==  Namelist  ==!
213      READ  ( numnat_ref, namtrc_adv, IOSTAT = ios, ERR = 901)
214901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namtrc_adv in reference namelist' )
215      READ  ( numnat_cfg, namtrc_adv, IOSTAT = ios, ERR = 902 )
216902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namtrc_adv in configuration namelist' )
217      IF(lwm) WRITE ( numont, namtrc_adv )
218      !
219      IF(lwp) THEN                           ! Namelist print
220         WRITE(numout,*)
221         WRITE(numout,*) 'trc_adv_ini : choice/control of the tracer advection scheme'
222         WRITE(numout,*) '~~~~~~~~~~~'
223         WRITE(numout,*) '   Namelist namtrc_adv : chose a advection scheme for tracers'
224         WRITE(numout,*) '      No advection on passive tracers           ln_trcadv_OFF = ', ln_trcadv_OFF
225         WRITE(numout,*) '      centered scheme                           ln_trcadv_cen = ', ln_trcadv_cen
226         WRITE(numout,*) '            horizontal 2nd/4th order               nn_cen_h   = ', nn_fct_h
227         WRITE(numout,*) '            vertical   2nd/4th order               nn_cen_v   = ', nn_fct_v
228         WRITE(numout,*) '      Flux Corrected Transport scheme           ln_trcadv_fct = ', ln_trcadv_fct
229         WRITE(numout,*) '            horizontal 2nd/4th order               nn_fct_h   = ', nn_fct_h
230         WRITE(numout,*) '            vertical   2nd/4th order               nn_fct_v   = ', nn_fct_v
231         WRITE(numout,*) '      MUSCL scheme                              ln_trcadv_mus = ', ln_trcadv_mus
232         WRITE(numout,*) '            + upstream scheme near river mouths    ln_mus_ups = ', ln_mus_ups
233         WRITE(numout,*) '      UBS scheme                                ln_trcadv_ubs = ', ln_trcadv_ubs
234         WRITE(numout,*) '            vertical   2nd/4th order               nn_ubs_v   = ', nn_ubs_v
235         WRITE(numout,*) '      QUICKEST scheme                           ln_trcadv_qck = ', ln_trcadv_qck
236      ENDIF
237      !
238      !                                !==  Parameter control & set nadv ==!
239      ioptio = 0
240      IF( ln_trcadv_OFF ) THEN   ;   ioptio = ioptio + 1   ;   nadv = np_NO_adv   ;   ENDIF
241      IF( ln_trcadv_cen ) THEN   ;   ioptio = ioptio + 1   ;   nadv = np_CEN      ;   ENDIF
242      IF( ln_trcadv_fct ) THEN   ;   ioptio = ioptio + 1   ;   nadv = np_FCT      ;   ENDIF
243      IF( ln_trcadv_mus ) THEN   ;   ioptio = ioptio + 1   ;   nadv = np_MUS      ;   ENDIF
244      IF( ln_trcadv_ubs ) THEN   ;   ioptio = ioptio + 1   ;   nadv = np_UBS      ;   ENDIF
245      IF( ln_trcadv_qck ) THEN   ;   ioptio = ioptio + 1   ;   nadv = np_QCK      ;   ENDIF
246      !
247      IF( ioptio /= 1 )   CALL ctl_stop( 'trc_adv_ini: Choose ONE advection option in namelist namtrc_adv' )
248      !
249      IF( ln_trcadv_cen .AND. ( nn_cen_h /= 2 .AND. nn_cen_h /= 4 )   &
250                        .AND. ( nn_cen_v /= 2 .AND. nn_cen_v /= 4 )   ) THEN
251        CALL ctl_stop( 'trc_adv_ini: CEN scheme, choose 2nd or 4th order' )
252      ENDIF
253      IF( ln_trcadv_fct .AND. ( nn_fct_h /= 2 .AND. nn_fct_h /= 4 )   &
254                        .AND. ( nn_fct_v /= 2 .AND. nn_fct_v /= 4 )   ) THEN
255        CALL ctl_stop( 'trc_adv_ini: FCT scheme, choose 2nd or 4th order' )
256      ENDIF
257      IF( ln_trcadv_ubs .AND. ( nn_ubs_v /= 2 .AND. nn_ubs_v /= 4 )   ) THEN
258        CALL ctl_stop( 'trc_adv_ini: UBS scheme, choose 2nd or 4th order' )
259      ENDIF
260      IF( ln_trcadv_ubs .AND. nn_ubs_v == 4 ) THEN
261         CALL ctl_warn( 'trc_adv_ini: UBS scheme, only 2nd FCT scheme available on the vertical. It will be used' )
262      ENDIF
263      IF( ln_isfcav ) THEN                                                       ! ice-shelf cavities
264         IF(  ln_trcadv_cen .AND. nn_cen_v == 4    .OR.   &                            ! NO 4th order with ISF
265            & ln_trcadv_fct .AND. nn_fct_v == 4   )   CALL ctl_stop( 'tra_adv_ini: 4th order COMPACT scheme not allowed with ISF' )
266      ENDIF
267      !
268      !                                !==  Print the choice  ==! 
269      IF(lwp) THEN
270         WRITE(numout,*)
271         SELECT CASE ( nadv )
272         CASE( np_NO_adv  )   ;   WRITE(numout,*) '      ===>>   NO passive tracer advection'
273         CASE( np_CEN     )   ;   WRITE(numout,*) '      ===>>   CEN      scheme is used. Horizontal order: ', nn_cen_h,   &
274            &                                                                     ' Vertical   order: ', nn_cen_v
275         CASE( np_FCT     )   ;   WRITE(numout,*) '      ===>>   FCT      scheme is used. Horizontal order: ', nn_fct_h,   &
276            &                                                                      ' Vertical   order: ', nn_fct_v
277         CASE( np_MUS     )   ;   WRITE(numout,*) '      ===>>   MUSCL    scheme is used'
278         CASE( np_UBS     )   ;   WRITE(numout,*) '      ===>>   UBS      scheme is used'
279         CASE( np_QCK     )   ;   WRITE(numout,*) '      ===>>   QUICKEST scheme is used'
280         END SELECT
281      ENDIF
282      !
283   END SUBROUTINE trc_adv_ini
284   
285#endif
286
287  !!======================================================================
288END MODULE trcadv
Note: See TracBrowser for help on using the repository browser.