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

source: trunk/NEMOGCM/NEMO/OPA_SRC/TRA/traadv_fct.F90 @ 6505

Last change on this file since 6505 was 6140, checked in by timgraham, 8 years ago

Merge of branches/2015/dev_merge_2015 back into trunk. Merge excludes NEMOGCM/TOOLS/OBSTOOLS/ for now due to issues with the change of file type. Will sort these manually with further commits.

Branch merged as follows:
In the working copy of branch ran:
svn merge svn+ssh://forge.ipsl.jussieu.fr/ipsl/forge/projets/nemo/svn/trunk@HEAD
Small conflicts due to bug fixes applied to trunk since the dev_merge_2015 was copied. Bug fixes were applied to the branch as well so these were easy to resolve.
Branch committed at this stage

In working copy run:
svn switch svn+ssh://forge.ipsl.jussieu.fr/ipsl/forge/projets/nemo/svn/trunk
to switch working copy

Run:
svn merge --reintegrate svn+ssh://forge.ipsl.jussieu.fr/ipsl/forge/projets/nemo/svn/branches/2015/dev_merge_2015
to merge the branch into the trunk and then commit - no conflicts at this stage.

  • Property svn:keywords set to Id
File size: 40.2 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   !!  tra_adv_fct_zts: update the tracer trend with a 3D advective trends using a 2nd order FCT scheme
12   !!                   with sub-time-stepping in the vertical direction
13   !!  nonosc         : compute monotonic tracer fluxes by a non-oscillatory algorithm
14   !!  interp_4th_cpt : 4th order compact scheme for the vertical component of the advection
[3]15   !!----------------------------------------------------------------------
[3625]16   USE oce            ! ocean dynamics and active tracers
17   USE dom_oce        ! ocean space and time domain
[4990]18   USE trc_oce        ! share passive tracers/Ocean variables
19   USE trd_oce        ! trends: ocean variables
[3625]20   USE trdtra         ! tracers trends
[4990]21   USE diaptr         ! poleward transport diagnostics
22   !
[5770]23   USE in_out_manager ! I/O manager
[3625]24   USE lib_mpp        ! MPP library
25   USE lbclnk         ! ocean lateral boundary condition (or mpp link)
[5770]26   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 
[3625]27   USE wrk_nemo       ! Memory Allocation
28   USE timing         ! Timing
[3]29
30   IMPLICIT NONE
31   PRIVATE
32
[5770]33   PUBLIC   tra_adv_fct        ! routine called by traadv.F90
34   PUBLIC   tra_adv_fct_zts    ! routine called by traadv.F90
35   PUBLIC   interp_4th_cpt     ! routine called by traadv_cen.F90
[3]36
[5770]37   LOGICAL  ::   l_trd   ! flag to compute trends
38   REAL(wp) ::   r1_6 = 1._wp / 6._wp   ! =1/6
[2528]39
[3]40   !! * Substitutions
41#  include "vectopt_loop_substitute.h90"
42   !!----------------------------------------------------------------------
[5770]43   !! NEMO/OPA 3.7 , NEMO Consortium (2014)
[1152]44   !! $Id$
[2528]45   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
[3]46   !!----------------------------------------------------------------------
47CONTAINS
48
[5770]49   SUBROUTINE tra_adv_fct( kt, kit000, cdtype, p2dt, pun, pvn, pwn,       &
50      &                                              ptb, ptn, pta, kjpt, kn_fct_h, kn_fct_v )
[3]51      !!----------------------------------------------------------------------
[5770]52      !!                  ***  ROUTINE tra_adv_fct  ***
[3]53      !!
[6140]54      !! **  Purpose :   Compute the now trend due to total advection of tracers
55      !!               and add it to the general trend of tracer equations
[3]56      !!
[5770]57      !! **  Method  : - 2nd or 4th FCT scheme on the horizontal direction
58      !!               (choice through the value of kn_fct)
[6140]59      !!               - on the vertical the 4th order is a compact scheme
[5770]60      !!               - corrected flux (monotonic correction)
[3]61      !!
[6140]62      !! ** Action : - update pta  with the now advective tracer trends
63      !!             - send trends to trdtra module for further diagnostcs (l_trdtra=T)
64      !!             - htr_adv, str_adv : poleward advective heat and salt transport (ln_diaptr=T)
[503]65      !!----------------------------------------------------------------------
[2528]66      INTEGER                              , INTENT(in   ) ::   kt              ! ocean time-step index
[3294]67      INTEGER                              , INTENT(in   ) ::   kit000          ! first time step index
[2528]68      CHARACTER(len=3)                     , INTENT(in   ) ::   cdtype          ! =TRA or TRC (tracer indicator)
69      INTEGER                              , INTENT(in   ) ::   kjpt            ! number of tracers
[5770]70      INTEGER                              , INTENT(in   ) ::   kn_fct_h        ! order of the FCT scheme (=2 or 4)
71      INTEGER                              , INTENT(in   ) ::   kn_fct_v        ! order of the FCT scheme (=2 or 4)
[6140]72      REAL(wp)                             , INTENT(in   ) ::   p2dt            ! tracer time-step
[2528]73      REAL(wp), DIMENSION(jpi,jpj,jpk     ), INTENT(in   ) ::   pun, pvn, pwn   ! 3 ocean velocity components
74      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(in   ) ::   ptb, ptn        ! before and now tracer fields
75      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::   pta             ! tracer trend
[2715]76      !
[5770]77      INTEGER  ::   ji, jj, jk, jn                           ! dummy loop indices 
[6140]78      REAL(wp) ::   ztra                                     ! local scalar
[5770]79      REAL(wp) ::   zfp_ui, zfp_vj, zfp_wk, zC2t_u, zC4t_u   !   -      -
80      REAL(wp) ::   zfm_ui, zfm_vj, zfm_wk, zC2t_v, zC4t_v   !   -      -
81      REAL(wp), POINTER, DIMENSION(:,:,:) ::   zwi, zwx, zwy, zwz, ztu, ztv, zltu, zltv, ztw
82      REAL(wp), POINTER, DIMENSION(:,:,:) ::   ztrdx, ztrdy, ztrdz
[3]83      !!----------------------------------------------------------------------
[3294]84      !
[5770]85      IF( nn_timing == 1 )  CALL timing_start('tra_adv_fct')
[3294]86      !
[5770]87      CALL wrk_alloc( jpi,jpj,jpk,   zwi, zwx, zwy, zwz, ztu, ztv, zltu, zltv, ztw )
[3294]88      !
89      IF( kt == kit000 )  THEN
[2528]90         IF(lwp) WRITE(numout,*)
[5770]91         IF(lwp) WRITE(numout,*) 'tra_adv_fct : FCT advection scheme on ', cdtype
[2528]92         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
[3]93      ENDIF
[2528]94      !
[5770]95      l_trd = .FALSE.
96      IF( ( cdtype == 'TRA' .AND. l_trdtra ) .OR. ( cdtype == 'TRC' .AND. l_trdtrc ) )   l_trd = .TRUE.
97      !
[2528]98      IF( l_trd )  THEN
[3294]99         CALL wrk_alloc( jpi, jpj, jpk, ztrdx, ztrdy, ztrdz )
[5770]100         ztrdx(:,:,:) = 0._wp   ;    ztrdy(:,:,:) = 0._wp   ;   ztrdz(:,:,:) = 0._wp
[3294]101      ENDIF
[2528]102      !
[6140]103      !                          ! surface & bottom value : flux set to zero one for all
104      zwz(:,:, 1 ) = 0._wp           
[5770]105      zwx(:,:,jpk) = 0._wp   ;   zwy(:,:,jpk) = 0._wp    ;    zwz(:,:,jpk) = 0._wp
[2528]106      !
[5770]107      zwi(:,:,:) = 0._wp       
[6140]108      !
109      DO jn = 1, kjpt            !==  loop over the tracers  ==!
[5770]110         !
111         !        !==  upstream advection with initial mass fluxes & intermediate update  ==!
112         !                    !* upstream tracer flux in the i and j direction
[2528]113         DO jk = 1, jpkm1
114            DO jj = 1, jpjm1
115               DO ji = 1, fs_jpim1   ! vector opt.
116                  ! upstream scheme
117                  zfp_ui = pun(ji,jj,jk) + ABS( pun(ji,jj,jk) )
118                  zfm_ui = pun(ji,jj,jk) - ABS( pun(ji,jj,jk) )
119                  zfp_vj = pvn(ji,jj,jk) + ABS( pvn(ji,jj,jk) )
120                  zfm_vj = pvn(ji,jj,jk) - ABS( pvn(ji,jj,jk) )
121                  zwx(ji,jj,jk) = 0.5 * ( zfp_ui * ptb(ji,jj,jk,jn) + zfm_ui * ptb(ji+1,jj  ,jk,jn) )
122                  zwy(ji,jj,jk) = 0.5 * ( zfp_vj * ptb(ji,jj,jk,jn) + zfm_vj * ptb(ji  ,jj+1,jk,jn) )
123               END DO
[3]124            END DO
125         END DO
[5770]126         !                    !* upstream tracer flux in the k direction *!
[6140]127         DO jk = 2, jpkm1        ! Interior value ( multiplied by wmask)
[4990]128            DO jj = 1, jpj
129               DO ji = 1, jpi
[2528]130                  zfp_wk = pwn(ji,jj,jk) + ABS( pwn(ji,jj,jk) )
131                  zfm_wk = pwn(ji,jj,jk) - ABS( pwn(ji,jj,jk) )
[5120]132                  zwz(ji,jj,jk) = 0.5 * ( zfp_wk * ptb(ji,jj,jk,jn) + zfm_wk * ptb(ji,jj,jk-1,jn) ) * wmask(ji,jj,jk)
[2528]133               END DO
[3]134            END DO
135         END DO
[6140]136         IF( ln_linssh ) THEN    ! top ocean value (only in linear free surface as zwz has been w-masked)
[5770]137            IF( ln_isfcav ) THEN             ! top of the ice-shelf cavities and at the ocean surface
[5120]138               DO jj = 1, jpj
139                  DO ji = 1, jpi
140                     zwz(ji,jj, mikt(ji,jj) ) = pwn(ji,jj,mikt(ji,jj)) * ptb(ji,jj,mikt(ji,jj),jn)   ! linear free surface
141                  END DO
142               END DO   
[5770]143            ELSE                             ! no cavities: only at the ocean surface
144               zwz(:,:,1) = pwn(:,:,1) * ptb(:,:,1,jn)
145            ENDIF
[5120]146         ENDIF
[5770]147         !               
148         DO jk = 1, jpkm1     !* trend and after field with monotonic scheme
[216]149            DO jj = 2, jpjm1
150               DO ji = fs_2, fs_jpim1   ! vector opt.
[2528]151                  ! total intermediate advective trends
[5770]152                  ztra = - (  zwx(ji,jj,jk) - zwx(ji-1,jj  ,jk  )   &
153                     &      + zwy(ji,jj,jk) - zwy(ji  ,jj-1,jk  )   &
[6140]154                     &      + zwz(ji,jj,jk) - zwz(ji  ,jj  ,jk+1) ) * r1_e1e2t(ji,jj) / e3t_n(ji,jj,jk)
[2528]155                  ! update and guess with monotonic sheme
[5770]156!!gm why tmask added in the two following lines ???    the mask is done in tranxt !
[4990]157                  pta(ji,jj,jk,jn) =   pta(ji,jj,jk,jn)         + ztra   * tmask(ji,jj,jk)
[6140]158                  zwi(ji,jj,jk)    = ( ptb(ji,jj,jk,jn) + p2dt * ztra ) * tmask(ji,jj,jk)
[216]159               END DO
160            END DO
161         END DO
[5770]162         CALL lbc_lnk( zwi, 'T', 1. )  ! Lateral boundary conditions on zwi  (unchanged sign)
163         !               
164         IF( l_trd )  THEN             ! trend diagnostics (contribution of upstream fluxes)
[2528]165            ztrdx(:,:,:) = zwx(:,:,:)   ;    ztrdy(:,:,:) = zwy(:,:,:)  ;   ztrdz(:,:,:) = zwz(:,:,:)
166         END IF
[5770]167         !                             ! "Poleward" heat and salt transports (contribution of upstream fluxes)
[5147]168         IF( cdtype == 'TRA' .AND. ln_diaptr ) THEN 
169           IF( jn == jp_tem )  htr_adv(:) = ptr_sj( zwy(:,:,:) )
170           IF( jn == jp_sal )  str_adv(:) = ptr_sj( zwy(:,:,:) )
[2528]171         ENDIF
[5770]172         !
173         !        !==  anti-diffusive flux : high order minus low order  ==!
174         !
[6140]175         SELECT CASE( kn_fct_h )    !* horizontal anti-diffusive fluxes
[5770]176         !
[6140]177         CASE(  2  )                   !- 2nd order centered
[5770]178            DO jk = 1, jpkm1
179               DO jj = 1, jpjm1
180                  DO ji = 1, fs_jpim1   ! vector opt.
181                     zwx(ji,jj,jk) = 0.5_wp * pun(ji,jj,jk) * ( ptn(ji,jj,jk,jn) + ptn(ji+1,jj,jk,jn) ) - zwx(ji,jj,jk)
182                     zwy(ji,jj,jk) = 0.5_wp * pvn(ji,jj,jk) * ( ptn(ji,jj,jk,jn) + ptn(ji,jj+1,jk,jn) ) - zwy(ji,jj,jk)
183                  END DO
[503]184               END DO
185            END DO
[5770]186            !
[6140]187         CASE(  4  )                   !- 4th order centered
188            zltu(:,:,jpk) = 0._wp            ! Bottom value : flux set to zero
[5770]189            zltv(:,:,jpk) = 0._wp
[6140]190            DO jk = 1, jpkm1                 ! Laplacian
191               DO jj = 1, jpjm1                    ! 1st derivative (gradient)
[5770]192                  DO ji = 1, fs_jpim1   ! vector opt.
193                     ztu(ji,jj,jk) = ( ptn(ji+1,jj  ,jk,jn) - ptn(ji,jj,jk,jn) ) * umask(ji,jj,jk)
194                     ztv(ji,jj,jk) = ( ptn(ji  ,jj+1,jk,jn) - ptn(ji,jj,jk,jn) ) * vmask(ji,jj,jk)
195                  END DO
[503]196               END DO
[6140]197               DO jj = 2, jpjm1                    ! 2nd derivative * 1/ 6
[5770]198                  DO ji = fs_2, fs_jpim1   ! vector opt.
199                     zltu(ji,jj,jk) = (  ztu(ji,jj,jk) + ztu(ji-1,jj,jk)  ) * r1_6
200                     zltv(ji,jj,jk) = (  ztv(ji,jj,jk) + ztv(ji,jj-1,jk)  ) * r1_6
201                  END DO
202               END DO
[503]203            END DO
[5770]204            CALL lbc_lnk( zltu, 'T', 1. )   ;    CALL lbc_lnk( zltv, 'T', 1. )   ! Lateral boundary cond. (unchanged sgn)
205            !
[6140]206            DO jk = 1, jpkm1                 ! Horizontal advective fluxes
[5770]207               DO jj = 1, jpjm1
208                  DO ji = 1, fs_jpim1   ! vector opt.
209                     zC2t_u = ptn(ji,jj,jk,jn) + ptn(ji+1,jj  ,jk,jn)   ! 2 x C2 interpolation of T at u- & v-points
210                     zC2t_v = ptn(ji,jj,jk,jn) + ptn(ji  ,jj+1,jk,jn)
211                     !                                                  ! C4 minus upstream advective fluxes
212                     zwx(ji,jj,jk) =  0.5_wp * pun(ji,jj,jk) * ( zC2t_u + zltu(ji,jj,jk) - zltu(ji+1,jj,jk) ) - zwx(ji,jj,jk)
213                     zwy(ji,jj,jk) =  0.5_wp * pvn(ji,jj,jk) * ( zC2t_v + zltv(ji,jj,jk) - zltv(ji,jj+1,jk) ) - zwy(ji,jj,jk)
214                  END DO
[5120]215               END DO
[5770]216            END DO         
217            !
[6140]218         CASE(  41 )                   !- 4th order centered       ==>>   !!gm coding attempt   need to be tested
219            ztu(:,:,jpk) = 0._wp             ! Bottom value : flux set to zero
[5770]220            ztv(:,:,jpk) = 0._wp
[6140]221            DO jk = 1, jpkm1                 ! 1st derivative (gradient)
222               DO jj = 1, jpjm1
[5770]223                  DO ji = 1, fs_jpim1   ! vector opt.
224                     ztu(ji,jj,jk) = ( ptn(ji+1,jj  ,jk,jn) - ptn(ji,jj,jk,jn) ) * umask(ji,jj,jk)
225                     ztv(ji,jj,jk) = ( ptn(ji  ,jj+1,jk,jn) - ptn(ji,jj,jk,jn) ) * vmask(ji,jj,jk)
226                  END DO
227               END DO
[5120]228            END DO
[5770]229            CALL lbc_lnk( ztu, 'U', -1. )   ;    CALL lbc_lnk( ztv, 'V', -1. )   ! Lateral boundary cond. (unchanged sgn)
230            !
[6140]231            DO jk = 1, jpkm1                 ! Horizontal advective fluxes
[5770]232               DO jj = 2, jpjm1
233                  DO ji = 2, fs_jpim1   ! vector opt.
234                     zC2t_u = ptn(ji,jj,jk,jn) + ptn(ji+1,jj  ,jk,jn)   ! 2 x C2 interpolation of T at u- & v-points (x2)
235                     zC2t_v = ptn(ji,jj,jk,jn) + ptn(ji  ,jj+1,jk,jn)
236                     !                                                  ! C4 interpolation of T at u- & v-points (x2)
237                     zC4t_u =  zC2t_u + r1_6 * ( ztu(ji-1,jj  ,jk) - ztu(ji+1,jj  ,jk) )
238                     zC4t_v =  zC2t_v + r1_6 * ( ztv(ji  ,jj-1,jk) - ztv(ji  ,jj+1,jk) )
239                     !                                                  ! C4 minus upstream advective fluxes
240                     zwx(ji,jj,jk) =  0.5_wp * pun(ji,jj,jk) * zC4t_u - zwx(ji,jj,jk)
241                     zwy(ji,jj,jk) =  0.5_wp * pvn(ji,jj,jk) * zC4t_v - zwy(ji,jj,jk)
242                  END DO
243               END DO
244            END DO
245            !
246         END SELECT
[6140]247         !                     
248         SELECT CASE( kn_fct_v )    !* vertical anti-diffusive fluxes (w-masked interior values)
[5770]249         !
[6140]250         CASE(  2  )                   !- 2nd order centered
[5770]251            DO jk = 2, jpkm1   
252               DO jj = 2, jpjm1
253                  DO ji = fs_2, fs_jpim1
[6140]254                     zwz(ji,jj,jk) =  (  pwn(ji,jj,jk) * 0.5_wp * ( ptn(ji,jj,jk,jn) + ptn(ji,jj,jk-1,jn) )   &
255                        &              - zwz(ji,jj,jk)  ) * wmask(ji,jj,jk)
[5770]256                  END DO
257               END DO
258            END DO
259            !
[6140]260         CASE(  4  )                   !- 4th order COMPACT
261            CALL interp_4th_cpt( ptn(:,:,:,jn) , ztw )   ! zwt = COMPACT interpolation of T at w-point
[5770]262            DO jk = 2, jpkm1
263               DO jj = 2, jpjm1
264                  DO ji = fs_2, fs_jpim1
265                     zwz(ji,jj,jk) = ( pwn(ji,jj,jk) * ztw(ji,jj,jk) - zwz(ji,jj,jk) ) * wmask(ji,jj,jk)
266                  END DO
267               END DO
268            END DO
269            !
270         END SELECT
[6140]271         IF( ln_linssh ) THEN    ! top ocean value: high order = upstream  ==>>  zwz=0
272            zwz(:,:,1) = 0._wp   ! only ocean surface as interior zwz values have been w-masked
273         ENDIF
[5770]274         !
[2528]275         CALL lbc_lnk( zwx, 'U', -1. )   ;   CALL lbc_lnk( zwy, 'V', -1. )         ! Lateral bondary conditions
276         CALL lbc_lnk( zwz, 'W',  1. )
[6140]277         !
[5770]278         !        !==  monotonicity algorithm  ==!
279         !
[2528]280         CALL nonosc( ptb(:,:,:,jn), zwx, zwy, zwz, zwi, p2dt )
[6140]281         !
[5770]282         !        !==  final trend with corrected fluxes  ==!
283         !
[216]284         DO jk = 1, jpkm1
285            DO jj = 2, jpjm1
[2528]286               DO ji = fs_2, fs_jpim1   ! vector opt. 
[5770]287                  pta(ji,jj,jk,jn) = pta(ji,jj,jk,jn) - (  zwx(ji,jj,jk) - zwx(ji-1,jj  ,jk  )   &
288                     &                                   + zwy(ji,jj,jk) - zwy(ji  ,jj-1,jk  )   &
289                     &                                   + zwz(ji,jj,jk) - zwz(ji  ,jj  ,jk+1) ) &
[6140]290                     &                                * r1_e1e2t(ji,jj) / e3t_n(ji,jj,jk)
[216]291               END DO
292            END DO
293         END DO
[5770]294         !
[6140]295         IF( l_trd ) THEN     ! trend diagnostics (contribution of upstream fluxes)
[2528]296            ztrdx(:,:,:) = ztrdx(:,:,:) + zwx(:,:,:)  ! <<< Add to previously computed
297            ztrdy(:,:,:) = ztrdy(:,:,:) + zwy(:,:,:)  ! <<< Add to previously computed
298            ztrdz(:,:,:) = ztrdz(:,:,:) + zwz(:,:,:)  ! <<< Add to previously computed
[5770]299            !
300            CALL trd_tra( kt, cdtype, jn, jptra_xad, ztrdx, pun, ptn(:,:,:,jn) )
301            CALL trd_tra( kt, cdtype, jn, jptra_yad, ztrdy, pvn, ptn(:,:,:,jn) )
302            CALL trd_tra( kt, cdtype, jn, jptra_zad, ztrdz, pwn, ptn(:,:,:,jn) )
303            !
304            CALL wrk_dealloc( jpi,jpj,jpk,   ztrdx, ztrdy, ztrdz )
[2528]305         END IF
[6140]306         !                    ! "Poleward" heat and salt transports (contribution of upstream fluxes)
[5147]307         IF( cdtype == 'TRA' .AND. ln_diaptr ) THEN 
[6140]308           IF( jn == jp_tem )   htr_adv(:) = htr_adv(:) + ptr_sj( zwy(:,:,:) )
309           IF( jn == jp_sal )   str_adv(:) = str_adv(:) + ptr_sj( zwy(:,:,:) )
[2528]310         ENDIF
[503]311         !
[6140]312      END DO                     ! end of tracer loop
[503]313      !
[5770]314      CALL wrk_dealloc( jpi,jpj,jpk,    zwi, zwx, zwy, zwz, ztu, ztv, zltu, zltv, ztw )
[2528]315      !
[5770]316      IF( nn_timing == 1 )  CALL timing_stop('tra_adv_fct')
[2715]317      !
[5770]318   END SUBROUTINE tra_adv_fct
[3]319
[5737]320
[5770]321   SUBROUTINE tra_adv_fct_zts( kt, kit000, cdtype, p2dt, pun, pvn, pwn,      &
322      &                                                  ptb, ptn, pta, kjpt, kn_fct_zts )
[4990]323      !!----------------------------------------------------------------------
[5770]324      !!                  ***  ROUTINE tra_adv_fct_zts  ***
[4990]325      !!
326      !! **  Purpose :   Compute the now trend due to total advection of
327      !!       tracers and add it to the general trend of tracer equations
328      !!
329      !! **  Method  :   TVD ZTS scheme, i.e. 2nd order centered scheme with
330      !!       corrected flux (monotonic correction). This version use sub-
331      !!       timestepping for the vertical advection which increases stability
332      !!       when vertical metrics are small.
333      !!       note: - this advection scheme needs a leap-frog time scheme
334      !!
335      !! ** Action : - update (pta) with the now advective tracer trends
336      !!             - save the trends
337      !!----------------------------------------------------------------------
338      INTEGER                              , INTENT(in   ) ::   kt              ! ocean time-step index
339      INTEGER                              , INTENT(in   ) ::   kit000          ! first time step index
340      CHARACTER(len=3)                     , INTENT(in   ) ::   cdtype          ! =TRA or TRC (tracer indicator)
341      INTEGER                              , INTENT(in   ) ::   kjpt            ! number of tracers
[5770]342      INTEGER                              , INTENT(in   ) ::   kn_fct_zts      ! number of number of vertical sub-timesteps
[6140]343      REAL(wp)                             , INTENT(in   ) ::   p2dt            ! tracer time-step
[4990]344      REAL(wp), DIMENSION(jpi,jpj,jpk     ), INTENT(in   ) ::   pun, pvn, pwn   ! 3 ocean velocity components
345      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(in   ) ::   ptb, ptn        ! before and now tracer fields
346      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::   pta             ! tracer trend
347      !
348      REAL(wp), DIMENSION( jpk )                           ::   zts             ! length of sub-timestep for vertical advection
[6140]349      REAL(wp)                                             ::   zr_p2dt         ! reciprocal of tracer timestep
[4990]350      INTEGER  ::   ji, jj, jk, jl, jn       ! dummy loop indices 
351      INTEGER  ::   jtb, jtn, jta   ! sub timestep pointers for leap-frog/euler forward steps
352      INTEGER  ::   jtaken          ! toggle for collecting appropriate fluxes from sub timesteps
353      REAL(wp) ::   z_rzts          ! Fractional length of Euler forward sub-timestep for vertical advection
[6140]354      REAL(wp) ::   ztra            ! local scalar
[4990]355      REAL(wp) ::   zfp_ui, zfp_vj, zfp_wk   !   -      -
356      REAL(wp) ::   zfm_ui, zfm_vj, zfm_wk   !   -      -
[5770]357      REAL(wp), POINTER, DIMENSION(:,:  )   ::   zwx_sav , zwy_sav
358      REAL(wp), POINTER, DIMENSION(:,:,:)   ::   zwi, zwx, zwy, zwz, zhdiv, zwzts, zwz_sav
359      REAL(wp), POINTER, DIMENSION(:,:,:)   ::   ztrdx, ztrdy, ztrdz
360      REAL(wp), POINTER, DIMENSION(:,:,:,:) ::   ztrs
[4990]361      !!----------------------------------------------------------------------
362      !
[5770]363      IF( nn_timing == 1 )  CALL timing_start('tra_adv_fct_zts')
[4990]364      !
[5770]365      CALL wrk_alloc( jpi,jpj,             zwx_sav, zwy_sav )
366      CALL wrk_alloc( jpi,jpj, jpk,        zwx, zwy, zwz, zwi, zhdiv, zwzts, zwz_sav )
367      CALL wrk_alloc( jpi,jpj,jpk,kjpt+1,  ztrs )
[4990]368      !
369      IF( kt == kit000 )  THEN
370         IF(lwp) WRITE(numout,*)
[5770]371         IF(lwp) WRITE(numout,*) 'tra_adv_fct_zts : 2nd order FCT scheme with ', kn_fct_zts, ' vertical sub-timestep on ', cdtype
[4990]372         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
373      ENDIF
374      !
375      l_trd = .FALSE.
[5770]376      IF( ( cdtype == 'TRA' .AND. l_trdtra ) .OR. ( cdtype == 'TRC' .AND. l_trdtrc ) )   l_trd = .TRUE.
[4990]377      !
378      IF( l_trd )  THEN
[5770]379         CALL wrk_alloc( jpi,jpj,jpk,   ztrdx, ztrdy, ztrdz )
[4990]380         ztrdx(:,:,:) = 0._wp  ;    ztrdy(:,:,:) = 0._wp  ;   ztrdz(:,:,:) = 0._wp
381      ENDIF
382      !
383      zwi(:,:,:) = 0._wp
[5770]384      z_rzts = 1._wp / REAL( kn_fct_zts, wp )
[6140]385      zr_p2dt = 1._wp / p2dt
[4990]386      !
[6140]387      ! surface & Bottom value : flux set to zero for all tracers
388      zwz(:,:, 1 ) = 0._wp
389      zwx(:,:,jpk) = 0._wp   ;    zwz(:,:,jpk) = 0._wp
390      zwy(:,:,jpk) = 0._wp   ;    zwi(:,:,jpk) = 0._wp
391      !
[4990]392      !                                                          ! ===========
393      DO jn = 1, kjpt                                            ! tracer loop
394         !                                                       ! ===========
[6140]395         !
396         ! Upstream advection with initial mass fluxes & intermediate update
397         DO jk = 1, jpkm1        ! upstream tracer flux in the i and j direction
[4990]398            DO jj = 1, jpjm1
399               DO ji = 1, fs_jpim1   ! vector opt.
400                  ! upstream scheme
401                  zfp_ui = pun(ji,jj,jk) + ABS( pun(ji,jj,jk) )
402                  zfm_ui = pun(ji,jj,jk) - ABS( pun(ji,jj,jk) )
403                  zfp_vj = pvn(ji,jj,jk) + ABS( pvn(ji,jj,jk) )
404                  zfm_vj = pvn(ji,jj,jk) - ABS( pvn(ji,jj,jk) )
405                  zwx(ji,jj,jk) = 0.5_wp * ( zfp_ui * ptb(ji,jj,jk,jn) + zfm_ui * ptb(ji+1,jj  ,jk,jn) )
406                  zwy(ji,jj,jk) = 0.5_wp * ( zfp_vj * ptb(ji,jj,jk,jn) + zfm_vj * ptb(ji  ,jj+1,jk,jn) )
407               END DO
408            END DO
409         END DO
[6140]410         !                       ! upstream tracer flux in the k direction
411         DO jk = 2, jpkm1              ! Interior value
[4990]412            DO jj = 1, jpj
413               DO ji = 1, jpi
414                  zfp_wk = pwn(ji,jj,jk) + ABS( pwn(ji,jj,jk) )
415                  zfm_wk = pwn(ji,jj,jk) - ABS( pwn(ji,jj,jk) )
[5770]416                  zwz(ji,jj,jk) = 0.5_wp * ( zfp_wk * ptb(ji,jj,jk,jn) + zfm_wk * ptb(ji,jj,jk-1,jn) ) * wmask(ji,jj,jk)
[4990]417               END DO
418            END DO
419         END DO
[6140]420         IF( ln_linssh ) THEN          ! top value : linear free surface case only (as zwz is multiplied by wmask)
421            IF( ln_isfcav ) THEN             ! ice-shelf cavities: top value
[5120]422               DO jj = 1, jpj
423                  DO ji = 1, jpi
[5770]424                     zwz(ji,jj, mikt(ji,jj) ) = pwn(ji,jj,mikt(ji,jj)) * ptb(ji,jj,mikt(ji,jj),jn) 
[5120]425                  END DO
[5770]426               END DO   
[6140]427            ELSE                             ! no cavities, surface value
[5770]428               zwz(:,:,1) = pwn(:,:,1) * ptb(:,:,1,jn)
429            ENDIF
[5120]430         ENDIF
[5770]431         !
432         DO jk = 1, jpkm1         ! total advective trend
[4990]433            DO jj = 2, jpjm1
434               DO ji = fs_2, fs_jpim1   ! vector opt.
[6140]435                  !                             ! total intermediate advective trends
[5770]436                  ztra = - (  zwx(ji,jj,jk) - zwx(ji-1,jj  ,jk  )   &
437                     &      + zwy(ji,jj,jk) - zwy(ji  ,jj-1,jk  )   &
[6140]438                     &      + zwz(ji,jj,jk) - zwz(ji  ,jj  ,jk+1)   ) * r1_e1e2t(ji,jj) / e3t_n(ji,jj,jk)
439                  !                             ! update and guess with monotonic sheme
[4990]440                  pta(ji,jj,jk,jn) =   pta(ji,jj,jk,jn)         + ztra
[6140]441                  zwi(ji,jj,jk)    = ( ptb(ji,jj,jk,jn) + p2dt * ztra ) * tmask(ji,jj,jk)
[4990]442               END DO
443            END DO
444         END DO
[5770]445         !                           
446         CALL lbc_lnk( zwi, 'T', 1. )     ! Lateral boundary conditions on zwi  (unchanged sign)
447         !               
448         IF( l_trd )  THEN                ! trend diagnostics (contribution of upstream fluxes)
[4990]449            ztrdx(:,:,:) = zwx(:,:,:)   ;    ztrdy(:,:,:) = zwy(:,:,:)  ;   ztrdz(:,:,:) = zwz(:,:,:)
450         END IF
[5770]451         !                                ! "Poleward" heat and salt transports (contribution of upstream fluxes)
[5147]452         IF( cdtype == 'TRA' .AND. ln_diaptr ) THEN 
453           IF( jn == jp_tem )  htr_adv(:) = ptr_sj( zwy(:,:,:) )
454           IF( jn == jp_sal )  str_adv(:) = ptr_sj( zwy(:,:,:) )
[4990]455         ENDIF
456
[5770]457         ! 3. anti-diffusive flux : high order minus low order
458         ! ---------------------------------------------------
[4990]459
[5770]460         DO jk = 1, jpkm1                    !* horizontal anti-diffusive fluxes
461            !
[4990]462            DO jj = 1, jpjm1
463               DO ji = 1, fs_jpim1   ! vector opt.
464                  zwx_sav(ji,jj) = zwx(ji,jj,jk)
465                  zwy_sav(ji,jj) = zwy(ji,jj,jk)
[5770]466                  !
[4990]467                  zwx(ji,jj,jk) = 0.5_wp * pun(ji,jj,jk) * ( ptn(ji,jj,jk,jn) + ptn(ji+1,jj,jk,jn) )
468                  zwy(ji,jj,jk) = 0.5_wp * pvn(ji,jj,jk) * ( ptn(ji,jj,jk,jn) + ptn(ji,jj+1,jk,jn) )
469               END DO
470            END DO
[5770]471            !
472            DO jj = 2, jpjm1                    ! partial horizontal divergence
[4990]473               DO ji = fs_2, fs_jpim1
474                  zhdiv(ji,jj,jk) = (  zwx(ji,jj,jk) - zwx(ji-1,jj  ,jk)   &
475                     &               + zwy(ji,jj,jk) - zwy(ji  ,jj-1,jk)  )
476               END DO
477            END DO
[5770]478            !
[4990]479            DO jj = 1, jpjm1
480               DO ji = 1, fs_jpim1   ! vector opt.
[5770]481                  zwx(ji,jj,jk) = zwx(ji,jj,jk) - zwx_sav(ji,jj)
482                  zwy(ji,jj,jk) = zwy(ji,jj,jk) - zwy_sav(ji,jj)
[4990]483               END DO
484            END DO
485         END DO
[6140]486         !
[5770]487         !                                !* vertical anti-diffusive flux
488         zwz_sav(:,:,:)   = zwz(:,:,:)
489         ztrs   (:,:,:,1) = ptb(:,:,:,jn)
490         zwzts  (:,:,:)   = 0._wp
[4990]491         !
[5770]492         DO jl = 1, kn_fct_zts                  ! Start of sub timestepping loop
493            !
494            IF( jl == 1 ) THEN                        ! Euler forward to kick things off
495               jtb = 1   ;   jtn = 1   ;   jta = 2
[6140]496               zts(:) = p2dt * z_rzts
[5770]497               jtaken = MOD( kn_fct_zts + 1 , 2)            ! Toggle to collect every second flux
498               !                                            ! starting at jl =1 if kn_fct_zts is odd;
499               !                                            ! starting at jl =2 otherwise
500            ELSEIF( jl == 2 ) THEN                    ! First leapfrog step
501               jtb = 1   ;   jtn = 2   ;   jta = 3
[6140]502               zts(:) = 2._wp * p2dt * z_rzts
[5770]503            ELSE                                      ! Shuffle pointers for subsequent leapfrog steps
504               jtb = MOD(jtb,3) + 1
505               jtn = MOD(jtn,3) + 1
506               jta = MOD(jta,3) + 1
[4990]507            ENDIF
[5770]508            DO jk = 2, jpkm1                          ! interior value
[4990]509               DO jj = 2, jpjm1
510                  DO ji = fs_2, fs_jpim1
[5770]511                     zwz(ji,jj,jk) = 0.5_wp * pwn(ji,jj,jk) * ( ztrs(ji,jj,jk,jtn) + ztrs(ji,jj,jk-1,jtn) ) * wmask(ji,jj,jk)
512                     IF( jtaken == 0 )   zwzts(ji,jj,jk) = zwzts(ji,jj,jk) + zwz(ji,jj,jk) * zts(jk)    ! Accumulate time-weighted vertcal flux
[4990]513                  END DO
514               END DO
515            END DO
[6140]516            IF( ln_linssh ) THEN                    ! top value (only in linear free surface case)
[5770]517               IF( ln_isfcav ) THEN                      ! ice-shelf cavities
518                  DO jj = 1, jpj
519                     DO ji = 1, jpi
520                        zwz(ji,jj, mikt(ji,jj) ) = pwn(ji,jj,mikt(ji,jj)) * ptb(ji,jj,mikt(ji,jj),jn)   ! linear free surface
521                     END DO
522                  END DO   
[6140]523               ELSE                                      ! no ocean cavities
[5770]524                  zwz(:,:,1) = pwn(:,:,1) * ptb(:,:,1,jn)
525               ENDIF
526            ENDIF
527            !
[4990]528            jtaken = MOD( jtaken + 1 , 2 )
[5770]529            !
530            DO jk = 2, jpkm1                             ! total advective trends
[4990]531               DO jj = 2, jpjm1
532                  DO ji = fs_2, fs_jpim1
[5770]533                     ztrs(ji,jj,jk,jta) = ztrs(ji,jj,jk,jtb)                                                 &
534                        &               - zts(jk) * (  zhdiv(ji,jj,jk) + zwz(ji,jj,jk) - zwz(ji,jj,jk+1) )   &
[6140]535                        &                         * r1_e1e2t(ji,jj) / e3t_n(ji,jj,jk)
[4990]536                  END DO
537               END DO
538            END DO
[5770]539            !
[4990]540         END DO
541
542         DO jk = 2, jpkm1          ! Anti-diffusive vertical flux using average flux from the sub-timestepping
543            DO jj = 2, jpjm1
544               DO ji = fs_2, fs_jpim1
[6140]545                  zwz(ji,jj,jk) = ( zwzts(ji,jj,jk) * zr_p2dt - zwz_sav(ji,jj,jk) ) * wmask(ji,jj,jk)
[4990]546               END DO
547            END DO
548         END DO
549         CALL lbc_lnk( zwx, 'U', -1. )   ;   CALL lbc_lnk( zwy, 'V', -1. )         ! Lateral bondary conditions
550         CALL lbc_lnk( zwz, 'W',  1. )
551
552         ! 4. monotonicity algorithm
553         ! -------------------------
554         CALL nonosc( ptb(:,:,:,jn), zwx, zwy, zwz, zwi, p2dt )
555
556
557         ! 5. final trend with corrected fluxes
558         ! ------------------------------------
559         DO jk = 1, jpkm1
560            DO jj = 2, jpjm1
561               DO ji = fs_2, fs_jpim1   ! vector opt. 
[5770]562                  pta(ji,jj,jk,jn) = pta(ji,jj,jk,jn) + (   zwy(ji,jj,jk) - zwy(ji  ,jj-1,jk  )       &
563                     &                                    + zwz(ji,jj,jk) - zwz(ji  ,jj  ,jk+1)   )   &
[6140]564                     &                                * r1_e1e2t(ji,jj) / e3t_n(ji,jj,jk)
[4990]565               END DO
566            END DO
567         END DO
568
569         !                                 ! trend diagnostics (contribution of upstream fluxes)
570         IF( l_trd )  THEN
571            ztrdx(:,:,:) = ztrdx(:,:,:) + zwx(:,:,:)  ! <<< Add to previously computed
572            ztrdy(:,:,:) = ztrdy(:,:,:) + zwy(:,:,:)  ! <<< Add to previously computed
573            ztrdz(:,:,:) = ztrdz(:,:,:) + zwz(:,:,:)  ! <<< Add to previously computed
[5770]574            !
[4990]575            CALL trd_tra( kt, cdtype, jn, jptra_xad, ztrdx, pun, ptn(:,:,:,jn) )   
576            CALL trd_tra( kt, cdtype, jn, jptra_yad, ztrdy, pvn, ptn(:,:,:,jn) ) 
577            CALL trd_tra( kt, cdtype, jn, jptra_zad, ztrdz, pwn, ptn(:,:,:,jn) ) 
[5770]578            !
579            CALL wrk_dealloc( jpi,jpj,jpk,   ztrdx, ztrdy, ztrdz )
[4990]580         END IF
581         !                                 ! "Poleward" heat and salt transports (contribution of upstream fluxes)
[5147]582         IF( cdtype == 'TRA' .AND. ln_diaptr ) THEN 
583           IF( jn == jp_tem )  htr_adv(:) = ptr_sj( zwy(:,:,:) ) + htr_adv(:)
584           IF( jn == jp_sal )  str_adv(:) = ptr_sj( zwy(:,:,:) ) + str_adv(:)
[4990]585         ENDIF
586         !
587      END DO
588      !
[5770]589      CALL wrk_alloc( jpi,jpj,             zwx_sav, zwy_sav )
590      CALL wrk_alloc( jpi,jpj, jpk,        zwx, zwy, zwz, zwi, zhdiv, zwzts, zwz_sav )
591      CALL wrk_alloc( jpi,jpj,jpk,kjpt+1,  ztrs )
[4990]592      !
[5770]593      IF( nn_timing == 1 )  CALL timing_stop('tra_adv_fct_zts')
[4990]594      !
[5770]595   END SUBROUTINE tra_adv_fct_zts
[4990]596
[5737]597
[2528]598   SUBROUTINE nonosc( pbef, paa, pbb, pcc, paft, p2dt )
[3]599      !!---------------------------------------------------------------------
600      !!                    ***  ROUTINE nonosc  ***
601      !!     
602      !! **  Purpose :   compute monotonic tracer fluxes from the upstream
603      !!       scheme and the before field by a nonoscillatory algorithm
604      !!
605      !! **  Method  :   ... ???
606      !!       warning : pbef and paft must be masked, but the boundaries
607      !!       conditions on the fluxes are not necessary zalezak (1979)
608      !!       drange (1995) multi-dimensional forward-in-time and upstream-
609      !!       in-space based differencing for fluid
610      !!----------------------------------------------------------------------
[6140]611      REAL(wp)                         , INTENT(in   ) ::   p2dt            ! tracer time-step
[2528]612      REAL(wp), DIMENSION (jpi,jpj,jpk), INTENT(in   ) ::   pbef, paft      ! before & after field
613      REAL(wp), DIMENSION (jpi,jpj,jpk), INTENT(inout) ::   paa, pbb, pcc   ! monotonic fluxes in the 3 directions
[2715]614      !
[4990]615      INTEGER  ::   ji, jj, jk   ! dummy loop indices
616      INTEGER  ::   ikm1         ! local integer
[6140]617      REAL(wp) ::   zpos, zneg, zbt, za, zb, zc, zbig, zrtrn    ! local scalars
[2715]618      REAL(wp) ::   zau, zbu, zcu, zav, zbv, zcv, zup, zdo            !   -      -
[3294]619      REAL(wp), POINTER, DIMENSION(:,:,:) :: zbetup, zbetdo, zbup, zbdo
[3]620      !!----------------------------------------------------------------------
[3294]621      !
622      IF( nn_timing == 1 )  CALL timing_start('nonosc')
623      !
624      CALL wrk_alloc( jpi, jpj, jpk, zbetup, zbetdo, zbup, zbdo )
625      !
[2715]626      zbig  = 1.e+40_wp
627      zrtrn = 1.e-15_wp
[4990]628      zbetup(:,:,:) = 0._wp   ;   zbetdo(:,:,:) = 0._wp
[785]629
[3]630      ! Search local extrema
631      ! --------------------
[785]632      ! max/min of pbef & paft with large negative/positive value (-/+zbig) inside land
[4990]633      zbup = MAX( pbef * tmask - zbig * ( 1._wp - tmask ),   &
634         &        paft * tmask - zbig * ( 1._wp - tmask )  )
635      zbdo = MIN( pbef * tmask + zbig * ( 1._wp - tmask ),   &
636         &        paft * tmask + zbig * ( 1._wp - tmask )  )
[785]637
[5120]638      DO jk = 1, jpkm1
639         ikm1 = MAX(jk-1,1)
640         DO jj = 2, jpjm1
641            DO ji = fs_2, fs_jpim1   ! vector opt.
642
[785]643               ! search maximum in neighbourhood
644               zup = MAX(  zbup(ji  ,jj  ,jk  ),   &
645                  &        zbup(ji-1,jj  ,jk  ), zbup(ji+1,jj  ,jk  ),   &
646                  &        zbup(ji  ,jj-1,jk  ), zbup(ji  ,jj+1,jk  ),   &
647                  &        zbup(ji  ,jj  ,ikm1), zbup(ji  ,jj  ,jk+1)  )
[3]648
[785]649               ! search minimum in neighbourhood
650               zdo = MIN(  zbdo(ji  ,jj  ,jk  ),   &
651                  &        zbdo(ji-1,jj  ,jk  ), zbdo(ji+1,jj  ,jk  ),   &
652                  &        zbdo(ji  ,jj-1,jk  ), zbdo(ji  ,jj+1,jk  ),   &
653                  &        zbdo(ji  ,jj  ,ikm1), zbdo(ji  ,jj  ,jk+1)  )
[3]654
[785]655               ! positive part of the flux
[3]656               zpos = MAX( 0., paa(ji-1,jj  ,jk  ) ) - MIN( 0., paa(ji  ,jj  ,jk  ) )   &
657                  & + MAX( 0., pbb(ji  ,jj-1,jk  ) ) - MIN( 0., pbb(ji  ,jj  ,jk  ) )   &
658                  & + MAX( 0., pcc(ji  ,jj  ,jk+1) ) - MIN( 0., pcc(ji  ,jj  ,jk  ) )
[785]659
660               ! negative part of the flux
[3]661               zneg = MAX( 0., paa(ji  ,jj  ,jk  ) ) - MIN( 0., paa(ji-1,jj  ,jk  ) )   &
662                  & + MAX( 0., pbb(ji  ,jj  ,jk  ) ) - MIN( 0., pbb(ji  ,jj-1,jk  ) )   &
663                  & + MAX( 0., pcc(ji  ,jj  ,jk  ) ) - MIN( 0., pcc(ji  ,jj  ,jk+1) )
[785]664
[3]665               ! up & down beta terms
[6140]666               zbt = e1e2t(ji,jj) * e3t_n(ji,jj,jk) / p2dt
[785]667               zbetup(ji,jj,jk) = ( zup            - paft(ji,jj,jk) ) / ( zpos + zrtrn ) * zbt
668               zbetdo(ji,jj,jk) = ( paft(ji,jj,jk) - zdo            ) / ( zneg + zrtrn ) * zbt
[3]669            END DO
670         END DO
671      END DO
[2528]672      CALL lbc_lnk( zbetup, 'T', 1. )   ;   CALL lbc_lnk( zbetdo, 'T', 1. )   ! lateral boundary cond. (unchanged sign)
[3]673
[237]674      ! 3. monotonic flux in the i & j direction (paa & pbb)
675      ! ----------------------------------------
[3]676      DO jk = 1, jpkm1
677         DO jj = 2, jpjm1
678            DO ji = fs_2, fs_jpim1   ! vector opt.
[4990]679               zau = MIN( 1._wp, zbetdo(ji,jj,jk), zbetup(ji+1,jj,jk) )
680               zbu = MIN( 1._wp, zbetup(ji,jj,jk), zbetdo(ji+1,jj,jk) )
[785]681               zcu =       ( 0.5  + SIGN( 0.5 , paa(ji,jj,jk) ) )
[4990]682               paa(ji,jj,jk) = paa(ji,jj,jk) * ( zcu * zau + ( 1._wp - zcu) * zbu )
[3]683
[4990]684               zav = MIN( 1._wp, zbetdo(ji,jj,jk), zbetup(ji,jj+1,jk) )
685               zbv = MIN( 1._wp, zbetup(ji,jj,jk), zbetdo(ji,jj+1,jk) )
[785]686               zcv =       ( 0.5  + SIGN( 0.5 , pbb(ji,jj,jk) ) )
[4990]687               pbb(ji,jj,jk) = pbb(ji,jj,jk) * ( zcv * zav + ( 1._wp - zcv) * zbv )
[3]688
689      ! monotonic flux in the k direction, i.e. pcc
690      ! -------------------------------------------
[785]691               za = MIN( 1., zbetdo(ji,jj,jk+1), zbetup(ji,jj,jk) )
692               zb = MIN( 1., zbetup(ji,jj,jk+1), zbetdo(ji,jj,jk) )
693               zc =       ( 0.5  + SIGN( 0.5 , pcc(ji,jj,jk+1) ) )
[4990]694               pcc(ji,jj,jk+1) = pcc(ji,jj,jk+1) * ( zc * za + ( 1._wp - zc) * zb )
[3]695            END DO
696         END DO
697      END DO
[2528]698      CALL lbc_lnk( paa, 'U', -1. )   ;   CALL lbc_lnk( pbb, 'V', -1. )   ! lateral boundary condition (changed sign)
[503]699      !
[3294]700      CALL wrk_dealloc( jpi, jpj, jpk, zbetup, zbetdo, zbup, zbdo )
[2715]701      !
[3294]702      IF( nn_timing == 1 )  CALL timing_stop('nonosc')
703      !
[3]704   END SUBROUTINE nonosc
705
[5770]706
707   SUBROUTINE interp_4th_cpt( pt_in, pt_out )
708      !!----------------------------------------------------------------------
709      !!                  ***  ROUTINE interp_4th_cpt  ***
710      !!
711      !! **  Purpose :   Compute the interpolation of tracer at w-point
712      !!
713      !! **  Method  :   4th order compact interpolation
714      !!----------------------------------------------------------------------
715      REAL(wp),DIMENSION(jpi,jpj,jpk), INTENT(in   ) ::   pt_in    ! now tracer fields
716      REAL(wp),DIMENSION(jpi,jpj,jpk), INTENT(  out) ::   pt_out   ! now tracer field interpolated at w-pts
717      !
718      INTEGER :: ji, jj, jk   ! dummy loop integers
719      REAL(wp),DIMENSION(jpi,jpj,jpk) :: zwd, zwi, zws, zwrm, zwt
720      !!----------------------------------------------------------------------
721     
722      DO jk = 3, jpkm1        !==  build the three diagonal matrix  ==!
723         DO jj = 1, jpj
724            DO ji = 1, jpi
725               zwd (ji,jj,jk) = 4._wp
726               zwi (ji,jj,jk) = 1._wp
727               zws (ji,jj,jk) = 1._wp
728               zwrm(ji,jj,jk) = 3._wp * ( pt_in(ji,jj,jk-1) + pt_in(ji,jj,jk) )
729               !
730               IF( tmask(ji,jj,jk+1) == 0._wp) THEN   ! Switch to second order centered at bottom
731                  zwd (ji,jj,jk) = 1._wp
732                  zwi (ji,jj,jk) = 0._wp
733                  zws (ji,jj,jk) = 0._wp
734                  zwrm(ji,jj,jk) = 0.5 * ( pt_in(ji,jj,jk-1) + pt_in(ji,jj,jk) )   
735               ENDIF
736            END DO
737         END DO
738      END DO
739      !
740      jk=2                                            ! Switch to second order centered at top
741      DO jj=1,jpj
742         DO ji=1,jpi
743            zwd (ji,jj,jk) = 1._wp
744            zwi (ji,jj,jk) = 0._wp
745            zws (ji,jj,jk) = 0._wp
746            zwrm(ji,jj,jk) = 0.5 * ( pt_in(ji,jj,jk-1) + pt_in(ji,jj,jk) )
747         END DO
748      END DO   
749      !
750      !                       !==  tridiagonal solve  ==!
751      DO jj = 1, jpj                ! first recurrence
752         DO ji = 1, jpi
753            zwt(ji,jj,2) = zwd(ji,jj,2)
754         END DO
755      END DO
756      DO jk = 3, jpkm1
757         DO jj = 1, jpj
758            DO ji = 1, jpi
759               zwt(ji,jj,jk) = zwd(ji,jj,jk) - zwi(ji,jj,jk) * zws(ji,jj,jk-1) /zwt(ji,jj,jk-1)
760            END DO
761         END DO
762      END DO
763      !
764      DO jj = 1, jpj                ! second recurrence:    Zk = Yk - Ik / Tk-1  Zk-1
765         DO ji = 1, jpi
766            pt_out(ji,jj,2) = zwrm(ji,jj,2)
767         END DO
768      END DO
769      DO jk = 3, jpkm1
770         DO jj = 1, jpj
771            DO ji = 1, jpi
772               pt_out(ji,jj,jk) = zwrm(ji,jj,jk) - zwi(ji,jj,jk) / zwt(ji,jj,jk-1) *pt_out(ji,jj,jk-1)             
773            END DO
774         END DO
775      END DO
776
777      DO jj = 1, jpj                ! third recurrence: Xk = (Zk - Sk Xk+1 ) / Tk
778         DO ji = 1, jpi
779            pt_out(ji,jj,jpkm1) = pt_out(ji,jj,jpkm1) / zwt(ji,jj,jpkm1)
780         END DO
781      END DO
782      DO jk = jpk-2, 2, -1
783         DO jj = 1, jpj
784            DO ji = 1, jpi
785               pt_out(ji,jj,jk) = ( pt_out(ji,jj,jk) - zws(ji,jj,jk) * pt_out(ji,jj,jk+1) ) / zwt(ji,jj,jk)
786            END DO
787         END DO
788      END DO
789      !   
790   END SUBROUTINE interp_4th_cpt
791   
[3]792   !!======================================================================
[5770]793END MODULE traadv_fct
Note: See TracBrowser for help on using the repository browser.