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_fct.F90 in NEMO/branches/2021/dev_r14116_HPC-10_mcastril_Mixed_Precision_implementation/src/OCE/TRA – NEMO

source: NEMO/branches/2021/dev_r14116_HPC-10_mcastril_Mixed_Precision_implementation/src/OCE/TRA/traadv_fct.F90 @ 14986

Last change on this file since 14986 was 14986, checked in by sparonuz, 3 years ago

Merge trunk -r14984:HEAD

  • Property svn:keywords set to Id
File size: 55.4 KB
RevLine 
[5770]1MODULE traadv_fct
[3]2   !!==============================================================================
[5770]3   !!                       ***  MODULE  traadv_fct  ***
4   !! Ocean  tracers:  horizontal & vertical advective trend (2nd/4th order Flux Corrected Transport method)
[3]5   !!==============================================================================
[5770]6   !! History :  3.7  !  2015-09  (L. Debreu, G. Madec)  original code (inspired from traadv_tvd.F90)
[503]7   !!----------------------------------------------------------------------
[3]8
9   !!----------------------------------------------------------------------
[5770]10   !!  tra_adv_fct    : update the tracer trend with a 3D advective trends using a 2nd or 4th order FCT scheme
11   !!                   with sub-time-stepping in the vertical direction
[14072]12   !!  nonosc         : compute monotonic tracer fluxes by a non-oscillatory algorithm
[5770]13   !!  interp_4th_cpt : 4th order compact scheme for the vertical component of the advection
[3]14   !!----------------------------------------------------------------------
[3625]15   USE oce            ! ocean dynamics and active tracers
16   USE dom_oce        ! ocean space and time domain
[4990]17   USE trc_oce        ! share passive tracers/Ocean variables
18   USE trd_oce        ! trends: ocean variables
[3625]19   USE trdtra         ! tracers trends
[4990]20   USE diaptr         ! poleward transport diagnostics
[7646]21   USE diaar5         ! AR5 diagnostics
[12489]22   USE phycst  , ONLY : rho0_rcp
[11407]23   USE zdf_oce , ONLY : ln_zad_Aimp
[4990]24   !
[5770]25   USE in_out_manager ! I/O manager
[14072]26   USE iom            !
[3625]27   USE lib_mpp        ! MPP library
[14072]28   USE lbclnk         ! ocean lateral boundary condition (or mpp link)
29   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined)
[3]30
31   IMPLICIT NONE
32   PRIVATE
33
[9019]34   PUBLIC   tra_adv_fct        ! called by traadv.F90
35   PUBLIC   interp_4th_cpt     ! called by traadv_cen.F90
[3]36
[5770]37   LOGICAL  ::   l_trd   ! flag to compute trends
[7646]38   LOGICAL  ::   l_ptr   ! flag to compute poleward transport
39   LOGICAL  ::   l_hst   ! flag to compute heat/salt transport
[5770]40   REAL(wp) ::   r1_6 = 1._wp / 6._wp   ! =1/6
[2528]41
[7646]42   !                                        ! tridiag solver associated indices:
43   INTEGER, PARAMETER ::   np_NH   = 0   ! Neumann homogeneous boundary condition
44   INTEGER, PARAMETER ::   np_CEN2 = 1   ! 2nd order centered  boundary condition
45
[3]46   !! * Substitutions
[12377]47#  include "do_loop_substitute.h90"
[13237]48#  include "domzgr_substitute.h90"
[14219]49#  include "single_precision_substitute.h90"
[3]50   !!----------------------------------------------------------------------
[9598]51   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
[1152]52   !! $Id$
[10068]53   !! Software governed by the CeCILL license (see ./LICENSE)
[3]54   !!----------------------------------------------------------------------
55CONTAINS
56
[12377]57   SUBROUTINE tra_adv_fct( kt, kit000, cdtype, p2dt, pU, pV, pW,       &
58      &                    Kbb, Kmm, pt, kjpt, Krhs, kn_fct_h, kn_fct_v )
[3]59      !!----------------------------------------------------------------------
[5770]60      !!                  ***  ROUTINE tra_adv_fct  ***
[14072]61      !!
[6140]62      !! **  Purpose :   Compute the now trend due to total advection of tracers
63      !!               and add it to the general trend of tracer equations
[3]64      !!
[5770]65      !! **  Method  : - 2nd or 4th FCT scheme on the horizontal direction
66      !!               (choice through the value of kn_fct)
[14072]67      !!               - on the vertical the 4th order is a compact scheme
68      !!               - corrected flux (monotonic correction)
[3]69      !!
[12377]70      !! ** Action : - update pt(:,:,:,:,Krhs)  with the now advective tracer trends
[9019]71      !!             - send trends to trdtra module for further diagnostics (l_trdtra=T)
[12377]72      !!             - poleward advective heat and salt transport (ln_diaptr=T)
[503]73      !!----------------------------------------------------------------------
[12377]74      INTEGER                                  , INTENT(in   ) ::   kt              ! ocean time-step index
75      INTEGER                                  , INTENT(in   ) ::   Kbb, Kmm, Krhs  ! ocean time level indices
76      INTEGER                                  , INTENT(in   ) ::   kit000          ! first time step index
77      CHARACTER(len=3)                         , INTENT(in   ) ::   cdtype          ! =TRA or TRC (tracer indicator)
78      INTEGER                                  , INTENT(in   ) ::   kjpt            ! number of tracers
79      INTEGER                                  , INTENT(in   ) ::   kn_fct_h        ! order of the FCT scheme (=2 or 4)
80      INTEGER                                  , INTENT(in   ) ::   kn_fct_v        ! order of the FCT scheme (=2 or 4)
81      REAL(wp)                                 , INTENT(in   ) ::   p2dt            ! tracer time-step
[14986]82      ! TEMP: [tiling] This can be A2D(nn_hls) after all lbc_lnks removed in the nn_hls = 2 case
[12377]83      REAL(wp), DIMENSION(jpi,jpj,jpk         ), INTENT(in   ) ::   pU, pV, pW      ! 3 ocean volume flux components
[14219]84      REAL(dp), DIMENSION(jpi,jpj,jpk,kjpt,jpt), INTENT(inout) ::   pt              ! tracers and RHS of tracer equation
[2715]85      !
[13982]86      INTEGER  ::   ji, jj, jk, jn                           ! dummy loop indices
[6140]87      REAL(wp) ::   ztra                                     ! local scalar
[5770]88      REAL(wp) ::   zfp_ui, zfp_vj, zfp_wk, zC2t_u, zC4t_u   !   -      -
89      REAL(wp) ::   zfm_ui, zfm_vj, zfm_wk, zC2t_v, zC4t_v   !   -      -
[14644]90      REAL(dp), DIMENSION(A2D(nn_hls),jpk)        ::   zwx, zwy, zwz
91      REAL(wp), DIMENSION(A2D(nn_hls),jpk)        ::   zwi, ztu, ztv, zltu, zltv, ztw
[9019]92      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::   ztrdx, ztrdy, ztrdz, zptry
[11407]93      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::   zwinf, zwdia, zwsup
94      LOGICAL  ::   ll_zAimp                                 ! flag to apply adaptive implicit vertical advection
[3]95      !!----------------------------------------------------------------------
[3294]96      !
[14986]97#if defined key_loop_fusion
98      CALL tra_adv_fct_lf ( kt, nit000, cdtype, p2dt, pU, pV, pW, Kbb, Kmm, pt, kjpt, Krhs, kn_fct_h, kn_fct_v )
99#else
100      IF( .NOT. l_istiled .OR. ntile == 1 )  THEN                       ! Do only on the first tile
[13982]101         IF( kt == kit000 )  THEN
102            IF(lwp) WRITE(numout,*)
103            IF(lwp) WRITE(numout,*) 'tra_adv_fct : FCT advection scheme on ', cdtype
104            IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
105         ENDIF
106         !
107         l_trd = .FALSE.            ! set local switches
108         l_hst = .FALSE.
109         l_ptr = .FALSE.
110         ll_zAimp = .FALSE.
111         IF( ( cdtype == 'TRA' .AND. l_trdtra  ) .OR. ( cdtype =='TRC' .AND. l_trdtrc ) )      l_trd = .TRUE.
112         IF(   cdtype == 'TRA' .AND. ( iom_use( 'sophtadv' ) .OR. iom_use( 'sophtadv' ) ) )    l_ptr = .TRUE.
113         IF(   cdtype == 'TRA' .AND. ( iom_use("uadv_heattr") .OR. iom_use("vadv_heattr") .OR.  &
114            &                         iom_use("uadv_salttr") .OR. iom_use("vadv_salttr")  ) )  l_hst = .TRUE.
115         !
[3]116      ENDIF
[13982]117
[13226]118      !! -- init to 0
119      zwi(:,:,:) = 0._wp
120      zwx(:,:,:) = 0._wp
121      zwy(:,:,:) = 0._wp
122      zwz(:,:,:) = 0._wp
123      ztu(:,:,:) = 0._wp
124      ztv(:,:,:) = 0._wp
125      zltu(:,:,:) = 0._wp
126      zltv(:,:,:) = 0._wp
127      ztw(:,:,:) = 0._wp
[2528]128      !
[7646]129      IF( l_trd .OR. l_hst )  THEN
[13982]130         ALLOCATE( ztrdx(A2D(nn_hls),jpk), ztrdy(A2D(nn_hls),jpk), ztrdz(A2D(nn_hls),jpk) )
[7753]131         ztrdx(:,:,:) = 0._wp   ;    ztrdy(:,:,:) = 0._wp   ;   ztrdz(:,:,:) = 0._wp
[3294]132      ENDIF
[2528]133      !
[13982]134      IF( l_ptr ) THEN
135         ALLOCATE( zptry(A2D(nn_hls),jpk) )
[7753]136         zptry(:,:,:) = 0._wp
[7646]137      ENDIF
[2528]138      !
[11407]139      ! If adaptive vertical advection, check if it is needed on this PE at this time
140      IF( ln_zad_Aimp ) THEN
[14986]141         IF( MAXVAL( ABS( wi(A2D(1),:) ) ) > 0._wp ) ll_zAimp = .TRUE.
[11407]142      END IF
143      ! If active adaptive vertical advection, build tridiagonal matrix
144      IF( ll_zAimp ) THEN
[13982]145         ALLOCATE(zwdia(A2D(nn_hls),jpk), zwinf(A2D(nn_hls),jpk), zwsup(A2D(nn_hls),jpk))
146         DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 1, jpkm1 )
[13237]147            zwdia(ji,jj,jk) =  1._wp + p2dt * ( MAX( wi(ji,jj,jk) , 0._wp ) - MIN( wi(ji,jj,jk+1) , 0._wp ) )   &
148            &                               / e3t(ji,jj,jk,Krhs)
[12377]149            zwinf(ji,jj,jk) =  p2dt * MIN( wi(ji,jj,jk  ) , 0._wp ) / e3t(ji,jj,jk,Krhs)
150            zwsup(ji,jj,jk) = -p2dt * MAX( wi(ji,jj,jk+1) , 0._wp ) / e3t(ji,jj,jk,Krhs)
151         END_3D
[11407]152      END IF
153      !
[6140]154      DO jn = 1, kjpt            !==  loop over the tracers  ==!
[5770]155         !
156         !        !==  upstream advection with initial mass fluxes & intermediate update  ==!
[14072]157         !                    !* upstream tracer flux in the i and j direction
[13982]158         DO_3D( nn_hls, nn_hls-1, nn_hls, nn_hls-1, 1, jpkm1 )
[12377]159            ! upstream scheme
160            zfp_ui = pU(ji,jj,jk) + ABS( pU(ji,jj,jk) )
161            zfm_ui = pU(ji,jj,jk) - ABS( pU(ji,jj,jk) )
162            zfp_vj = pV(ji,jj,jk) + ABS( pV(ji,jj,jk) )
163            zfm_vj = pV(ji,jj,jk) - ABS( pV(ji,jj,jk) )
164            zwx(ji,jj,jk) = 0.5 * ( zfp_ui * pt(ji,jj,jk,jn,Kbb) + zfm_ui * pt(ji+1,jj  ,jk,jn,Kbb) )
165            zwy(ji,jj,jk) = 0.5 * ( zfp_vj * pt(ji,jj,jk,jn,Kbb) + zfm_vj * pt(ji  ,jj+1,jk,jn,Kbb) )
166         END_3D
[13497]167         !                               !* upstream tracer flux in the k direction *!
[14986]168         DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )      ! Interior value ( multiplied by wmask)
[12377]169            zfp_wk = pW(ji,jj,jk) + ABS( pW(ji,jj,jk) )
170            zfm_wk = pW(ji,jj,jk) - ABS( pW(ji,jj,jk) )
171            zwz(ji,jj,jk) = 0.5 * ( zfp_wk * pt(ji,jj,jk,jn,Kbb) + zfm_wk * pt(ji,jj,jk-1,jn,Kbb) ) * wmask(ji,jj,jk)
172         END_3D
[13497]173         IF( ln_linssh ) THEN               ! top ocean value (only in linear free surface as zwz has been w-masked)
174            IF( ln_isfcav ) THEN                        ! top of the ice-shelf cavities and at the ocean surface
[14986]175               DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )
[14072]176                  zwz(ji,jj, mikt(ji,jj) ) = pW(ji,jj,mikt(ji,jj)) * pt(ji,jj,mikt(ji,jj),jn,Kbb)   ! linear free surface
[12377]177               END_2D
[13497]178            ELSE                                        ! no cavities: only at the ocean surface
[14986]179               DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )
[13286]180                  zwz(ji,jj,1) = pW(ji,jj,1) * pt(ji,jj,1,jn,Kbb)
181               END_2D
[5770]182            ENDIF
[5120]183         ENDIF
[14072]184         !
[13982]185         DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 1, jpkm1 )   !* trend and after field with monotonic scheme
[13497]186            !                               ! total intermediate advective trends
[12377]187            ztra = - (  zwx(ji,jj,jk) - zwx(ji-1,jj  ,jk  )   &
188               &      + zwy(ji,jj,jk) - zwy(ji  ,jj-1,jk  )   &
189               &      + zwz(ji,jj,jk) - zwz(ji  ,jj  ,jk+1) ) * r1_e1e2t(ji,jj)
[13497]190            !                               ! update and guess with monotonic sheme
[13237]191            pt(ji,jj,jk,jn,Krhs) =                   pt(ji,jj,jk,jn,Krhs) +       ztra   &
192               &                                  / e3t(ji,jj,jk,Kmm ) * tmask(ji,jj,jk)
193            zwi(ji,jj,jk)    = ( e3t(ji,jj,jk,Kbb) * pt(ji,jj,jk,jn,Kbb) + p2dt * ztra ) &
194               &                                  / e3t(ji,jj,jk,Krhs) * tmask(ji,jj,jk)
[12377]195         END_3D
[14072]196
[11407]197         IF ( ll_zAimp ) THEN
198            CALL tridia_solver( zwdia, zwsup, zwinf, zwi, zwi , 0 )
199            !
[11411]200            ztw(:,:,1) = 0._wp ; ztw(:,:,jpk) = 0._wp ;
[13982]201            DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )       ! Interior value ( multiplied by wmask)
[12377]202               zfp_wk = wi(ji,jj,jk) + ABS( wi(ji,jj,jk) )
203               zfm_wk = wi(ji,jj,jk) - ABS( wi(ji,jj,jk) )
204               ztw(ji,jj,jk) =  0.5 * e1e2t(ji,jj) * ( zfp_wk * zwi(ji,jj,jk) + zfm_wk * zwi(ji,jj,jk-1) ) * wmask(ji,jj,jk)
205               zwz(ji,jj,jk) = zwz(ji,jj,jk) + ztw(ji,jj,jk) ! update vertical fluxes
206            END_3D
[13295]207            DO_3D( 0, 0, 0, 0, 1, jpkm1 )
[12377]208               pt(ji,jj,jk,jn,Krhs) = pt(ji,jj,jk,jn,Krhs) - ( ztw(ji,jj,jk) - ztw(ji  ,jj  ,jk+1) ) &
209                  &                                        * r1_e1e2t(ji,jj) / e3t(ji,jj,jk,Kmm)
210            END_3D
[11407]211            !
212         END IF
[13982]213         !
[7646]214         IF( l_trd .OR. l_hst )  THEN             ! trend diagnostics (contribution of upstream fluxes)
[7753]215            ztrdx(:,:,:) = zwx(:,:,:)   ;   ztrdy(:,:,:) = zwy(:,:,:)   ;   ztrdz(:,:,:) = zwz(:,:,:)
[2528]216         END IF
[5770]217         !                             ! "Poleward" heat and salt transports (contribution of upstream fluxes)
[14072]218         IF( l_ptr )   zptry(:,:,:) = zwy(:,:,:)
[5770]219         !
220         !        !==  anti-diffusive flux : high order minus low order  ==!
221         !
[6140]222         SELECT CASE( kn_fct_h )    !* horizontal anti-diffusive fluxes
[5770]223         !
[6140]224         CASE(  2  )                   !- 2nd order centered
[13982]225            DO_3D( nn_hls, nn_hls-1, nn_hls, nn_hls-1, 1, jpkm1 )
[12377]226               zwx(ji,jj,jk) = 0.5_wp * pU(ji,jj,jk) * ( pt(ji,jj,jk,jn,Kmm) + pt(ji+1,jj,jk,jn,Kmm) ) - zwx(ji,jj,jk)
227               zwy(ji,jj,jk) = 0.5_wp * pV(ji,jj,jk) * ( pt(ji,jj,jk,jn,Kmm) + pt(ji,jj+1,jk,jn,Kmm) ) - zwy(ji,jj,jk)
228            END_3D
[5770]229            !
[6140]230         CASE(  4  )                   !- 4th order centered
[7753]231            zltu(:,:,jpk) = 0._wp            ! Bottom value : flux set to zero
232            zltv(:,:,jpk) = 0._wp
[6140]233            DO jk = 1, jpkm1                 ! Laplacian
[13497]234               DO_2D( 1, 0, 1, 0 )                 ! 1st derivative (gradient)
[12377]235                  ztu(ji,jj,jk) = ( pt(ji+1,jj  ,jk,jn,Kmm) - pt(ji,jj,jk,jn,Kmm) ) * umask(ji,jj,jk)
236                  ztv(ji,jj,jk) = ( pt(ji  ,jj+1,jk,jn,Kmm) - pt(ji,jj,jk,jn,Kmm) ) * vmask(ji,jj,jk)
237               END_2D
[13497]238               DO_2D( 0, 0, 0, 0 )                 ! 2nd derivative * 1/ 6
[12377]239                  zltu(ji,jj,jk) = (  ztu(ji,jj,jk) + ztu(ji-1,jj,jk)  ) * r1_6
240                  zltv(ji,jj,jk) = (  ztv(ji,jj,jk) + ztv(ji,jj-1,jk)  ) * r1_6
241               END_2D
[503]242            END DO
[14986]243            ! NOTE [ comm_cleanup ] : need to change sign to ensure halo 1 - halo 2 compatibility
244            CALL lbc_lnk( 'traadv_fct', zltu, 'T', -1.0_wp , zltv, 'T', -1.0_wp, ld4only= .TRUE. )   ! Lateral boundary cond. (unchanged sgn)
[5770]245            !
[14986]246            DO_3D( nn_hls, nn_hls-1, nn_hls, nn_hls-1, 1, jpkm1 )
[12377]247               zC2t_u = pt(ji,jj,jk,jn,Kmm) + pt(ji+1,jj  ,jk,jn,Kmm)   ! 2 x C2 interpolation of T at u- & v-points
248               zC2t_v = pt(ji,jj,jk,jn,Kmm) + pt(ji  ,jj+1,jk,jn,Kmm)
[13982]249               !                                                        ! C4 minus upstream advective fluxes
[14986]250               ! round brackets added to fix the order of floating point operations
251               ! needed to ensure halo 1 - halo 2 compatibility
252               zwx(ji,jj,jk) =  0.5_wp * pU(ji,jj,jk) * ( zC2t_u + ( zltu(ji,jj,jk) - zltu(ji+1,jj,jk)   &
253                             &                                     )                                     & ! bracket for halo 1 - halo 2 compatibility
254                             &                          ) - zwx(ji,jj,jk)
255               zwy(ji,jj,jk) =  0.5_wp * pV(ji,jj,jk) * ( zC2t_v + ( zltv(ji,jj,jk) - zltv(ji,jj+1,jk)   &
256                             &                                     )                                     & ! bracket for halo 1 - halo 2 compatibility
257                             &                          ) - zwy(ji,jj,jk)
[12377]258            END_3D
[5770]259            !
[6140]260         CASE(  41 )                   !- 4th order centered       ==>>   !!gm coding attempt   need to be tested
[7753]261            ztu(:,:,jpk) = 0._wp             ! Bottom value : flux set to zero
262            ztv(:,:,jpk) = 0._wp
[14986]263            DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 1, jpkm1 )    ! 1st derivative (gradient)
[12377]264               ztu(ji,jj,jk) = ( pt(ji+1,jj  ,jk,jn,Kmm) - pt(ji,jj,jk,jn,Kmm) ) * umask(ji,jj,jk)
265               ztv(ji,jj,jk) = ( pt(ji  ,jj+1,jk,jn,Kmm) - pt(ji,jj,jk,jn,Kmm) ) * vmask(ji,jj,jk)
266            END_3D
[14986]267            IF (nn_hls==1) CALL lbc_lnk( 'traadv_fct', ztu, 'U', -1.0_wp , ztv, 'V', -1.0_wp, ld4only= .TRUE. )   ! Lateral boundary cond. (unchanged sgn)
[5770]268            !
[13497]269            DO_3D( 0, 0, 0, 0, 1, jpkm1 )    ! Horizontal advective fluxes
[12377]270               zC2t_u = pt(ji,jj,jk,jn,Kmm) + pt(ji+1,jj  ,jk,jn,Kmm)   ! 2 x C2 interpolation of T at u- & v-points (x2)
271               zC2t_v = pt(ji,jj,jk,jn,Kmm) + pt(ji  ,jj+1,jk,jn,Kmm)
272               !                                                  ! C4 interpolation of T at u- & v-points (x2)
273               zC4t_u =  zC2t_u + r1_6 * ( ztu(ji-1,jj  ,jk) - ztu(ji+1,jj  ,jk) )
274               zC4t_v =  zC2t_v + r1_6 * ( ztv(ji  ,jj-1,jk) - ztv(ji  ,jj+1,jk) )
[14072]275               !                                                  ! C4 minus upstream advective fluxes
[12377]276               zwx(ji,jj,jk) =  0.5_wp * pU(ji,jj,jk) * zC4t_u - zwx(ji,jj,jk)
277               zwy(ji,jj,jk) =  0.5_wp * pV(ji,jj,jk) * zC4t_v - zwy(ji,jj,jk)
278            END_3D
[14986]279            IF (nn_hls==2) CALL lbc_lnk( 'traadv_fct', zwx, 'U', -1.0_wp , zwy, 'V', -1.0_wp )   ! Lateral boundary cond. (unchanged sgn)
[5770]280            !
281         END SELECT
[14072]282         !
[6140]283         SELECT CASE( kn_fct_v )    !* vertical anti-diffusive fluxes (w-masked interior values)
[5770]284         !
[6140]285         CASE(  2  )                   !- 2nd order centered
[13982]286            DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )
[12377]287               zwz(ji,jj,jk) =  (  pW(ji,jj,jk) * 0.5_wp * ( pt(ji,jj,jk,jn,Kmm) + pt(ji,jj,jk-1,jn,Kmm) )   &
288                  &              - zwz(ji,jj,jk)  ) * wmask(ji,jj,jk)
289            END_3D
[5770]290            !
[6140]291         CASE(  4  )                   !- 4th order COMPACT
[14219]292            CALL interp_4th_cpt( CASTWP(pt(:,:,:,jn,Kmm)) , ztw )   ! zwt = COMPACT interpolation of T at w-point
[13982]293            DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )
[12377]294               zwz(ji,jj,jk) = ( pW(ji,jj,jk) * ztw(ji,jj,jk) - zwz(ji,jj,jk) ) * wmask(ji,jj,jk)
295            END_3D
[5770]296            !
297         END SELECT
[6140]298         IF( ln_linssh ) THEN    ! top ocean value: high order = upstream  ==>>  zwz=0
[7753]299            zwz(:,:,1) = 0._wp   ! only ocean surface as interior zwz values have been w-masked
[6140]300         ENDIF
[13982]301         !
[14986]302         IF (nn_hls==1) THEN
[14648]303            CALL lbc_lnk( 'traadv_fct', zwi, 'T', 1.0_wp )
304            CALL lbc_lnk( 'traadv_fct', zwx, 'U', -1.0_wp , zwy, 'V', -1.0_wp, zwz, 'T', 1.0_wp )
[13982]305         ELSE
306            CALL lbc_lnk( 'traadv_fct', zwi, 'T', 1.0_wp)
307         END IF
308         !
[11407]309         IF ( ll_zAimp ) THEN
[13982]310            DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 1, jpkm1 )    !* trend and after field with monotonic scheme
[13497]311               !                                                ! total intermediate advective trends
[12377]312               ztra = - (  zwx(ji,jj,jk) - zwx(ji-1,jj  ,jk  )   &
313                  &      + zwy(ji,jj,jk) - zwy(ji  ,jj-1,jk  )   &
314                  &      + zwz(ji,jj,jk) - zwz(ji  ,jj  ,jk+1) ) * r1_e1e2t(ji,jj)
[13982]315               ztw(ji,jj,jk) = zwi(ji,jj,jk) + p2dt * ztra / e3t(ji,jj,jk,Krhs) * tmask(ji,jj,jk)
[12377]316            END_3D
[11407]317            !
[11411]318            CALL tridia_solver( zwdia, zwsup, zwinf, ztw, ztw , 0 )
[11407]319            !
[13982]320            DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )       ! Interior value ( multiplied by wmask)
[12377]321               zfp_wk = wi(ji,jj,jk) + ABS( wi(ji,jj,jk) )
322               zfm_wk = wi(ji,jj,jk) - ABS( wi(ji,jj,jk) )
[13982]323               zwz(ji,jj,jk) = zwz(ji,jj,jk) + 0.5 * e1e2t(ji,jj) * ( zfp_wk * ztw(ji,jj,jk) + zfm_wk * ztw(ji,jj,jk-1) ) * wmask(ji,jj,jk)
[12377]324            END_3D
[11407]325         END IF
[11411]326         !
[5770]327         !        !==  monotonicity algorithm  ==!
328         !
[14219]329         CALL nonosc( Kmm, CASTWP(pt(:,:,:,jn,Kbb)), zwx, zwy, zwz, zwi, p2dt )
[6140]330         !
[5770]331         !        !==  final trend with corrected fluxes  ==!
332         !
[13295]333         DO_3D( 0, 0, 0, 0, 1, jpkm1 )
[12377]334            ztra = - (  zwx(ji,jj,jk) - zwx(ji-1,jj  ,jk  )   &
335               &      + zwy(ji,jj,jk) - zwy(ji  ,jj-1,jk  )   &
336               &      + zwz(ji,jj,jk) - zwz(ji  ,jj  ,jk+1) ) * r1_e1e2t(ji,jj)
337            pt(ji,jj,jk,jn,Krhs) = pt(ji,jj,jk,jn,Krhs) + ztra / e3t(ji,jj,jk,Kmm)
338            zwi(ji,jj,jk) = zwi(ji,jj,jk) + p2dt * ztra / e3t(ji,jj,jk,Krhs) * tmask(ji,jj,jk)
339         END_3D
[5770]340         !
[11407]341         IF ( ll_zAimp ) THEN
342            !
[11411]343            ztw(:,:,1) = 0._wp ; ztw(:,:,jpk) = 0._wp
[13497]344            DO_3D( 0, 0, 0, 0, 2, jpkm1 )      ! Interior value ( multiplied by wmask)
[12377]345               zfp_wk = wi(ji,jj,jk) + ABS( wi(ji,jj,jk) )
346               zfm_wk = wi(ji,jj,jk) - ABS( wi(ji,jj,jk) )
347               ztw(ji,jj,jk) = - 0.5 * e1e2t(ji,jj) * ( zfp_wk * zwi(ji,jj,jk) + zfm_wk * zwi(ji,jj,jk-1) ) * wmask(ji,jj,jk)
348               zwz(ji,jj,jk) = zwz(ji,jj,jk) + ztw(ji,jj,jk) ! Update vertical fluxes for trend diagnostic
349            END_3D
[13295]350            DO_3D( 0, 0, 0, 0, 1, jpkm1 )
[12377]351               pt(ji,jj,jk,jn,Krhs) = pt(ji,jj,jk,jn,Krhs) - ( ztw(ji,jj,jk) - ztw(ji  ,jj  ,jk+1) ) &
352                  &                                        * r1_e1e2t(ji,jj) / e3t(ji,jj,jk,Kmm)
353            END_3D
[13982]354         END IF
[9019]355         IF( l_trd .OR. l_hst ) THEN   ! trend diagnostics // heat/salt transport
[13982]356            ztrdx(:,:,:) = ztrdx(:,:,:) + zwx(:,:,:)  ! <<< add anti-diffusive fluxes
[9019]357            ztrdy(:,:,:) = ztrdy(:,:,:) + zwy(:,:,:)  !     to upstream fluxes
358            ztrdz(:,:,:) = ztrdz(:,:,:) + zwz(:,:,:)  !
[5770]359            !
[9019]360            IF( l_trd ) THEN              ! trend diagnostics
[14219]361               CALL trd_tra( kt, Kmm, Krhs, cdtype, jn, jptra_xad, ztrdx, pU, CASTWP(pt(:,:,:,jn,Kmm)) )
362               CALL trd_tra( kt, Kmm, Krhs, cdtype, jn, jptra_yad, ztrdy, pV, CASTWP(pt(:,:,:,jn,Kmm)) )
363               CALL trd_tra( kt, Kmm, Krhs, cdtype, jn, jptra_zad, ztrdz, pW, CASTWP(pt(:,:,:,jn,Kmm)) )
[9019]364            ENDIF
365            !                             ! heat/salt transport
366            IF( l_hst )   CALL dia_ar5_hst( jn, 'adv', ztrdx(:,:,:), ztrdy(:,:,:) )
[5770]367            !
[9019]368         ENDIF
369         IF( l_ptr ) THEN              ! "Poleward" transports
370            zptry(:,:,:) = zptry(:,:,:) + zwy(:,:,:)  ! <<< add anti-diffusive fluxes
[7646]371            CALL dia_ptr_hst( jn, 'adv', zptry(:,:,:) )
[2528]372         ENDIF
[503]373         !
[6140]374      END DO                     ! end of tracer loop
[503]375      !
[11407]376      IF ( ll_zAimp ) THEN
377         DEALLOCATE( zwdia, zwinf, zwsup )
378      ENDIF
[13982]379      IF( l_trd .OR. l_hst ) THEN
[10024]380         DEALLOCATE( ztrdx, ztrdy, ztrdz )
381      ENDIF
[14072]382      IF( l_ptr ) THEN
[10024]383         DEALLOCATE( zptry )
384      ENDIF
385      !
[14986]386#endif
[5770]387   END SUBROUTINE tra_adv_fct
[3]388
[5737]389
[12377]390   SUBROUTINE nonosc( Kmm, pbef, paa, pbb, pcc, paft, p2dt )
[3]391      !!---------------------------------------------------------------------
392      !!                    ***  ROUTINE nonosc  ***
393      !!
[14072]394      !! **  Purpose :   compute monotonic tracer fluxes from the upstream
395      !!       scheme and the before field by a nonoscillatory algorithm
396      !!
[3]397      !! **  Method  :   ... ???
398      !!       warning : pbef and paft must be masked, but the boundaries
399      !!       conditions on the fluxes are not necessary zalezak (1979)
400      !!       drange (1995) multi-dimensional forward-in-time and upstream-
401      !!       in-space based differencing for fluid
402      !!----------------------------------------------------------------------
[13982]403      INTEGER                         , INTENT(in   ) ::   Kmm             ! time level index
404      REAL(wp)                        , INTENT(in   ) ::   p2dt            ! tracer time-step
405      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in   ) ::   pbef            ! before field
406      REAL(wp), DIMENSION(A2D(nn_hls)    ,jpk), INTENT(in   ) ::   paft            ! after field
[14219]407      REAL(dp), DIMENSION(A2D(nn_hls)    ,jpk), INTENT(inout) ::   paa, pbb, pcc   ! monotonic fluxes in the 3 directions
[2715]408      !
[4990]409      INTEGER  ::   ji, jj, jk   ! dummy loop indices
410      INTEGER  ::   ikm1         ! local integer
[13226]411      REAL(dp) ::   zpos, zneg, zbt, za, zb, zc, zbig, zrtrn    ! local scalars
412      REAL(dp) ::   zau, zbu, zcu, zav, zbv, zcv, zup, zdo            !   -      -
[13982]413      REAL(dp), DIMENSION(A2D(nn_hls),jpk) :: zbetup, zbetdo, zbup, zbdo
[3]414      !!----------------------------------------------------------------------
[3294]415      !
[13226]416      zbig  = 1.e+40_dp
417      zrtrn = 1.e-15_dp
418      zbetup(:,:,:) = 0._dp   ;   zbetdo(:,:,:) = 0._dp
[785]419
[3]420      ! Search local extrema
421      ! --------------------
[785]422      ! max/min of pbef & paft with large negative/positive value (-/+zbig) inside land
[13982]423      DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, 1, jpk )
424         zbup(ji,jj,jk) = MAX( pbef(ji,jj,jk) * tmask(ji,jj,jk) - zbig * ( 1._wp - tmask(ji,jj,jk) ),   &
425            &                  paft(ji,jj,jk) * tmask(ji,jj,jk) - zbig * ( 1._wp - tmask(ji,jj,jk) )  )
426         zbdo(ji,jj,jk) = MIN( pbef(ji,jj,jk) * tmask(ji,jj,jk) + zbig * ( 1._wp - tmask(ji,jj,jk) ),   &
427            &                  paft(ji,jj,jk) * tmask(ji,jj,jk) + zbig * ( 1._wp - tmask(ji,jj,jk) )  )
428      END_3D
[785]429
[5120]430      DO jk = 1, jpkm1
431         ikm1 = MAX(jk-1,1)
[13982]432         DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )
[5120]433
[12377]434            ! search maximum in neighbourhood
435            zup = MAX(  zbup(ji  ,jj  ,jk  ),   &
436               &        zbup(ji-1,jj  ,jk  ), zbup(ji+1,jj  ,jk  ),   &
437               &        zbup(ji  ,jj-1,jk  ), zbup(ji  ,jj+1,jk  ),   &
438               &        zbup(ji  ,jj  ,ikm1), zbup(ji  ,jj  ,jk+1)  )
[3]439
[12377]440            ! search minimum in neighbourhood
441            zdo = MIN(  zbdo(ji  ,jj  ,jk  ),   &
442               &        zbdo(ji-1,jj  ,jk  ), zbdo(ji+1,jj  ,jk  ),   &
443               &        zbdo(ji  ,jj-1,jk  ), zbdo(ji  ,jj+1,jk  ),   &
444               &        zbdo(ji  ,jj  ,ikm1), zbdo(ji  ,jj  ,jk+1)  )
[3]445
[12377]446            ! positive part of the flux
447            zpos = MAX( 0., paa(ji-1,jj  ,jk  ) ) - MIN( 0., paa(ji  ,jj  ,jk  ) )   &
448               & + MAX( 0., pbb(ji  ,jj-1,jk  ) ) - MIN( 0., pbb(ji  ,jj  ,jk  ) )   &
449               & + MAX( 0., pcc(ji  ,jj  ,jk+1) ) - MIN( 0., pcc(ji  ,jj  ,jk  ) )
[785]450
[12377]451            ! negative part of the flux
452            zneg = MAX( 0., paa(ji  ,jj  ,jk  ) ) - MIN( 0., paa(ji-1,jj  ,jk  ) )   &
453               & + MAX( 0., pbb(ji  ,jj  ,jk  ) ) - MIN( 0., pbb(ji  ,jj-1,jk  ) )   &
454               & + MAX( 0., pcc(ji  ,jj  ,jk  ) ) - MIN( 0., pcc(ji  ,jj  ,jk+1) )
[785]455
[12377]456            ! up & down beta terms
457            zbt = e1e2t(ji,jj) * e3t(ji,jj,jk,Kmm) / p2dt
458            zbetup(ji,jj,jk) = ( zup            - paft(ji,jj,jk) ) / ( zpos + zrtrn ) * zbt
459            zbetdo(ji,jj,jk) = ( paft(ji,jj,jk) - zdo            ) / ( zneg + zrtrn ) * zbt
460         END_2D
[3]461      END DO
[14986]462      IF (nn_hls==1) CALL lbc_lnk( 'traadv_fct', zbetup, 'T', 1.0_wp , zbetdo, 'T', 1.0_wp, ld4only= .TRUE. )   ! lateral boundary cond. (unchanged sign)
[3]463
[237]464      ! 3. monotonic flux in the i & j direction (paa & pbb)
465      ! ----------------------------------------
[13982]466      DO_3D( 1, 0, 1, 0, 1, jpkm1 )
[12377]467         zau = MIN( 1._wp, zbetdo(ji,jj,jk), zbetup(ji+1,jj,jk) )
468         zbu = MIN( 1._wp, zbetup(ji,jj,jk), zbetdo(ji+1,jj,jk) )
[14495]469         zcu =       ( 0.5  + SIGN( 0.5_wp , CASTWP(paa(ji,jj,jk)) ) )
[12377]470         paa(ji,jj,jk) = paa(ji,jj,jk) * ( zcu * zau + ( 1._wp - zcu) * zbu )
[3]471
[12377]472         zav = MIN( 1._wp, zbetdo(ji,jj,jk), zbetup(ji,jj+1,jk) )
473         zbv = MIN( 1._wp, zbetup(ji,jj,jk), zbetdo(ji,jj+1,jk) )
[14495]474         zcv =       ( 0.5  + SIGN( 0.5_wp , CASTWP(pbb(ji,jj,jk)) ) )
[12377]475         pbb(ji,jj,jk) = pbb(ji,jj,jk) * ( zcv * zav + ( 1._wp - zcv) * zbv )
[3]476
[13497]477      ! monotonic flux in the k direction, i.e. pcc
478      ! -------------------------------------------
[12377]479         za = MIN( 1., zbetdo(ji,jj,jk+1), zbetup(ji,jj,jk) )
480         zb = MIN( 1., zbetup(ji,jj,jk+1), zbetdo(ji,jj,jk) )
[14495]481         zc =       ( 0.5  + SIGN( 0.5_wp , CASTWP(pcc(ji,jj,jk+1)) ) )
[12377]482         pcc(ji,jj,jk+1) = pcc(ji,jj,jk+1) * ( zc * za + ( 1._wp - zc) * zb )
483      END_3D
[503]484      !
[3]485   END SUBROUTINE nonosc
486
[5770]487
[7646]488   SUBROUTINE interp_4th_cpt_org( pt_in, pt_out )
[5770]489      !!----------------------------------------------------------------------
[7646]490      !!                  ***  ROUTINE interp_4th_cpt_org  ***
[14072]491      !!
[5770]492      !! **  Purpose :   Compute the interpolation of tracer at w-point
493      !!
494      !! **  Method  :   4th order compact interpolation
495      !!----------------------------------------------------------------------
496      REAL(wp),DIMENSION(jpi,jpj,jpk), INTENT(in   ) ::   pt_in    ! now tracer fields
497      REAL(wp),DIMENSION(jpi,jpj,jpk), INTENT(  out) ::   pt_out   ! now tracer field interpolated at w-pts
498      !
499      INTEGER :: ji, jj, jk   ! dummy loop integers
500      REAL(wp),DIMENSION(jpi,jpj,jpk) :: zwd, zwi, zws, zwrm, zwt
501      !!----------------------------------------------------------------------
[14072]502
[13497]503      DO_3D( 1, 1, 1, 1, 3, jpkm1 )       !==  build the three diagonal matrix  ==!
[12377]504         zwd (ji,jj,jk) = 4._wp
505         zwi (ji,jj,jk) = 1._wp
506         zws (ji,jj,jk) = 1._wp
507         zwrm(ji,jj,jk) = 3._wp * ( pt_in(ji,jj,jk-1) + pt_in(ji,jj,jk) )
508         !
509         IF( tmask(ji,jj,jk+1) == 0._wp) THEN   ! Switch to second order centered at bottom
[5770]510            zwd (ji,jj,jk) = 1._wp
511            zwi (ji,jj,jk) = 0._wp
512            zws (ji,jj,jk) = 0._wp
[14072]513            zwrm(ji,jj,jk) = 0.5 * ( pt_in(ji,jj,jk-1) + pt_in(ji,jj,jk) )
[12377]514         ENDIF
515      END_3D
[5770]516      !
[13497]517      jk = 2                                    ! Switch to second order centered at top
[13295]518      DO_2D( 1, 1, 1, 1 )
[12377]519         zwd (ji,jj,jk) = 1._wp
520         zwi (ji,jj,jk) = 0._wp
521         zws (ji,jj,jk) = 0._wp
522         zwrm(ji,jj,jk) = 0.5 * ( pt_in(ji,jj,jk-1) + pt_in(ji,jj,jk) )
523      END_2D
524      !
[5770]525      !                       !==  tridiagonal solve  ==!
[13497]526      DO_2D( 1, 1, 1, 1 )           ! first recurrence
[12377]527         zwt(ji,jj,2) = zwd(ji,jj,2)
528      END_2D
[13295]529      DO_3D( 1, 1, 1, 1, 3, jpkm1 )
[12377]530         zwt(ji,jj,jk) = zwd(ji,jj,jk) - zwi(ji,jj,jk) * zws(ji,jj,jk-1) /zwt(ji,jj,jk-1)
531      END_3D
[5770]532      !
[13497]533      DO_2D( 1, 1, 1, 1 )           ! second recurrence:    Zk = Yk - Ik / Tk-1  Zk-1
[12377]534         pt_out(ji,jj,2) = zwrm(ji,jj,2)
535      END_2D
[13295]536      DO_3D( 1, 1, 1, 1, 3, jpkm1 )
[14072]537         pt_out(ji,jj,jk) = zwrm(ji,jj,jk) - zwi(ji,jj,jk) / zwt(ji,jj,jk-1) *pt_out(ji,jj,jk-1)
[12377]538      END_3D
[5770]539
[13497]540      DO_2D( 1, 1, 1, 1 )           ! third recurrence: Xk = (Zk - Sk Xk+1 ) / Tk
[12377]541         pt_out(ji,jj,jpkm1) = pt_out(ji,jj,jpkm1) / zwt(ji,jj,jpkm1)
542      END_2D
[13295]543      DO_3DS( 1, 1, 1, 1, jpk-2, 2, -1 )
[12377]544         pt_out(ji,jj,jk) = ( pt_out(ji,jj,jk) - zws(ji,jj,jk) * pt_out(ji,jj,jk+1) ) / zwt(ji,jj,jk)
545      END_3D
[14072]546      !
[7646]547   END SUBROUTINE interp_4th_cpt_org
548
[14072]549
[7646]550   SUBROUTINE interp_4th_cpt( pt_in, pt_out )
551      !!----------------------------------------------------------------------
552      !!                  ***  ROUTINE interp_4th_cpt  ***
[14072]553      !!
[7646]554      !! **  Purpose :   Compute the interpolation of tracer at w-point
555      !!
556      !! **  Method  :   4th order compact interpolation
557      !!----------------------------------------------------------------------
558      REAL(wp),DIMENSION(jpi,jpj,jpk), INTENT(in   ) ::   pt_in    ! field at t-point
[13982]559      REAL(wp),DIMENSION(A2D(nn_hls)    ,jpk), INTENT(  out) ::   pt_out   ! field interpolated at w-point
[7646]560      !
561      INTEGER ::   ji, jj, jk   ! dummy loop integers
562      INTEGER ::   ikt, ikb     ! local integers
[13982]563      REAL(wp),DIMENSION(A2D(nn_hls),jpk) :: zwd, zwi, zws, zwrm, zwt
[7646]564      !!----------------------------------------------------------------------
565      !
566      !                      !==  build the three diagonal matrix & the RHS  ==!
567      !
[13982]568      DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 3, jpkm1 )    ! interior (from jk=3 to jpk-1)
[12377]569         zwd (ji,jj,jk) = 3._wp * wmask(ji,jj,jk) + 1._wp                 !       diagonal
570         zwi (ji,jj,jk) =         wmask(ji,jj,jk)                         ! lower diagonal
571         zws (ji,jj,jk) =         wmask(ji,jj,jk)                         ! upper diagonal
572         zwrm(ji,jj,jk) = 3._wp * wmask(ji,jj,jk)                     &   ! RHS
573            &           *       ( pt_in(ji,jj,jk) + pt_in(ji,jj,jk-1) )
574      END_3D
[7646]575      !
576!!gm
577!      SELECT CASE( kbc )               !* boundary condition
578!      CASE( np_NH   )   ! Neumann homogeneous at top & bottom
579!      CASE( np_CEN2 )   ! 2nd order centered  at top & bottom
580!      END SELECT
[14072]581!!gm
[7646]582      !
[9901]583      IF ( ln_isfcav ) THEN            ! set level two values which may not be set in ISF case
584         zwd(:,:,2) = 1._wp  ;  zwi(:,:,2) = 0._wp  ;  zws(:,:,2) = 0._wp  ;  zwrm(:,:,2) = 0._wp
585      END IF
586      !
[13982]587      DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )              ! 2nd order centered at top & bottom
[12377]588         ikt = mikt(ji,jj) + 1            ! w-point below the 1st  wet point
589         ikb = MAX(mbkt(ji,jj), 2)        !     -   above the last wet point
590         !
591         zwd (ji,jj,ikt) = 1._wp          ! top
592         zwi (ji,jj,ikt) = 0._wp
593         zws (ji,jj,ikt) = 0._wp
594         zwrm(ji,jj,ikt) = 0.5_wp * ( pt_in(ji,jj,ikt-1) + pt_in(ji,jj,ikt) )
595         !
596         zwd (ji,jj,ikb) = 1._wp          ! bottom
597         zwi (ji,jj,ikb) = 0._wp
598         zws (ji,jj,ikb) = 0._wp
[14072]599         zwrm(ji,jj,ikb) = 0.5_wp * ( pt_in(ji,jj,ikb-1) + pt_in(ji,jj,ikb) )
[12377]600      END_2D
[7646]601      !
602      !                       !==  tridiagonal solver  ==!
603      !
[13982]604      DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )           !* 1st recurrence:   Tk = Dk - Ik Sk-1 / Tk-1
[12377]605         zwt(ji,jj,2) = zwd(ji,jj,2)
606      END_2D
[13982]607      DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 3, jpkm1 )
[12377]608         zwt(ji,jj,jk) = zwd(ji,jj,jk) - zwi(ji,jj,jk) * zws(ji,jj,jk-1) /zwt(ji,jj,jk-1)
609      END_3D
[7646]610      !
[13982]611      DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )           !* 2nd recurrence:    Zk = Yk - Ik / Tk-1  Zk-1
[12377]612         pt_out(ji,jj,2) = zwrm(ji,jj,2)
613      END_2D
[13982]614      DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 3, jpkm1 )
[14072]615         pt_out(ji,jj,jk) = zwrm(ji,jj,jk) - zwi(ji,jj,jk) / zwt(ji,jj,jk-1) *pt_out(ji,jj,jk-1)
[12377]616      END_3D
[7646]617
[13982]618      DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )           !* 3d recurrence:    Xk = (Zk - Sk Xk+1 ) / Tk
[12377]619         pt_out(ji,jj,jpkm1) = pt_out(ji,jj,jpkm1) / zwt(ji,jj,jpkm1)
620      END_2D
[13982]621      DO_3DS( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, jpk-2, 2, -1 )
[12377]622         pt_out(ji,jj,jk) = ( pt_out(ji,jj,jk) - zws(ji,jj,jk) * pt_out(ji,jj,jk+1) ) / zwt(ji,jj,jk)
623      END_3D
[14072]624      !
[5770]625   END SUBROUTINE interp_4th_cpt
[7646]626
627
628   SUBROUTINE tridia_solver( pD, pU, pL, pRHS, pt_out , klev )
629      !!----------------------------------------------------------------------
630      !!                  ***  ROUTINE tridia_solver  ***
[14072]631      !!
[7646]632      !! **  Purpose :   solve a symmetric 3diagonal system
633      !!
634      !! **  Method  :   solve M.t_out = RHS(t)  where M is a tri diagonal matrix ( jpk*jpk )
[14072]635      !!
[7646]636      !!             ( D_1 U_1  0   0   0  )( t_1 )   ( RHS_1 )
637      !!             ( L_2 D_2 U_2  0   0  )( t_2 )   ( RHS_2 )
638      !!             (  0  L_3 D_3 U_3  0  )( t_3 ) = ( RHS_3 )
639      !!             (        ...          )( ... )   ( ...  )
640      !!             (  0   0   0  L_k D_k )( t_k )   ( RHS_k )
[14072]641      !!
[7646]642      !!        M is decomposed in the product of an upper and lower triangular matrix.
[14072]643      !!        The tri-diagonals matrix is given as input 3D arrays:   pD, pU, pL
[7646]644      !!        (i.e. the Diagonal, the Upper diagonal, and the Lower diagonal).
645      !!        The solution is pta.
646      !!        The 3d array zwt is used as a work space array.
647      !!----------------------------------------------------------------------
[14649]648      REAL(wp),DIMENSION(A2D(nn_hls),jpk), INTENT(in   ) ::   pD, pU, pL    ! 3-diagonal matrix
[13982]649      REAL(wp),DIMENSION(A2D(nn_hls),jpk), INTENT(in   ) ::   pRHS          ! Right-Hand-Side
650      REAL(wp),DIMENSION(A2D(nn_hls),jpk), INTENT(  out) ::   pt_out        !!gm field at level=F(klev)
651      INTEGER                    , INTENT(in   ) ::   klev          ! =1 pt_out at w-level
652      !                                                             ! =0 pt at t-level
[7646]653      INTEGER ::   ji, jj, jk   ! dummy loop integers
654      INTEGER ::   kstart       ! local indices
[13982]655      REAL(wp),DIMENSION(A2D(nn_hls),jpk) ::   zwt   ! 3D work array
[7646]656      !!----------------------------------------------------------------------
657      !
658      kstart =  1  + klev
659      !
[13982]660      DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )                         !* 1st recurrence:   Tk = Dk - Ik Sk-1 / Tk-1
[12377]661         zwt(ji,jj,kstart) = pD(ji,jj,kstart)
662      END_2D
[13982]663      DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, kstart+1, jpkm1 )
[12377]664         zwt(ji,jj,jk) = pD(ji,jj,jk) - pL(ji,jj,jk) * pU(ji,jj,jk-1) /zwt(ji,jj,jk-1)
665      END_3D
[7646]666      !
[13982]667      DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )                        !* 2nd recurrence:    Zk = Yk - Ik / Tk-1  Zk-1
[12377]668         pt_out(ji,jj,kstart) = pRHS(ji,jj,kstart)
669      END_2D
[13982]670      DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, kstart+1, jpkm1 )
[14072]671         pt_out(ji,jj,jk) = pRHS(ji,jj,jk) - pL(ji,jj,jk) / zwt(ji,jj,jk-1) *pt_out(ji,jj,jk-1)
[12377]672      END_3D
[7646]673
[13982]674      DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )                       !* 3d recurrence:    Xk = (Zk - Sk Xk+1 ) / Tk
[12377]675         pt_out(ji,jj,jpkm1) = pt_out(ji,jj,jpkm1) / zwt(ji,jj,jpkm1)
676      END_2D
[13982]677      DO_3DS( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, jpk-2, kstart, -1 )
[12377]678         pt_out(ji,jj,jk) = ( pt_out(ji,jj,jk) - pU(ji,jj,jk) * pt_out(ji,jj,jk+1) ) / zwt(ji,jj,jk)
679      END_3D
[7646]680      !
681   END SUBROUTINE tridia_solver
682
[14986]683#if defined key_loop_fusion
684#define tracer_flux_i(out,zfp,zfm,ji,jj,jk) \
685        zfp = pU(ji,jj,jk) + ABS( pU(ji,jj,jk) ) ; \
686        zfm = pU(ji,jj,jk) - ABS( pU(ji,jj,jk) ) ; \
687        out = 0.5 * ( zfp * pt(ji,jj,jk,jn,Kbb) + zfm * pt(ji+1,jj,jk,jn,Kbb) )
688
689#define tracer_flux_j(out,zfp,zfm,ji,jj,jk) \
690        zfp = pV(ji,jj,jk) + ABS( pV(ji,jj,jk) ) ; \
691        zfm = pV(ji,jj,jk) - ABS( pV(ji,jj,jk) ) ; \
692        out = 0.5 * ( zfp * pt(ji,jj,jk,jn,Kbb) + zfm * pt(ji,jj+1,jk,jn,Kbb) )
693
694   SUBROUTINE tra_adv_fct_lf( kt, kit000, cdtype, p2dt, pU, pV, pW,       &
695      &                    Kbb, Kmm, pt, kjpt, Krhs, kn_fct_h, kn_fct_v )
696      !!----------------------------------------------------------------------
697      !!                  ***  ROUTINE tra_adv_fct  ***
698      !!
699      !! **  Purpose :   Compute the now trend due to total advection of tracers
700      !!               and add it to the general trend of tracer equations
701      !!
702      !! **  Method  : - 2nd or 4th FCT scheme on the horizontal direction
703      !!               (choice through the value of kn_fct)
704      !!               - on the vertical the 4th order is a compact scheme
705      !!               - corrected flux (monotonic correction)
706      !!
707      !! ** Action : - update pt(:,:,:,:,Krhs)  with the now advective tracer trends
708      !!             - send trends to trdtra module for further diagnostics (l_trdtra=T)
709      !!             - poleward advective heat and salt transport (ln_diaptr=T)
710      !!----------------------------------------------------------------------
711      INTEGER                                  , INTENT(in   ) ::   kt              ! ocean time-step index
712      INTEGER                                  , INTENT(in   ) ::   Kbb, Kmm, Krhs  ! ocean time level indices
713      INTEGER                                  , INTENT(in   ) ::   kit000          ! first time step index
714      CHARACTER(len=3)                         , INTENT(in   ) ::   cdtype          ! =TRA or TRC (tracer indicator)
715      INTEGER                                  , INTENT(in   ) ::   kjpt            ! number of tracers
716      INTEGER                                  , INTENT(in   ) ::   kn_fct_h        ! order of the FCT scheme (=2 or 4)
717      INTEGER                                  , INTENT(in   ) ::   kn_fct_v        ! order of the FCT scheme (=2 or 4)
718      REAL(wp)                                 , INTENT(in   ) ::   p2dt            ! tracer time-step
719      REAL(wp), DIMENSION(jpi,jpj,jpk         ), INTENT(in   ) ::   pU, pV, pW      ! 3 ocean volume flux components
720      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt,jpt), INTENT(inout) ::   pt              ! tracers and RHS of tracer equation
721      !
722      INTEGER  ::   ji, jj, jk, jn                           ! dummy loop indices
723      REAL(wp) ::   ztra                                     ! local scalar
724      REAL(wp) ::   zwx_im1, zfp_ui, zfp_ui_m1, zfp_vj, zfp_vj_m1, zfp_wk, zC2t_u, zC4t_u   !   -      -
725      REAL(wp) ::   zwy_jm1, zfm_ui, zfm_ui_m1, zfm_vj, zfm_vj_m1, zfm_wk, zC2t_v, zC4t_v   !   -      -
726      REAL(wp) ::   ztu, ztv, ztu_im1, ztu_ip1, ztv_jm1, ztv_jp1
727      REAL(wp), DIMENSION(jpi,jpj,jpk)        ::   zwi, zwx_3d, zwy_3d, zwz, ztw, zltu_3d, zltv_3d
728      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::   ztrdx, ztrdy, ztrdz, zptry
729      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::   zwinf, zwdia, zwsup
730      LOGICAL  ::   ll_zAimp                                 ! flag to apply adaptive implicit vertical advection
731      !!----------------------------------------------------------------------
732      !
733      IF( kt == kit000 )  THEN
734         IF(lwp) WRITE(numout,*)
735         IF(lwp) WRITE(numout,*) 'tra_adv_fct_lf : FCT advection scheme on ', cdtype
736         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
737      ENDIF
738      !! -- init to 0
739      zwx_3d(:,:,:) = 0._wp
740      zwy_3d(:,:,:) = 0._wp
741      zwz(:,:,:) = 0._wp
742      zwi(:,:,:) = 0._wp
743      !
744      l_trd = .FALSE.            ! set local switches
745      l_hst = .FALSE.
746      l_ptr = .FALSE.
747      ll_zAimp = .FALSE.
748      IF( ( cdtype == 'TRA' .AND. l_trdtra  ) .OR. ( cdtype =='TRC' .AND. l_trdtrc ) )      l_trd = .TRUE.
749      IF(   cdtype == 'TRA' .AND. ( iom_use( 'sophtadv' ) .OR. iom_use( 'sophtadv' ) ) )    l_ptr = .TRUE.
750      IF(   cdtype == 'TRA' .AND. ( iom_use("uadv_heattr") .OR. iom_use("vadv_heattr") .OR.  &
751         &                         iom_use("uadv_salttr") .OR. iom_use("vadv_salttr")  ) )  l_hst = .TRUE.
752      !
753      IF( l_trd .OR. l_hst )  THEN
754         ALLOCATE( ztrdx(jpi,jpj,jpk), ztrdy(jpi,jpj,jpk), ztrdz(jpi,jpj,jpk) )
755         ztrdx(:,:,:) = 0._wp   ;    ztrdy(:,:,:) = 0._wp   ;   ztrdz(:,:,:) = 0._wp
756      ENDIF
757      !
758      IF( l_ptr ) THEN
759         ALLOCATE( zptry(jpi,jpj,jpk) )
760         zptry(:,:,:) = 0._wp
761      ENDIF
762      !
763      ! If adaptive vertical advection, check if it is needed on this PE at this time
764      IF( ln_zad_Aimp ) THEN
765         IF( MAXVAL( ABS( wi(:,:,:) ) ) > 0._wp ) ll_zAimp = .TRUE.
766      END IF
767      ! If active adaptive vertical advection, build tridiagonal matrix
768      IF( ll_zAimp ) THEN
769         ALLOCATE(zwdia(jpi,jpj,jpk), zwinf(jpi,jpj,jpk),zwsup(jpi,jpj,jpk))
770         DO_3D( 1, 1, 1, 1, 1, jpkm1 )
771            zwdia(ji,jj,jk) =  1._wp + p2dt * ( MAX( wi(ji,jj,jk) , 0._wp ) - MIN( wi(ji,jj,jk+1) , 0._wp ) )   &
772            &                               / e3t(ji,jj,jk,Krhs)
773            zwinf(ji,jj,jk) =  p2dt * MIN( wi(ji,jj,jk  ) , 0._wp ) / e3t(ji,jj,jk,Krhs)
774            zwsup(ji,jj,jk) = -p2dt * MAX( wi(ji,jj,jk+1) , 0._wp ) / e3t(ji,jj,jk,Krhs)
775         END_3D
776      END IF
777      !
778      DO jn = 1, kjpt            !==  loop over the tracers  ==!
779         !
780         !        !==  upstream advection with initial mass fluxes & intermediate update  ==!
781         !                               !* upstream tracer flux in the k direction *!
782         DO_3D( 1, 1, 1, 1, 2, jpkm1 )      ! Interior value ( multiplied by wmask)
783            zfp_wk = pW(ji,jj,jk) + ABS( pW(ji,jj,jk) )
784            zfm_wk = pW(ji,jj,jk) - ABS( pW(ji,jj,jk) )
785            zwz(ji,jj,jk) = 0.5 * ( zfp_wk * pt(ji,jj,jk,jn,Kbb) + zfm_wk * pt(ji,jj,jk-1,jn,Kbb) ) * wmask(ji,jj,jk)
786         END_3D
787         IF( ln_linssh ) THEN               ! top ocean value (only in linear free surface as zwz has been w-masked)
788            IF( ln_isfcav ) THEN                        ! top of the ice-shelf cavities and at the ocean surface
789               DO_2D( 1, 1, 1, 1 )
790                  zwz(ji,jj, mikt(ji,jj) ) = pW(ji,jj,mikt(ji,jj)) * pt(ji,jj,mikt(ji,jj),jn,Kbb)   ! linear free surface
791               END_2D
792            ELSE                                        ! no cavities: only at the ocean surface
793               DO_2D( 1, 1, 1, 1 )
794                  zwz(ji,jj,1) = pW(ji,jj,1) * pt(ji,jj,1,jn,Kbb)
795               END_2D
796            ENDIF
797         ENDIF
798         !
799         !                    !* upstream tracer flux in the i and j direction
800         DO jk = 1, jpkm1
801            DO jj = 1, jpj-1
802               tracer_flux_i(zwx_3d(1,jj,jk),zfp_ui,zfm_ui,1,jj,jk)
803               tracer_flux_j(zwy_3d(1,jj,jk),zfp_vj,zfm_vj,1,jj,jk)
804            END DO
805            DO ji = 1, jpi-1
806               tracer_flux_i(zwx_3d(ji,1,jk),zfp_ui,zfm_ui,ji,1,jk)
807               tracer_flux_j(zwy_3d(ji,1,jk),zfp_vj,zfm_vj,ji,1,jk)
808            END DO
809            DO_2D( 1, 1, 1, 1 )
810               tracer_flux_i(zwx_3d(ji,jj,jk),zfp_ui,zfm_ui,ji,jj,jk)
811               tracer_flux_i(zwx_im1,zfp_ui_m1,zfm_ui_m1,ji-1,jj,jk)
812               tracer_flux_j(zwy_3d(ji,jj,jk),zfp_vj,zfm_vj,ji,jj,jk)
813               tracer_flux_j(zwy_jm1,zfp_vj_m1,zfm_vj_m1,ji,jj-1,jk)
814               ztra = - ( zwx_3d(ji,jj,jk) - zwx_im1 + zwy_3d(ji,jj,jk) - zwy_jm1 + zwz(ji,jj,jk) - zwz(ji,jj,jk+1) ) * r1_e1e2t(ji,jj)
815               !                               ! update and guess with monotonic sheme
816               pt(ji,jj,jk,jn,Krhs) =                   pt(ji,jj,jk,jn,Krhs) +       ztra   &
817                  &                                  / e3t(ji,jj,jk,Kmm ) * tmask(ji,jj,jk)
818               zwi(ji,jj,jk)    = ( e3t(ji,jj,jk,Kbb) * pt(ji,jj,jk,jn,Kbb) + p2dt * ztra ) &
819                  &                                  / e3t(ji,jj,jk,Krhs) * tmask(ji,jj,jk)
820            END_2D
821         END DO
822
823         IF ( ll_zAimp ) THEN
824            CALL tridia_solver( zwdia, zwsup, zwinf, zwi, zwi , 0 )
825            !
826            ztw(:,:,1) = 0._wp ; ztw(:,:,jpk) = 0._wp ;
827            DO_3D( 1, 1, 1, 1, 2, jpkm1 )       ! Interior value ( multiplied by wmask)
828               zfp_wk = wi(ji,jj,jk) + ABS( wi(ji,jj,jk) )
829               zfm_wk = wi(ji,jj,jk) - ABS( wi(ji,jj,jk) )
830               ztw(ji,jj,jk) =  0.5 * e1e2t(ji,jj) * ( zfp_wk * zwi(ji,jj,jk) + zfm_wk * zwi(ji,jj,jk-1) ) * wmask(ji,jj,jk)
831               zwz(ji,jj,jk) = zwz(ji,jj,jk) + ztw(ji,jj,jk) ! update vertical fluxes
832            END_3D
833            DO_3D( 0, 0, 0, 0, 1, jpkm1 )
834               pt(ji,jj,jk,jn,Krhs) = pt(ji,jj,jk,jn,Krhs) - ( ztw(ji,jj,jk) - ztw(ji  ,jj  ,jk+1) ) &
835                  &                                        * r1_e1e2t(ji,jj) / e3t(ji,jj,jk,Kmm)
836            END_3D
837            !
838         END IF
839         !
840         IF( l_trd .OR. l_hst )  THEN             ! trend diagnostics (contribution of upstream fluxes)
841            ztrdx(:,:,:) = zwx_3d(:,:,:)   ;   ztrdy(:,:,:) = zwy_3d(:,:,:)   ;   ztrdz(:,:,:) = zwz(:,:,:)
842         END IF
843         !                             ! "Poleward" heat and salt transports (contribution of upstream fluxes)
844         IF( l_ptr )   zptry(:,:,:) = zwy_3d(:,:,:)
845         !
846         !        !==  anti-diffusive flux : high order minus low order  ==!
847         !
848         SELECT CASE( kn_fct_h )    !* horizontal anti-diffusive fluxes
849         !
850         CASE(  2  )                   !- 2nd order centered
851            DO_3D( 2, 1, 2, 1, 1, jpkm1 )
852               zwx_3d(ji,jj,jk) = 0.5_wp * pU(ji,jj,jk) * ( pt(ji,jj,jk,jn,Kmm) + pt(ji+1,jj,jk,jn,Kmm) ) - zwx_3d(ji,jj,jk)
853               zwy_3d(ji,jj,jk) = 0.5_wp * pV(ji,jj,jk) * ( pt(ji,jj,jk,jn,Kmm) + pt(ji,jj+1,jk,jn,Kmm) ) - zwy_3d(ji,jj,jk)
854            END_3D
855            !
856         CASE(  4  )                   !- 4th order centered
857            zltu_3d(:,:,jpk) = 0._wp            ! Bottom value : flux set to zero
858            zltv_3d(:,:,jpk) = 0._wp
859            !                          ! Laplacian
860            DO_3D( 0, 0, 0, 0, 1, jpkm1 )                 ! 2nd derivative * 1/ 6
861                  !             ! 1st derivative (gradient)
862                  ztu = ( pt(ji+1,jj,jk,jn,Kmm) - pt(ji,jj,jk,jn,Kmm) ) * umask(ji,jj,jk)
863                  ztu_im1 = ( pt(ji,jj,jk,jn,Kmm) - pt(ji-1,jj,jk,jn,Kmm) ) * umask(ji-1,jj,jk)
864                  ztv = ( pt(ji,jj+1,jk,jn,Kmm) - pt(ji,jj,jk,jn,Kmm) ) * vmask(ji,jj,jk)
865                  ztv_jm1 = ( pt(ji,jj,jk,jn,Kmm) - pt(ji,jj-1,jk,jn,Kmm) ) * vmask(ji,jj-1,jk)
866                  !             ! 2nd derivative * 1/ 6
867                  zltu_3d(ji,jj,jk) = (  ztu + ztu_im1  ) * r1_6
868                  zltv_3d(ji,jj,jk) = (  ztv + ztv_jm1  ) * r1_6
869               END_2D
870            END DO
871            ! NOTE [ comm_cleanup ] : need to change sign to ensure halo 1 - halo 2 compatibility
872            CALL lbc_lnk( 'traadv_fct', zltu_3d, 'T', -1.0_wp , zltv_3d, 'T', -1.0_wp )   ! Lateral boundary cond. (unchanged sgn)
873            !
874            DO_3D( 2, 1, 2, 1, 1, jpkm1 )
875               zC2t_u = pt(ji,jj,jk,jn,Kmm) + pt(ji+1,jj  ,jk,jn,Kmm)   ! 2 x C2 interpolation of T at u- & v-points
876               zC2t_v = pt(ji,jj,jk,jn,Kmm) + pt(ji  ,jj+1,jk,jn,Kmm)
877               !                                                        ! C4 minus upstream advective fluxes
878               ! round brackets added to fix the order of floating point operations
879               ! needed to ensure halo 1 - halo 2 compatibility
880               zwx_3d(ji,jj,jk) =  0.5_wp * pU(ji,jj,jk) * ( zC2t_u + ( zltu_3d(ji,jj,jk) - zltu_3d(ji+1,jj,jk)   &
881                             &                                        )                                           & ! bracket for halo 1 - halo 2 compatibility
882                             &                             ) - zwx_3d(ji,jj,jk)
883               zwy_3d(ji,jj,jk) =  0.5_wp * pV(ji,jj,jk) * ( zC2t_v + ( zltv_3d(ji,jj,jk) - zltv_3d(ji,jj+1,jk)   &
884                             &                                        )                                           & ! bracket for halo 1 - halo 2 compatibility
885                             &                             ) - zwy_3d(ji,jj,jk)
886            END_3D
887            !
888         CASE(  41 )                   !- 4th order centered       ==>>   !!gm coding attempt   need to be tested
889            DO_3D( 0, 0, 0, 0, 1, jpkm1 )    ! Horizontal advective fluxes
890               ztu_im1 = ( pt(ji  ,jj  ,jk,jn,Kmm) - pt(ji-1,jj,jk,jn,Kmm) ) * umask(ji-1,jj,jk)
891               ztu_ip1 = ( pt(ji+2,jj  ,jk,jn,Kmm) - pt(ji+1,jj,jk,jn,Kmm) ) * umask(ji+1,jj,jk)
892
893               ztv_jm1 = ( pt(ji,jj  ,jk,jn,Kmm) - pt(ji,jj-1,jk,jn,Kmm) ) * vmask(ji,jj-1,jk)
894               ztv_jp1 = ( pt(ji,jj+2,jk,jn,Kmm) - pt(ji,jj+1,jk,jn,Kmm) ) * vmask(ji,jj+1,jk)
895               zC2t_u = pt(ji,jj,jk,jn,Kmm) + pt(ji+1,jj  ,jk,jn,Kmm)   ! 2 x C2 interpolation of T at u- & v-points (x2)
896               zC2t_v = pt(ji,jj,jk,jn,Kmm) + pt(ji  ,jj+1,jk,jn,Kmm)
897               !                                                  ! C4 interpolation of T at u- & v-points (x2)
898               zC4t_u =  zC2t_u + r1_6 * ( ztu_im1 - ztu_ip1 )
899               zC4t_v =  zC2t_v + r1_6 * ( ztv_jm1 - ztv_jp1 )
900               !                                                  ! C4 minus upstream advective fluxes
901               zwx_3d(ji,jj,jk) =  0.5_wp * pU(ji,jj,jk) * zC4t_u - zwx_3d(ji,jj,jk)
902               zwy_3d(ji,jj,jk) =  0.5_wp * pV(ji,jj,jk) * zC4t_v - zwy_3d(ji,jj,jk)
903            END_3D
904            CALL lbc_lnk( 'traadv_fct', zwx_3d, 'U', -1.0_wp , zwy_3d, 'V', -1.0_wp )   ! Lateral boundary cond. (unchanged sgn)
905            !
906         END SELECT
907         !
908         SELECT CASE( kn_fct_v )    !* vertical anti-diffusive fluxes (w-masked interior values)
909         !
910         CASE(  2  )                   !- 2nd order centered
911            DO_3D( 1, 1, 1, 1, 2, jpkm1 )
912               zwz(ji,jj,jk) =  (  pW(ji,jj,jk) * 0.5_wp * ( pt(ji,jj,jk,jn,Kmm) + pt(ji,jj,jk-1,jn,Kmm) )   &
913                  &              - zwz(ji,jj,jk)  ) * wmask(ji,jj,jk)
914            END_3D
915            !
916         CASE(  4  )                   !- 4th order COMPACT
917            CALL interp_4th_cpt( pt(:,:,:,jn,Kmm) , ztw )   ! zwt = COMPACT interpolation of T at w-point
918            DO_3D( 1, 1, 1, 1, 2, jpkm1 )
919               zwz(ji,jj,jk) = ( pW(ji,jj,jk) * ztw(ji,jj,jk) - zwz(ji,jj,jk) ) * wmask(ji,jj,jk)
920            END_3D
921            !
922         END SELECT
923         IF( ln_linssh ) THEN    ! top ocean value: high order = upstream  ==>>  zwz=0
924            zwz(:,:,1) = 0._wp   ! only ocean surface as interior zwz values have been w-masked
925         ENDIF
926         !
927         CALL lbc_lnk( 'traadv_fct', zwi, 'T', 1.0_wp)
928         !
929         IF ( ll_zAimp ) THEN
930            DO_3D( 1, 1, 1, 1, 1, jpkm1 )    !* trend and after field with monotonic scheme
931               !                                                ! total intermediate advective trends
932               ztra = - (  zwx_3d(ji,jj,jk) - zwx_3d(ji-1,jj  ,jk  )   &
933                  &      + zwy_3d(ji,jj,jk) - zwy_3d(ji  ,jj-1,jk  )   &
934                  &      + zwz(ji,jj,jk) - zwz(ji  ,jj  ,jk+1) ) * r1_e1e2t(ji,jj)
935               ztw(ji,jj,jk) = zwi(ji,jj,jk) + p2dt * ztra / e3t(ji,jj,jk,Krhs) * tmask(ji,jj,jk)
936            END_3D
937            !
938            CALL tridia_solver( zwdia, zwsup, zwinf, ztw, ztw , 0 )
939            !
940            DO_3D( 1, 1, 1, 1, 2, jpkm1 )       ! Interior value ( multiplied by wmask)
941               zfp_wk = wi(ji,jj,jk) + ABS( wi(ji,jj,jk) )
942               zfm_wk = wi(ji,jj,jk) - ABS( wi(ji,jj,jk) )
943               zwz(ji,jj,jk) = zwz(ji,jj,jk) + 0.5 * e1e2t(ji,jj) * ( zfp_wk * ztw(ji,jj,jk) + zfm_wk * ztw(ji,jj,jk-1) ) * wmask(ji,jj,jk)
944            END_3D
945         END IF
946         !
947         !        !==  monotonicity algorithm  ==!
948         !
949         CALL nonosc( Kmm, pt(:,:,:,jn,Kbb), zwx_3d, zwy_3d, zwz, zwi, p2dt )
950         !
951         !        !==  final trend with corrected fluxes  ==!
952         !
953         DO_3D( 0, 0, 0, 0, 1, jpkm1 )
954            ztra = - (  zwx_3d(ji,jj,jk) - zwx_3d(ji-1,jj  ,jk  )   &
955               &      + zwy_3d(ji,jj,jk) - zwy_3d(ji  ,jj-1,jk  )   &
956               &      + zwz(ji,jj,jk) - zwz(ji  ,jj  ,jk+1) ) * r1_e1e2t(ji,jj)
957            pt(ji,jj,jk,jn,Krhs) = pt(ji,jj,jk,jn,Krhs) + ztra / e3t(ji,jj,jk,Kmm)
958            zwi(ji,jj,jk) = zwi(ji,jj,jk) + p2dt * ztra / e3t(ji,jj,jk,Krhs) * tmask(ji,jj,jk)
959         END_3D
960         !
961         IF ( ll_zAimp ) THEN
962            !
963            ztw(:,:,1) = 0._wp ; ztw(:,:,jpk) = 0._wp
964            DO_3D( 0, 0, 0, 0, 2, jpkm1 )      ! Interior value ( multiplied by wmask)
965               zfp_wk = wi(ji,jj,jk) + ABS( wi(ji,jj,jk) )
966               zfm_wk = wi(ji,jj,jk) - ABS( wi(ji,jj,jk) )
967               ztw(ji,jj,jk) = - 0.5 * e1e2t(ji,jj) * ( zfp_wk * zwi(ji,jj,jk) + zfm_wk * zwi(ji,jj,jk-1) ) * wmask(ji,jj,jk)
968               zwz(ji,jj,jk) = zwz(ji,jj,jk) + ztw(ji,jj,jk) ! Update vertical fluxes for trend diagnostic
969            END_3D
970            DO_3D( 0, 0, 0, 0, 1, jpkm1 )
971               pt(ji,jj,jk,jn,Krhs) = pt(ji,jj,jk,jn,Krhs) - ( ztw(ji,jj,jk) - ztw(ji  ,jj  ,jk+1) ) &
972                  &                                        * r1_e1e2t(ji,jj) / e3t(ji,jj,jk,Kmm)
973            END_3D
974         END IF
975         ! NOT TESTED - NEED l_trd OR l_hst TRUE
976         IF( l_trd .OR. l_hst ) THEN   ! trend diagnostics // heat/salt transport
977            ztrdx(:,:,:) = ztrdx(:,:,:) + zwx_3d(:,:,:)  ! <<< add anti-diffusive fluxes
978            ztrdy(:,:,:) = ztrdy(:,:,:) + zwy_3d(:,:,:)  !     to upstream fluxes
979            ztrdz(:,:,:) = ztrdz(:,:,:) + zwz(:,:,:)  !
980            !
981            IF( l_trd ) THEN              ! trend diagnostics
982               CALL trd_tra( kt, Kmm, Krhs, cdtype, jn, jptra_xad, ztrdx, pU, pt(:,:,:,jn,Kmm) )
983               CALL trd_tra( kt, Kmm, Krhs, cdtype, jn, jptra_yad, ztrdy, pV, pt(:,:,:,jn,Kmm) )
984               CALL trd_tra( kt, Kmm, Krhs, cdtype, jn, jptra_zad, ztrdz, pW, pt(:,:,:,jn,Kmm) )
985            ENDIF
986            !                             ! heat/salt transport
987            IF( l_hst )   CALL dia_ar5_hst( jn, 'adv', ztrdx(:,:,:), ztrdy(:,:,:) )
988            !
989         ENDIF
990         ! NOT TESTED - NEED l_ptr TRUE
991         IF( l_ptr ) THEN              ! "Poleward" transports
992            zptry(:,:,:) = zptry(:,:,:) + zwy_3d(:,:,:)  ! <<< add anti-diffusive fluxes
993            CALL dia_ptr_hst( jn, 'adv', zptry(:,:,:) )
994         ENDIF
995         !
996      END DO                     ! end of tracer loop
997      !
998      IF ( ll_zAimp ) THEN
999         DEALLOCATE( zwdia, zwinf, zwsup )
1000      ENDIF
1001      IF( l_trd .OR. l_hst ) THEN
1002         DEALLOCATE( ztrdx, ztrdy, ztrdz )
1003      ENDIF
1004      IF( l_ptr ) THEN
1005         DEALLOCATE( zptry )
1006      ENDIF
1007      !
1008   END SUBROUTINE tra_adv_fct_lf
1009#endif
[3]1010   !!======================================================================
[5770]1011END MODULE traadv_fct
Note: See TracBrowser for help on using the repository browser.