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/2019/dev_r11470_HPC_12_mpi3/src/OCE/TRA – NEMO

source: NEMO/branches/2019/dev_r11470_HPC_12_mpi3/src/OCE/TRA/traadv_fct.F90 @ 11940

Last change on this file since 11940 was 11940, checked in by mocavero, 4 years ago

Add MPI3 neighbourhood collectives halo exchange in LBC and call it in tracer advection FCT scheme #2011

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