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

source: NEMO/trunk/src/OCE/TRA/traadv_fct.F90 @ 10068

Last change on this file since 10068 was 10068, checked in by nicolasmartin, 6 years ago

First part of modifications to have a common default header : fix typos and SVN keywords properties

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