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/2020/dev_r13648_ASINTER-04_laurent_bulk_ice/src/OCE/TRA – NEMO

source: NEMO/branches/2020/dev_r13648_ASINTER-04_laurent_bulk_ice/src/OCE/TRA/traadv_fct.F90 @ 14017

Last change on this file since 14017 was 14017, checked in by laurent, 4 years ago

Keep up with trunk revision 13999

  • Property svn:keywords set to Id
File size: 36.5 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
[14017]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
[14017]26   USE iom            !
[3625]27   USE lib_mpp        ! MPP library
[14017]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
[14017]36   PUBLIC   tridia_solver      ! called by traadv_fct_lf.F90
37   PUBLIC   nonosc             ! called by traadv_fct_lf.F90 - key_agrif
[3]38
[5770]39   LOGICAL  ::   l_trd   ! flag to compute trends
[7646]40   LOGICAL  ::   l_ptr   ! flag to compute poleward transport
41   LOGICAL  ::   l_hst   ! flag to compute heat/salt transport
[5770]42   REAL(wp) ::   r1_6 = 1._wp / 6._wp   ! =1/6
[2528]43
[7646]44   !                                        ! tridiag solver associated indices:
45   INTEGER, PARAMETER ::   np_NH   = 0   ! Neumann homogeneous boundary condition
46   INTEGER, PARAMETER ::   np_CEN2 = 1   ! 2nd order centered  boundary condition
47
[3]48   !! * Substitutions
[12377]49#  include "do_loop_substitute.h90"
[13237]50#  include "domzgr_substitute.h90"
[3]51   !!----------------------------------------------------------------------
[9598]52   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
[1152]53   !! $Id$
[10068]54   !! Software governed by the CeCILL license (see ./LICENSE)
[3]55   !!----------------------------------------------------------------------
56CONTAINS
57
[12377]58   SUBROUTINE tra_adv_fct( kt, kit000, cdtype, p2dt, pU, pV, pW,       &
59      &                    Kbb, Kmm, pt, kjpt, Krhs, kn_fct_h, kn_fct_v )
[3]60      !!----------------------------------------------------------------------
[5770]61      !!                  ***  ROUTINE tra_adv_fct  ***
[14017]62      !!
[6140]63      !! **  Purpose :   Compute the now trend due to total advection of tracers
64      !!               and add it to the general trend of tracer equations
[3]65      !!
[5770]66      !! **  Method  : - 2nd or 4th FCT scheme on the horizontal direction
67      !!               (choice through the value of kn_fct)
[14017]68      !!               - on the vertical the 4th order is a compact scheme
69      !!               - corrected flux (monotonic correction)
[3]70      !!
[12377]71      !! ** Action : - update pt(:,:,:,:,Krhs)  with the now advective tracer trends
[9019]72      !!             - send trends to trdtra module for further diagnostics (l_trdtra=T)
[12377]73      !!             - poleward advective heat and salt transport (ln_diaptr=T)
[503]74      !!----------------------------------------------------------------------
[12377]75      INTEGER                                  , INTENT(in   ) ::   kt              ! ocean time-step index
76      INTEGER                                  , INTENT(in   ) ::   Kbb, Kmm, Krhs  ! ocean time level indices
77      INTEGER                                  , INTENT(in   ) ::   kit000          ! first time step index
78      CHARACTER(len=3)                         , INTENT(in   ) ::   cdtype          ! =TRA or TRC (tracer indicator)
79      INTEGER                                  , INTENT(in   ) ::   kjpt            ! number of tracers
80      INTEGER                                  , INTENT(in   ) ::   kn_fct_h        ! order of the FCT scheme (=2 or 4)
81      INTEGER                                  , INTENT(in   ) ::   kn_fct_v        ! order of the FCT scheme (=2 or 4)
82      REAL(wp)                                 , INTENT(in   ) ::   p2dt            ! tracer time-step
[14017]83      ! TEMP: [tiling] This can be A2D(nn_hls) if using XIOS (subdomain support)
84      ! NOTE: [tiling-comms-merge] These were changed to INTENT(inout) but they are not modified, so it is reverted
[12377]85      REAL(wp), DIMENSION(jpi,jpj,jpk         ), INTENT(in   ) ::   pU, pV, pW      ! 3 ocean volume flux components
86      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt,jpt), INTENT(inout) ::   pt              ! tracers and RHS of tracer equation
[2715]87      !
[14017]88      INTEGER  ::   ji, jj, jk, jn                           ! dummy loop indices
[6140]89      REAL(wp) ::   ztra                                     ! local scalar
[5770]90      REAL(wp) ::   zfp_ui, zfp_vj, zfp_wk, zC2t_u, zC4t_u   !   -      -
91      REAL(wp) ::   zfm_ui, zfm_vj, zfm_wk, zC2t_v, zC4t_v   !   -      -
[14017]92      REAL(wp), DIMENSION(A2D(nn_hls),jpk)        ::   zwi, zwx, zwy, zwz, ztu, ztv, zltu, zltv, ztw
[9019]93      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::   ztrdx, ztrdy, ztrdz, zptry
[11407]94      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::   zwinf, zwdia, zwsup
95      LOGICAL  ::   ll_zAimp                                 ! flag to apply adaptive implicit vertical advection
[3]96      !!----------------------------------------------------------------------
[3294]97      !
[14017]98      IF( ntile == 0 .OR. ntile == 1 )  THEN                       ! Do only on the first tile
99         IF( kt == kit000 )  THEN
100            IF(lwp) WRITE(numout,*)
101            IF(lwp) WRITE(numout,*) 'tra_adv_fct : FCT advection scheme on ', cdtype
102            IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
103         ENDIF
104         ! NOTE: [tiling-comms-merge] Bug fix- move array zeroing out of this IF block
105         !
106         l_trd = .FALSE.            ! set local switches
107         l_hst = .FALSE.
108         l_ptr = .FALSE.
109         ll_zAimp = .FALSE.
110         IF( ( cdtype == 'TRA' .AND. l_trdtra  ) .OR. ( cdtype =='TRC' .AND. l_trdtrc ) )      l_trd = .TRUE.
111         IF(   cdtype == 'TRA' .AND. ( iom_use( 'sophtadv' ) .OR. iom_use( 'sophtadv' ) ) )    l_ptr = .TRUE.
112         IF(   cdtype == 'TRA' .AND. ( iom_use("uadv_heattr") .OR. iom_use("vadv_heattr") .OR.  &
113            &                         iom_use("uadv_salttr") .OR. iom_use("vadv_salttr")  ) )  l_hst = .TRUE.
114         !
[3]115      ENDIF
[14017]116
[13226]117      !! -- init to 0
118      zwi(:,:,:) = 0._wp
119      zwx(:,:,:) = 0._wp
120      zwy(:,:,:) = 0._wp
121      zwz(:,:,:) = 0._wp
122      ztu(:,:,:) = 0._wp
123      ztv(:,:,:) = 0._wp
124      zltu(:,:,:) = 0._wp
125      zltv(:,:,:) = 0._wp
126      ztw(:,:,:) = 0._wp
[2528]127      !
[7646]128      IF( l_trd .OR. l_hst )  THEN
[14017]129         ALLOCATE( ztrdx(A2D(nn_hls),jpk), ztrdy(A2D(nn_hls),jpk), ztrdz(A2D(nn_hls),jpk) )
[7753]130         ztrdx(:,:,:) = 0._wp   ;    ztrdy(:,:,:) = 0._wp   ;   ztrdz(:,:,:) = 0._wp
[3294]131      ENDIF
[2528]132      !
[14017]133      IF( l_ptr ) THEN
134         ALLOCATE( zptry(A2D(nn_hls),jpk) )
[7753]135         zptry(:,:,:) = 0._wp
[7646]136      ENDIF
[2528]137      !
[11407]138      ! If adaptive vertical advection, check if it is needed on this PE at this time
139      IF( ln_zad_Aimp ) THEN
[14017]140         IF( MAXVAL( ABS( wi(A2D(nn_hls),:) ) ) > 0._wp ) ll_zAimp = .TRUE.
[11407]141      END IF
142      ! If active adaptive vertical advection, build tridiagonal matrix
143      IF( ll_zAimp ) THEN
[14017]144         ALLOCATE(zwdia(A2D(nn_hls),jpk), zwinf(A2D(nn_hls),jpk), zwsup(A2D(nn_hls),jpk))
145         DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 1, jpkm1 )
[13237]146            zwdia(ji,jj,jk) =  1._wp + p2dt * ( MAX( wi(ji,jj,jk) , 0._wp ) - MIN( wi(ji,jj,jk+1) , 0._wp ) )   &
147            &                               / e3t(ji,jj,jk,Krhs)
[12377]148            zwinf(ji,jj,jk) =  p2dt * MIN( wi(ji,jj,jk  ) , 0._wp ) / e3t(ji,jj,jk,Krhs)
149            zwsup(ji,jj,jk) = -p2dt * MAX( wi(ji,jj,jk+1) , 0._wp ) / e3t(ji,jj,jk,Krhs)
150         END_3D
[11407]151      END IF
152      !
[6140]153      DO jn = 1, kjpt            !==  loop over the tracers  ==!
[5770]154         !
155         !        !==  upstream advection with initial mass fluxes & intermediate update  ==!
[14017]156         !                    !* upstream tracer flux in the i and j direction
157         DO_3D( nn_hls, nn_hls-1, nn_hls, nn_hls-1, 1, jpkm1 )
[12377]158            ! upstream scheme
159            zfp_ui = pU(ji,jj,jk) + ABS( pU(ji,jj,jk) )
160            zfm_ui = pU(ji,jj,jk) - ABS( pU(ji,jj,jk) )
161            zfp_vj = pV(ji,jj,jk) + ABS( pV(ji,jj,jk) )
162            zfm_vj = pV(ji,jj,jk) - ABS( pV(ji,jj,jk) )
163            zwx(ji,jj,jk) = 0.5 * ( zfp_ui * pt(ji,jj,jk,jn,Kbb) + zfm_ui * pt(ji+1,jj  ,jk,jn,Kbb) )
164            zwy(ji,jj,jk) = 0.5 * ( zfp_vj * pt(ji,jj,jk,jn,Kbb) + zfm_vj * pt(ji  ,jj+1,jk,jn,Kbb) )
165         END_3D
[13497]166         !                               !* upstream tracer flux in the k direction *!
167         DO_3D( 1, 1, 1, 1, 2, jpkm1 )      ! Interior value ( multiplied by wmask)
[12377]168            zfp_wk = pW(ji,jj,jk) + ABS( pW(ji,jj,jk) )
169            zfm_wk = pW(ji,jj,jk) - ABS( pW(ji,jj,jk) )
170            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)
171         END_3D
[13497]172         IF( ln_linssh ) THEN               ! top ocean value (only in linear free surface as zwz has been w-masked)
173            IF( ln_isfcav ) THEN                        ! top of the ice-shelf cavities and at the ocean surface
[13295]174               DO_2D( 1, 1, 1, 1 )
[14017]175                  zwz(ji,jj, mikt(ji,jj) ) = pW(ji,jj,mikt(ji,jj)) * pt(ji,jj,mikt(ji,jj),jn,Kbb)   ! linear free surface
[12377]176               END_2D
[13497]177            ELSE                                        ! no cavities: only at the ocean surface
[13295]178               DO_2D( 1, 1, 1, 1 )
[13286]179                  zwz(ji,jj,1) = pW(ji,jj,1) * pt(ji,jj,1,jn,Kbb)
180               END_2D
[5770]181            ENDIF
[5120]182         ENDIF
[14017]183         !
184         DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 1, jpkm1 )   !* trend and after field with monotonic scheme
[13497]185            !                               ! total intermediate advective trends
[12377]186            ztra = - (  zwx(ji,jj,jk) - zwx(ji-1,jj  ,jk  )   &
187               &      + zwy(ji,jj,jk) - zwy(ji  ,jj-1,jk  )   &
188               &      + zwz(ji,jj,jk) - zwz(ji  ,jj  ,jk+1) ) * r1_e1e2t(ji,jj)
[13497]189            !                               ! update and guess with monotonic sheme
[13237]190            pt(ji,jj,jk,jn,Krhs) =                   pt(ji,jj,jk,jn,Krhs) +       ztra   &
191               &                                  / e3t(ji,jj,jk,Kmm ) * tmask(ji,jj,jk)
192            zwi(ji,jj,jk)    = ( e3t(ji,jj,jk,Kbb) * pt(ji,jj,jk,jn,Kbb) + p2dt * ztra ) &
193               &                                  / e3t(ji,jj,jk,Krhs) * tmask(ji,jj,jk)
[12377]194         END_3D
[14017]195
[11407]196         IF ( ll_zAimp ) THEN
197            CALL tridia_solver( zwdia, zwsup, zwinf, zwi, zwi , 0 )
198            !
[11411]199            ztw(:,:,1) = 0._wp ; ztw(:,:,jpk) = 0._wp ;
[14017]200            DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )       ! Interior value ( multiplied by wmask)
[12377]201               zfp_wk = wi(ji,jj,jk) + ABS( wi(ji,jj,jk) )
202               zfm_wk = wi(ji,jj,jk) - ABS( wi(ji,jj,jk) )
203               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)
204               zwz(ji,jj,jk) = zwz(ji,jj,jk) + ztw(ji,jj,jk) ! update vertical fluxes
205            END_3D
[13295]206            DO_3D( 0, 0, 0, 0, 1, jpkm1 )
[12377]207               pt(ji,jj,jk,jn,Krhs) = pt(ji,jj,jk,jn,Krhs) - ( ztw(ji,jj,jk) - ztw(ji  ,jj  ,jk+1) ) &
208                  &                                        * r1_e1e2t(ji,jj) / e3t(ji,jj,jk,Kmm)
209            END_3D
[11407]210            !
211         END IF
[14017]212         !
[7646]213         IF( l_trd .OR. l_hst )  THEN             ! trend diagnostics (contribution of upstream fluxes)
[7753]214            ztrdx(:,:,:) = zwx(:,:,:)   ;   ztrdy(:,:,:) = zwy(:,:,:)   ;   ztrdz(:,:,:) = zwz(:,:,:)
[2528]215         END IF
[5770]216         !                             ! "Poleward" heat and salt transports (contribution of upstream fluxes)
[14017]217         IF( l_ptr )   zptry(:,:,:) = zwy(:,:,:)
[5770]218         !
219         !        !==  anti-diffusive flux : high order minus low order  ==!
220         !
[6140]221         SELECT CASE( kn_fct_h )    !* horizontal anti-diffusive fluxes
[5770]222         !
[6140]223         CASE(  2  )                   !- 2nd order centered
[14017]224            DO_3D( nn_hls, nn_hls-1, nn_hls, nn_hls-1, 1, jpkm1 )
[12377]225               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)
226               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)
227            END_3D
[5770]228            !
[6140]229         CASE(  4  )                   !- 4th order centered
[7753]230            zltu(:,:,jpk) = 0._wp            ! Bottom value : flux set to zero
231            zltv(:,:,jpk) = 0._wp
[6140]232            DO jk = 1, jpkm1                 ! Laplacian
[13497]233               DO_2D( 1, 0, 1, 0 )                 ! 1st derivative (gradient)
[12377]234                  ztu(ji,jj,jk) = ( pt(ji+1,jj  ,jk,jn,Kmm) - pt(ji,jj,jk,jn,Kmm) ) * umask(ji,jj,jk)
235                  ztv(ji,jj,jk) = ( pt(ji  ,jj+1,jk,jn,Kmm) - pt(ji,jj,jk,jn,Kmm) ) * vmask(ji,jj,jk)
236               END_2D
[13497]237               DO_2D( 0, 0, 0, 0 )                 ! 2nd derivative * 1/ 6
[12377]238                  zltu(ji,jj,jk) = (  ztu(ji,jj,jk) + ztu(ji-1,jj,jk)  ) * r1_6
239                  zltv(ji,jj,jk) = (  ztv(ji,jj,jk) + ztv(ji,jj-1,jk)  ) * r1_6
240               END_2D
[503]241            END DO
[13226]242            CALL lbc_lnk_multi( 'traadv_fct', zltu, 'T', 1.0_wp , zltv, 'T', 1.0_wp )   ! Lateral boundary cond. (unchanged sgn)
[5770]243            !
[14017]244            DO_3D( 1, 0, 1, 0, 1, jpkm1 )
[12377]245               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
246               zC2t_v = pt(ji,jj,jk,jn,Kmm) + pt(ji  ,jj+1,jk,jn,Kmm)
[14017]247               !                                                        ! C4 minus upstream advective fluxes
[12377]248               zwx(ji,jj,jk) =  0.5_wp * pU(ji,jj,jk) * ( zC2t_u + zltu(ji,jj,jk) - zltu(ji+1,jj,jk) ) - zwx(ji,jj,jk)
249               zwy(ji,jj,jk) =  0.5_wp * pV(ji,jj,jk) * ( zC2t_v + zltv(ji,jj,jk) - zltv(ji,jj+1,jk) ) - zwy(ji,jj,jk)
250            END_3D
[14017]251            IF (nn_hls.EQ.2) CALL lbc_lnk_multi( 'traadv_fct', zwx, 'U', -1.0_wp, zwy, 'V', -1.0_wp )   ! Lateral boundary cond. (unchanged sgn)
[5770]252            !
[6140]253         CASE(  41 )                   !- 4th order centered       ==>>   !!gm coding attempt   need to be tested
[7753]254            ztu(:,:,jpk) = 0._wp             ! Bottom value : flux set to zero
255            ztv(:,:,jpk) = 0._wp
[14017]256            DO_3D( nn_hls, nn_hls-1, nn_hls, nn_hls-1, 1, jpkm1 )    ! 1st derivative (gradient)
[12377]257               ztu(ji,jj,jk) = ( pt(ji+1,jj  ,jk,jn,Kmm) - pt(ji,jj,jk,jn,Kmm) ) * umask(ji,jj,jk)
258               ztv(ji,jj,jk) = ( pt(ji  ,jj+1,jk,jn,Kmm) - pt(ji,jj,jk,jn,Kmm) ) * vmask(ji,jj,jk)
259            END_3D
[14017]260            IF (nn_hls.EQ.1) CALL lbc_lnk_multi( 'traadv_fct', ztu, 'U', -1.0_wp , ztv, 'V', -1.0_wp )   ! Lateral boundary cond. (unchanged sgn)
[5770]261            !
[14017]262            IF (nn_hls.EQ.1) CALL lbc_lnk_multi( 'traadv_fct', ztu, 'U', -1.0_wp , ztv, 'V', -1.0_wp )   ! Lateral boundary cond. (unchanged sgn)
263            !
[13497]264            DO_3D( 0, 0, 0, 0, 1, jpkm1 )    ! Horizontal advective fluxes
[12377]265               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)
266               zC2t_v = pt(ji,jj,jk,jn,Kmm) + pt(ji  ,jj+1,jk,jn,Kmm)
267               !                                                  ! C4 interpolation of T at u- & v-points (x2)
268               zC4t_u =  zC2t_u + r1_6 * ( ztu(ji-1,jj  ,jk) - ztu(ji+1,jj  ,jk) )
269               zC4t_v =  zC2t_v + r1_6 * ( ztv(ji  ,jj-1,jk) - ztv(ji  ,jj+1,jk) )
[14017]270               !                                                  ! C4 minus upstream advective fluxes
[12377]271               zwx(ji,jj,jk) =  0.5_wp * pU(ji,jj,jk) * zC4t_u - zwx(ji,jj,jk)
272               zwy(ji,jj,jk) =  0.5_wp * pV(ji,jj,jk) * zC4t_v - zwy(ji,jj,jk)
273            END_3D
[14017]274            IF (nn_hls.EQ.2) CALL lbc_lnk_multi( 'traadv_fct', zwx, 'U', -1.0_wp , zwy, 'V', -1.0_wp )   ! Lateral boundary cond. (unchanged sgn)
[5770]275            !
276         END SELECT
[14017]277         !
[6140]278         SELECT CASE( kn_fct_v )    !* vertical anti-diffusive fluxes (w-masked interior values)
[5770]279         !
[6140]280         CASE(  2  )                   !- 2nd order centered
[14017]281            DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )
[12377]282               zwz(ji,jj,jk) =  (  pW(ji,jj,jk) * 0.5_wp * ( pt(ji,jj,jk,jn,Kmm) + pt(ji,jj,jk-1,jn,Kmm) )   &
283                  &              - zwz(ji,jj,jk)  ) * wmask(ji,jj,jk)
284            END_3D
[5770]285            !
[6140]286         CASE(  4  )                   !- 4th order COMPACT
[12377]287            CALL interp_4th_cpt( pt(:,:,:,jn,Kmm) , ztw )   ! zwt = COMPACT interpolation of T at w-point
[14017]288            DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )
[12377]289               zwz(ji,jj,jk) = ( pW(ji,jj,jk) * ztw(ji,jj,jk) - zwz(ji,jj,jk) ) * wmask(ji,jj,jk)
290            END_3D
[5770]291            !
292         END SELECT
[6140]293         IF( ln_linssh ) THEN    ! top ocean value: high order = upstream  ==>>  zwz=0
[7753]294            zwz(:,:,1) = 0._wp   ! only ocean surface as interior zwz values have been w-masked
[6140]295         ENDIF
[14017]296         !
297         IF (nn_hls.EQ.1) THEN
298            CALL lbc_lnk_multi( 'traadv_fct', zwi, 'T', 1.0_wp, zwx, 'U', -1.0_wp , zwy, 'V', -1.0_wp, zwz, 'T', 1.0_wp )
299         ELSE
300            CALL lbc_lnk( 'traadv_fct', zwi, 'T', 1.0_wp)
301         END IF
302         !
303         IF (nn_hls.EQ.1) THEN
304            CALL lbc_lnk_multi( 'traadv_fct', zwi, 'T', 1.0_wp, zwx, 'U', -1.0_wp , zwy, 'V', -1.0_wp, zwz, 'T', 1.0_wp )
305         ELSE
306            CALL lbc_lnk( 'traadv_fct', zwi, 'T', 1.0_wp)
307         END IF
308         !
[11407]309         IF ( ll_zAimp ) THEN
[14017]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)
[14017]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            !
[14017]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) )
[14017]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         !
[12377]329         CALL nonosc( Kmm, 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
[14017]354         END IF
355         ! NOTE: [tiling-comms-merge] I tested this
356         ! NOT TESTED - NEED l_trd OR l_hst TRUE
[9019]357         IF( l_trd .OR. l_hst ) THEN   ! trend diagnostics // heat/salt transport
[14017]358            ztrdx(:,:,:) = ztrdx(:,:,:) + zwx(:,:,:)  ! <<< add anti-diffusive fluxes
[9019]359            ztrdy(:,:,:) = ztrdy(:,:,:) + zwy(:,:,:)  !     to upstream fluxes
360            ztrdz(:,:,:) = ztrdz(:,:,:) + zwz(:,:,:)  !
[5770]361            !
[9019]362            IF( l_trd ) THEN              ! trend diagnostics
[12377]363               CALL trd_tra( kt, Kmm, Krhs, cdtype, jn, jptra_xad, ztrdx, pU, pt(:,:,:,jn,Kmm) )
364               CALL trd_tra( kt, Kmm, Krhs, cdtype, jn, jptra_yad, ztrdy, pV, pt(:,:,:,jn,Kmm) )
365               CALL trd_tra( kt, Kmm, Krhs, cdtype, jn, jptra_zad, ztrdz, pW, pt(:,:,:,jn,Kmm) )
[9019]366            ENDIF
367            !                             ! heat/salt transport
368            IF( l_hst )   CALL dia_ar5_hst( jn, 'adv', ztrdx(:,:,:), ztrdy(:,:,:) )
[5770]369            !
[9019]370         ENDIF
[14017]371         ! NOTE: [tiling-comms-merge] I tested this
372         ! NOT TESTED - NEED l_ptr TRUE
[9019]373         IF( l_ptr ) THEN              ! "Poleward" transports
374            zptry(:,:,:) = zptry(:,:,:) + zwy(:,:,:)  ! <<< add anti-diffusive fluxes
[7646]375            CALL dia_ptr_hst( jn, 'adv', zptry(:,:,:) )
[2528]376         ENDIF
[503]377         !
[6140]378      END DO                     ! end of tracer loop
[503]379      !
[11407]380      IF ( ll_zAimp ) THEN
381         DEALLOCATE( zwdia, zwinf, zwsup )
382      ENDIF
[14017]383      IF( l_trd .OR. l_hst ) THEN
[10024]384         DEALLOCATE( ztrdx, ztrdy, ztrdz )
385      ENDIF
[14017]386      IF( l_ptr ) THEN
[10024]387         DEALLOCATE( zptry )
388      ENDIF
389      !
[5770]390   END SUBROUTINE tra_adv_fct
[3]391
[5737]392
[12377]393   SUBROUTINE nonosc( Kmm, pbef, paa, pbb, pcc, paft, p2dt )
[3]394      !!---------------------------------------------------------------------
395      !!                    ***  ROUTINE nonosc  ***
396      !!
[14017]397      !! **  Purpose :   compute monotonic tracer fluxes from the upstream
398      !!       scheme and the before field by a nonoscillatory algorithm
399      !!
[3]400      !! **  Method  :   ... ???
401      !!       warning : pbef and paft must be masked, but the boundaries
402      !!       conditions on the fluxes are not necessary zalezak (1979)
403      !!       drange (1995) multi-dimensional forward-in-time and upstream-
404      !!       in-space based differencing for fluid
405      !!----------------------------------------------------------------------
[14017]406      INTEGER                         , INTENT(in   ) ::   Kmm             ! time level index
407      REAL(wp)                        , INTENT(in   ) ::   p2dt            ! tracer time-step
408      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in   ) ::   pbef            ! before field
409      REAL(wp), DIMENSION(A2D(nn_hls)    ,jpk), INTENT(in   ) ::   paft            ! after field
410      REAL(wp), DIMENSION(A2D(nn_hls)    ,jpk), INTENT(inout) ::   paa, pbb, pcc   ! monotonic fluxes in the 3 directions
[2715]411      !
[4990]412      INTEGER  ::   ji, jj, jk   ! dummy loop indices
413      INTEGER  ::   ikm1         ! local integer
[13226]414      REAL(dp) ::   zpos, zneg, zbt, za, zb, zc, zbig, zrtrn    ! local scalars
415      REAL(dp) ::   zau, zbu, zcu, zav, zbv, zcv, zup, zdo            !   -      -
[14017]416      REAL(dp), DIMENSION(A2D(nn_hls),jpk) :: zbetup, zbetdo, zbup, zbdo
[3]417      !!----------------------------------------------------------------------
[3294]418      !
[13226]419      zbig  = 1.e+40_dp
420      zrtrn = 1.e-15_dp
421      zbetup(:,:,:) = 0._dp   ;   zbetdo(:,:,:) = 0._dp
[785]422
[3]423      ! Search local extrema
424      ! --------------------
[785]425      ! max/min of pbef & paft with large negative/positive value (-/+zbig) inside land
[14017]426      DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, 1, jpk )
427         zbup(ji,jj,jk) = MAX( pbef(ji,jj,jk) * tmask(ji,jj,jk) - zbig * ( 1._wp - tmask(ji,jj,jk) ),   &
428            &                  paft(ji,jj,jk) * tmask(ji,jj,jk) - zbig * ( 1._wp - tmask(ji,jj,jk) )  )
429         zbdo(ji,jj,jk) = MIN( pbef(ji,jj,jk) * tmask(ji,jj,jk) + zbig * ( 1._wp - tmask(ji,jj,jk) ),   &
430            &                  paft(ji,jj,jk) * tmask(ji,jj,jk) + zbig * ( 1._wp - tmask(ji,jj,jk) )  )
431      END_3D
[785]432
[5120]433      DO jk = 1, jpkm1
434         ikm1 = MAX(jk-1,1)
[14017]435         DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )
[5120]436
[12377]437            ! search maximum in neighbourhood
438            zup = MAX(  zbup(ji  ,jj  ,jk  ),   &
439               &        zbup(ji-1,jj  ,jk  ), zbup(ji+1,jj  ,jk  ),   &
440               &        zbup(ji  ,jj-1,jk  ), zbup(ji  ,jj+1,jk  ),   &
441               &        zbup(ji  ,jj  ,ikm1), zbup(ji  ,jj  ,jk+1)  )
[3]442
[12377]443            ! search minimum in neighbourhood
444            zdo = MIN(  zbdo(ji  ,jj  ,jk  ),   &
445               &        zbdo(ji-1,jj  ,jk  ), zbdo(ji+1,jj  ,jk  ),   &
446               &        zbdo(ji  ,jj-1,jk  ), zbdo(ji  ,jj+1,jk  ),   &
447               &        zbdo(ji  ,jj  ,ikm1), zbdo(ji  ,jj  ,jk+1)  )
[3]448
[12377]449            ! positive part of the flux
450            zpos = MAX( 0., paa(ji-1,jj  ,jk  ) ) - MIN( 0., paa(ji  ,jj  ,jk  ) )   &
451               & + MAX( 0., pbb(ji  ,jj-1,jk  ) ) - MIN( 0., pbb(ji  ,jj  ,jk  ) )   &
452               & + MAX( 0., pcc(ji  ,jj  ,jk+1) ) - MIN( 0., pcc(ji  ,jj  ,jk  ) )
[785]453
[12377]454            ! negative part of the flux
455            zneg = MAX( 0., paa(ji  ,jj  ,jk  ) ) - MIN( 0., paa(ji-1,jj  ,jk  ) )   &
456               & + MAX( 0., pbb(ji  ,jj  ,jk  ) ) - MIN( 0., pbb(ji  ,jj-1,jk  ) )   &
457               & + MAX( 0., pcc(ji  ,jj  ,jk  ) ) - MIN( 0., pcc(ji  ,jj  ,jk+1) )
[785]458
[12377]459            ! up & down beta terms
460            zbt = e1e2t(ji,jj) * e3t(ji,jj,jk,Kmm) / p2dt
461            zbetup(ji,jj,jk) = ( zup            - paft(ji,jj,jk) ) / ( zpos + zrtrn ) * zbt
462            zbetdo(ji,jj,jk) = ( paft(ji,jj,jk) - zdo            ) / ( zneg + zrtrn ) * zbt
463         END_2D
[3]464      END DO
[14017]465      IF (nn_hls.EQ.1) CALL lbc_lnk_multi( 'traadv_fct', zbetup, 'T', 1.0_wp , zbetdo, 'T', 1.0_wp )   ! lateral boundary cond. (unchanged sign)
[3]466
[237]467      ! 3. monotonic flux in the i & j direction (paa & pbb)
468      ! ----------------------------------------
[14017]469      DO_3D( 1, 0, 1, 0, 1, jpkm1 )
[12377]470         zau = MIN( 1._wp, zbetdo(ji,jj,jk), zbetup(ji+1,jj,jk) )
471         zbu = MIN( 1._wp, zbetup(ji,jj,jk), zbetdo(ji+1,jj,jk) )
[13226]472         zcu =       ( 0.5  + SIGN( 0.5_wp , paa(ji,jj,jk) ) )
[12377]473         paa(ji,jj,jk) = paa(ji,jj,jk) * ( zcu * zau + ( 1._wp - zcu) * zbu )
[3]474
[12377]475         zav = MIN( 1._wp, zbetdo(ji,jj,jk), zbetup(ji,jj+1,jk) )
476         zbv = MIN( 1._wp, zbetup(ji,jj,jk), zbetdo(ji,jj+1,jk) )
[13226]477         zcv =       ( 0.5  + SIGN( 0.5_wp , pbb(ji,jj,jk) ) )
[12377]478         pbb(ji,jj,jk) = pbb(ji,jj,jk) * ( zcv * zav + ( 1._wp - zcv) * zbv )
[3]479
[13497]480      ! monotonic flux in the k direction, i.e. pcc
481      ! -------------------------------------------
[12377]482         za = MIN( 1., zbetdo(ji,jj,jk+1), zbetup(ji,jj,jk) )
483         zb = MIN( 1., zbetup(ji,jj,jk+1), zbetdo(ji,jj,jk) )
[13226]484         zc =       ( 0.5  + SIGN( 0.5_wp , pcc(ji,jj,jk+1) ) )
[12377]485         pcc(ji,jj,jk+1) = pcc(ji,jj,jk+1) * ( zc * za + ( 1._wp - zc) * zb )
486      END_3D
[503]487      !
[3]488   END SUBROUTINE nonosc
489
[5770]490
[7646]491   SUBROUTINE interp_4th_cpt_org( pt_in, pt_out )
[5770]492      !!----------------------------------------------------------------------
[7646]493      !!                  ***  ROUTINE interp_4th_cpt_org  ***
[14017]494      !!
[5770]495      !! **  Purpose :   Compute the interpolation of tracer at w-point
496      !!
497      !! **  Method  :   4th order compact interpolation
498      !!----------------------------------------------------------------------
499      REAL(wp),DIMENSION(jpi,jpj,jpk), INTENT(in   ) ::   pt_in    ! now tracer fields
500      REAL(wp),DIMENSION(jpi,jpj,jpk), INTENT(  out) ::   pt_out   ! now tracer field interpolated at w-pts
501      !
502      INTEGER :: ji, jj, jk   ! dummy loop integers
503      REAL(wp),DIMENSION(jpi,jpj,jpk) :: zwd, zwi, zws, zwrm, zwt
504      !!----------------------------------------------------------------------
[14017]505
[13497]506      DO_3D( 1, 1, 1, 1, 3, jpkm1 )       !==  build the three diagonal matrix  ==!
[12377]507         zwd (ji,jj,jk) = 4._wp
508         zwi (ji,jj,jk) = 1._wp
509         zws (ji,jj,jk) = 1._wp
510         zwrm(ji,jj,jk) = 3._wp * ( pt_in(ji,jj,jk-1) + pt_in(ji,jj,jk) )
511         !
512         IF( tmask(ji,jj,jk+1) == 0._wp) THEN   ! Switch to second order centered at bottom
[5770]513            zwd (ji,jj,jk) = 1._wp
514            zwi (ji,jj,jk) = 0._wp
515            zws (ji,jj,jk) = 0._wp
[14017]516            zwrm(ji,jj,jk) = 0.5 * ( pt_in(ji,jj,jk-1) + pt_in(ji,jj,jk) )
[12377]517         ENDIF
518      END_3D
[5770]519      !
[13497]520      jk = 2                                    ! Switch to second order centered at top
[13295]521      DO_2D( 1, 1, 1, 1 )
[12377]522         zwd (ji,jj,jk) = 1._wp
523         zwi (ji,jj,jk) = 0._wp
524         zws (ji,jj,jk) = 0._wp
525         zwrm(ji,jj,jk) = 0.5 * ( pt_in(ji,jj,jk-1) + pt_in(ji,jj,jk) )
526      END_2D
527      !
[5770]528      !                       !==  tridiagonal solve  ==!
[13497]529      DO_2D( 1, 1, 1, 1 )           ! first recurrence
[12377]530         zwt(ji,jj,2) = zwd(ji,jj,2)
531      END_2D
[13295]532      DO_3D( 1, 1, 1, 1, 3, jpkm1 )
[12377]533         zwt(ji,jj,jk) = zwd(ji,jj,jk) - zwi(ji,jj,jk) * zws(ji,jj,jk-1) /zwt(ji,jj,jk-1)
534      END_3D
[5770]535      !
[13497]536      DO_2D( 1, 1, 1, 1 )           ! second recurrence:    Zk = Yk - Ik / Tk-1  Zk-1
[12377]537         pt_out(ji,jj,2) = zwrm(ji,jj,2)
538      END_2D
[13295]539      DO_3D( 1, 1, 1, 1, 3, jpkm1 )
[14017]540         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]541      END_3D
[5770]542
[13497]543      DO_2D( 1, 1, 1, 1 )           ! third recurrence: Xk = (Zk - Sk Xk+1 ) / Tk
[12377]544         pt_out(ji,jj,jpkm1) = pt_out(ji,jj,jpkm1) / zwt(ji,jj,jpkm1)
545      END_2D
[13295]546      DO_3DS( 1, 1, 1, 1, jpk-2, 2, -1 )
[12377]547         pt_out(ji,jj,jk) = ( pt_out(ji,jj,jk) - zws(ji,jj,jk) * pt_out(ji,jj,jk+1) ) / zwt(ji,jj,jk)
548      END_3D
[14017]549      !
[7646]550   END SUBROUTINE interp_4th_cpt_org
551
[14017]552
[7646]553   SUBROUTINE interp_4th_cpt( pt_in, pt_out )
554      !!----------------------------------------------------------------------
555      !!                  ***  ROUTINE interp_4th_cpt  ***
[14017]556      !!
[7646]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    ! field at t-point
[14017]562      REAL(wp),DIMENSION(A2D(nn_hls)    ,jpk), INTENT(  out) ::   pt_out   ! field interpolated at w-point
[7646]563      !
564      INTEGER ::   ji, jj, jk   ! dummy loop integers
565      INTEGER ::   ikt, ikb     ! local integers
[14017]566      REAL(wp),DIMENSION(A2D(nn_hls),jpk) :: zwd, zwi, zws, zwrm, zwt
[7646]567      !!----------------------------------------------------------------------
568      !
569      !                      !==  build the three diagonal matrix & the RHS  ==!
570      !
[14017]571      DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 3, jpkm1 )    ! interior (from jk=3 to jpk-1)
[12377]572         zwd (ji,jj,jk) = 3._wp * wmask(ji,jj,jk) + 1._wp                 !       diagonal
573         zwi (ji,jj,jk) =         wmask(ji,jj,jk)                         ! lower diagonal
574         zws (ji,jj,jk) =         wmask(ji,jj,jk)                         ! upper diagonal
575         zwrm(ji,jj,jk) = 3._wp * wmask(ji,jj,jk)                     &   ! RHS
576            &           *       ( pt_in(ji,jj,jk) + pt_in(ji,jj,jk-1) )
577      END_3D
[7646]578      !
579!!gm
580!      SELECT CASE( kbc )               !* boundary condition
581!      CASE( np_NH   )   ! Neumann homogeneous at top & bottom
582!      CASE( np_CEN2 )   ! 2nd order centered  at top & bottom
583!      END SELECT
[14017]584!!gm
[7646]585      !
[9901]586      IF ( ln_isfcav ) THEN            ! set level two values which may not be set in ISF case
587         zwd(:,:,2) = 1._wp  ;  zwi(:,:,2) = 0._wp  ;  zws(:,:,2) = 0._wp  ;  zwrm(:,:,2) = 0._wp
588      END IF
589      !
[14017]590      DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )              ! 2nd order centered at top & bottom
[12377]591         ikt = mikt(ji,jj) + 1            ! w-point below the 1st  wet point
592         ikb = MAX(mbkt(ji,jj), 2)        !     -   above the last wet point
593         !
594         zwd (ji,jj,ikt) = 1._wp          ! top
595         zwi (ji,jj,ikt) = 0._wp
596         zws (ji,jj,ikt) = 0._wp
597         zwrm(ji,jj,ikt) = 0.5_wp * ( pt_in(ji,jj,ikt-1) + pt_in(ji,jj,ikt) )
598         !
599         zwd (ji,jj,ikb) = 1._wp          ! bottom
600         zwi (ji,jj,ikb) = 0._wp
601         zws (ji,jj,ikb) = 0._wp
[14017]602         zwrm(ji,jj,ikb) = 0.5_wp * ( pt_in(ji,jj,ikb-1) + pt_in(ji,jj,ikb) )
[12377]603      END_2D
[7646]604      !
605      !                       !==  tridiagonal solver  ==!
606      !
[14017]607      DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )           !* 1st recurrence:   Tk = Dk - Ik Sk-1 / Tk-1
[12377]608         zwt(ji,jj,2) = zwd(ji,jj,2)
609      END_2D
[14017]610      DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 3, jpkm1 )
[12377]611         zwt(ji,jj,jk) = zwd(ji,jj,jk) - zwi(ji,jj,jk) * zws(ji,jj,jk-1) /zwt(ji,jj,jk-1)
612      END_3D
[7646]613      !
[14017]614      DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )           !* 2nd recurrence:    Zk = Yk - Ik / Tk-1  Zk-1
[12377]615         pt_out(ji,jj,2) = zwrm(ji,jj,2)
616      END_2D
[14017]617      DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 3, jpkm1 )
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)
[12377]619      END_3D
[7646]620
[14017]621      DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )           !* 3d recurrence:    Xk = (Zk - Sk Xk+1 ) / Tk
[12377]622         pt_out(ji,jj,jpkm1) = pt_out(ji,jj,jpkm1) / zwt(ji,jj,jpkm1)
623      END_2D
[14017]624      DO_3DS( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, jpk-2, 2, -1 )
[12377]625         pt_out(ji,jj,jk) = ( pt_out(ji,jj,jk) - zws(ji,jj,jk) * pt_out(ji,jj,jk+1) ) / zwt(ji,jj,jk)
626      END_3D
[14017]627      !
[5770]628   END SUBROUTINE interp_4th_cpt
[7646]629
630
631   SUBROUTINE tridia_solver( pD, pU, pL, pRHS, pt_out , klev )
632      !!----------------------------------------------------------------------
633      !!                  ***  ROUTINE tridia_solver  ***
[14017]634      !!
[7646]635      !! **  Purpose :   solve a symmetric 3diagonal system
636      !!
637      !! **  Method  :   solve M.t_out = RHS(t)  where M is a tri diagonal matrix ( jpk*jpk )
[14017]638      !!
[7646]639      !!             ( D_1 U_1  0   0   0  )( t_1 )   ( RHS_1 )
640      !!             ( L_2 D_2 U_2  0   0  )( t_2 )   ( RHS_2 )
641      !!             (  0  L_3 D_3 U_3  0  )( t_3 ) = ( RHS_3 )
642      !!             (        ...          )( ... )   ( ...  )
643      !!             (  0   0   0  L_k D_k )( t_k )   ( RHS_k )
[14017]644      !!
[7646]645      !!        M is decomposed in the product of an upper and lower triangular matrix.
[14017]646      !!        The tri-diagonals matrix is given as input 3D arrays:   pD, pU, pL
[7646]647      !!        (i.e. the Diagonal, the Upper diagonal, and the Lower diagonal).
648      !!        The solution is pta.
649      !!        The 3d array zwt is used as a work space array.
650      !!----------------------------------------------------------------------
[14017]651      REAL(wp),DIMENSION(A2D(nn_hls),jpk), INTENT(in   ) ::   pD, pU, PL    ! 3-diagonal matrix
652      REAL(wp),DIMENSION(A2D(nn_hls),jpk), INTENT(in   ) ::   pRHS          ! Right-Hand-Side
653      REAL(wp),DIMENSION(A2D(nn_hls),jpk), INTENT(  out) ::   pt_out        !!gm field at level=F(klev)
654      INTEGER                    , INTENT(in   ) ::   klev          ! =1 pt_out at w-level
655      !                                                             ! =0 pt at t-level
[7646]656      INTEGER ::   ji, jj, jk   ! dummy loop integers
657      INTEGER ::   kstart       ! local indices
[14017]658      REAL(wp),DIMENSION(A2D(nn_hls),jpk) ::   zwt   ! 3D work array
[7646]659      !!----------------------------------------------------------------------
660      !
661      kstart =  1  + klev
662      !
[14017]663      DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )                         !* 1st recurrence:   Tk = Dk - Ik Sk-1 / Tk-1
[12377]664         zwt(ji,jj,kstart) = pD(ji,jj,kstart)
665      END_2D
[14017]666      DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, kstart+1, jpkm1 )
[12377]667         zwt(ji,jj,jk) = pD(ji,jj,jk) - pL(ji,jj,jk) * pU(ji,jj,jk-1) /zwt(ji,jj,jk-1)
668      END_3D
[7646]669      !
[14017]670      DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )                        !* 2nd recurrence:    Zk = Yk - Ik / Tk-1  Zk-1
[12377]671         pt_out(ji,jj,kstart) = pRHS(ji,jj,kstart)
672      END_2D
[14017]673      DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, kstart+1, jpkm1 )
674         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]675      END_3D
[7646]676
[14017]677      DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )                       !* 3d recurrence:    Xk = (Zk - Sk Xk+1 ) / Tk
[12377]678         pt_out(ji,jj,jpkm1) = pt_out(ji,jj,jpkm1) / zwt(ji,jj,jpkm1)
679      END_2D
[14017]680      DO_3DS( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, jpk-2, kstart, -1 )
[12377]681         pt_out(ji,jj,jk) = ( pt_out(ji,jj,jk) - pU(ji,jj,jk) * pt_out(ji,jj,jk+1) ) / zwt(ji,jj,jk)
682      END_3D
[7646]683      !
684   END SUBROUTINE tridia_solver
685
[3]686   !!======================================================================
[5770]687END MODULE traadv_fct
Note: See TracBrowser for help on using the repository browser.