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 branches/2017/dev_merge_2017/NEMOGCM/NEMO/OPA_SRC/TRA – NEMO

source: branches/2017/dev_merge_2017/NEMOGCM/NEMO/OPA_SRC/TRA/traadv_fct.F90 @ 9106

Last change on this file since 9106 was 9094, checked in by cetlod, 7 years ago

Use of lbclnk_multi in subdir LDF & TRA

  • 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) 
[3625]29   USE timing         ! Timing
[3]30
31   IMPLICIT NONE
32   PRIVATE
33
[9019]34   PUBLIC   tra_adv_fct        ! called by traadv.F90
35   PUBLIC   interp_4th_cpt     ! called by traadv_cen.F90
[3]36
[5770]37   LOGICAL  ::   l_trd   ! flag to compute trends
[7646]38   LOGICAL  ::   l_ptr   ! flag to compute poleward transport
39   LOGICAL  ::   l_hst   ! flag to compute heat/salt transport
[5770]40   REAL(wp) ::   r1_6 = 1._wp / 6._wp   ! =1/6
[2528]41
[7646]42   !                                        ! tridiag solver associated indices:
43   INTEGER, PARAMETER ::   np_NH   = 0   ! Neumann homogeneous boundary condition
44   INTEGER, PARAMETER ::   np_CEN2 = 1   ! 2nd order centered  boundary condition
45
[3]46   !! * Substitutions
47#  include "vectopt_loop_substitute.h90"
48   !!----------------------------------------------------------------------
[9019]49   !! NEMO/OPA 4.0 , NEMO Consortium (2017)
[1152]50   !! $Id$
[2528]51   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
[3]52   !!----------------------------------------------------------------------
53CONTAINS
54
[5770]55   SUBROUTINE tra_adv_fct( kt, kit000, cdtype, p2dt, pun, pvn, pwn,       &
56      &                                              ptb, ptn, pta, kjpt, kn_fct_h, kn_fct_v )
[3]57      !!----------------------------------------------------------------------
[5770]58      !!                  ***  ROUTINE tra_adv_fct  ***
[3]59      !!
[6140]60      !! **  Purpose :   Compute the now trend due to total advection of tracers
61      !!               and add it to the general trend of tracer equations
[3]62      !!
[5770]63      !! **  Method  : - 2nd or 4th FCT scheme on the horizontal direction
64      !!               (choice through the value of kn_fct)
[6140]65      !!               - on the vertical the 4th order is a compact scheme
[5770]66      !!               - corrected flux (monotonic correction)
[3]67      !!
[6140]68      !! ** Action : - update pta  with the now advective tracer trends
[9019]69      !!             - send trends to trdtra module for further diagnostics (l_trdtra=T)
[6140]70      !!             - htr_adv, str_adv : poleward advective heat and salt transport (ln_diaptr=T)
[503]71      !!----------------------------------------------------------------------
[2528]72      INTEGER                              , INTENT(in   ) ::   kt              ! ocean time-step index
[3294]73      INTEGER                              , INTENT(in   ) ::   kit000          ! first time step index
[2528]74      CHARACTER(len=3)                     , INTENT(in   ) ::   cdtype          ! =TRA or TRC (tracer indicator)
75      INTEGER                              , INTENT(in   ) ::   kjpt            ! number of tracers
[5770]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)
[6140]78      REAL(wp)                             , INTENT(in   ) ::   p2dt            ! tracer time-step
[2528]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
[2715]82      !
[5770]83      INTEGER  ::   ji, jj, jk, jn                           ! dummy loop indices 
[6140]84      REAL(wp) ::   ztra                                     ! local scalar
[5770]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   !   -      -
[9019]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
[3]89      !!----------------------------------------------------------------------
[3294]90      !
[9019]91      IF( ln_timing )   CALL timing_start('tra_adv_fct')
[3294]92      !
93      IF( kt == kit000 )  THEN
[2528]94         IF(lwp) WRITE(numout,*)
[5770]95         IF(lwp) WRITE(numout,*) 'tra_adv_fct : FCT advection scheme on ', cdtype
[2528]96         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
[3]97      ENDIF
[2528]98      !
[9019]99      l_trd = .FALSE.            ! set local switches
[7646]100      l_hst = .FALSE.
101      l_ptr = .FALSE.
[9019]102      IF( ( cdtype =='TRA' .AND. l_trdtra  ) .OR. ( cdtype =='TRC' .AND. l_trdtrc ) )       l_trd = .TRUE.
103      IF(   cdtype =='TRA' .AND. ln_diaptr )                                                l_ptr = .TRUE. 
104      IF(   cdtype =='TRA' .AND. ( iom_use("uadv_heattr") .OR. iom_use("vadv_heattr") .OR.  &
105         &                         iom_use("uadv_salttr") .OR. iom_use("vadv_salttr")  ) )  l_hst = .TRUE.
[5770]106      !
[7646]107      IF( l_trd .OR. l_hst )  THEN
[9019]108         ALLOCATE( ztrdx(jpi,jpj,jpk), ztrdy(jpi,jpj,jpk), ztrdz(jpi,jpj,jpk) )
[7753]109         ztrdx(:,:,:) = 0._wp   ;    ztrdy(:,:,:) = 0._wp   ;   ztrdz(:,:,:) = 0._wp
[3294]110      ENDIF
[2528]111      !
[7646]112      IF( l_ptr ) THEN 
[9019]113         ALLOCATE( zptry(jpi,jpj,jpk) )
[7753]114         zptry(:,:,:) = 0._wp
[7646]115      ENDIF
[6140]116      !                          ! surface & bottom value : flux set to zero one for all
[7753]117      zwz(:,:, 1 ) = 0._wp           
118      zwx(:,:,jpk) = 0._wp   ;   zwy(:,:,jpk) = 0._wp    ;    zwz(:,:,jpk) = 0._wp
[2528]119      !
[7753]120      zwi(:,:,:) = 0._wp       
121      !
[6140]122      DO jn = 1, kjpt            !==  loop over the tracers  ==!
[5770]123         !
124         !        !==  upstream advection with initial mass fluxes & intermediate update  ==!
125         !                    !* upstream tracer flux in the i and j direction
[2528]126         DO jk = 1, jpkm1
127            DO jj = 1, jpjm1
128               DO ji = 1, fs_jpim1   ! vector opt.
129                  ! upstream scheme
130                  zfp_ui = pun(ji,jj,jk) + ABS( pun(ji,jj,jk) )
131                  zfm_ui = pun(ji,jj,jk) - ABS( pun(ji,jj,jk) )
132                  zfp_vj = pvn(ji,jj,jk) + ABS( pvn(ji,jj,jk) )
133                  zfm_vj = pvn(ji,jj,jk) - ABS( pvn(ji,jj,jk) )
134                  zwx(ji,jj,jk) = 0.5 * ( zfp_ui * ptb(ji,jj,jk,jn) + zfm_ui * ptb(ji+1,jj  ,jk,jn) )
135                  zwy(ji,jj,jk) = 0.5 * ( zfp_vj * ptb(ji,jj,jk,jn) + zfm_vj * ptb(ji  ,jj+1,jk,jn) )
136               END DO
[3]137            END DO
138         END DO
[5770]139         !                    !* upstream tracer flux in the k direction *!
[6140]140         DO jk = 2, jpkm1        ! Interior value ( multiplied by wmask)
[4990]141            DO jj = 1, jpj
142               DO ji = 1, jpi
[2528]143                  zfp_wk = pwn(ji,jj,jk) + ABS( pwn(ji,jj,jk) )
144                  zfm_wk = pwn(ji,jj,jk) - ABS( pwn(ji,jj,jk) )
[5120]145                  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]146               END DO
[3]147            END DO
148         END DO
[6140]149         IF( ln_linssh ) THEN    ! top ocean value (only in linear free surface as zwz has been w-masked)
[5770]150            IF( ln_isfcav ) THEN             ! top of the ice-shelf cavities and at the ocean surface
[5120]151               DO jj = 1, jpj
152                  DO ji = 1, jpi
153                     zwz(ji,jj, mikt(ji,jj) ) = pwn(ji,jj,mikt(ji,jj)) * ptb(ji,jj,mikt(ji,jj),jn)   ! linear free surface
154                  END DO
155               END DO   
[5770]156            ELSE                             ! no cavities: only at the ocean surface
[7753]157               zwz(:,:,1) = pwn(:,:,1) * ptb(:,:,1,jn)
[5770]158            ENDIF
[5120]159         ENDIF
[5770]160         !               
161         DO jk = 1, jpkm1     !* trend and after field with monotonic scheme
[216]162            DO jj = 2, jpjm1
163               DO ji = fs_2, fs_jpim1   ! vector opt.
[6691]164                  !                             ! total intermediate advective trends
[5770]165                  ztra = - (  zwx(ji,jj,jk) - zwx(ji-1,jj  ,jk  )   &
166                     &      + zwy(ji,jj,jk) - zwy(ji  ,jj-1,jk  )   &
[6691]167                     &      + zwz(ji,jj,jk) - zwz(ji  ,jj  ,jk+1) ) * r1_e1e2t(ji,jj)
168                  !                             ! update and guess with monotonic sheme
169                  pta(ji,jj,jk,jn) =                     pta(ji,jj,jk,jn) +        ztra   / e3t_n(ji,jj,jk) * tmask(ji,jj,jk)
170                  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]171               END DO
172            END DO
173         END DO
[5770]174         CALL lbc_lnk( zwi, 'T', 1. )  ! Lateral boundary conditions on zwi  (unchanged sign)
175         !               
[7646]176         IF( l_trd .OR. l_hst )  THEN             ! trend diagnostics (contribution of upstream fluxes)
[7753]177            ztrdx(:,:,:) = zwx(:,:,:)   ;   ztrdy(:,:,:) = zwy(:,:,:)   ;   ztrdz(:,:,:) = zwz(:,:,:)
[2528]178         END IF
[5770]179         !                             ! "Poleward" heat and salt transports (contribution of upstream fluxes)
[9019]180         IF( l_ptr )   zptry(:,:,:) = zwy(:,:,:) 
[5770]181         !
182         !        !==  anti-diffusive flux : high order minus low order  ==!
183         !
[6140]184         SELECT CASE( kn_fct_h )    !* horizontal anti-diffusive fluxes
[5770]185         !
[6140]186         CASE(  2  )                   !- 2nd order centered
[5770]187            DO jk = 1, jpkm1
188               DO jj = 1, jpjm1
189                  DO ji = 1, fs_jpim1   ! vector opt.
190                     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)
191                     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)
192                  END DO
[503]193               END DO
194            END DO
[5770]195            !
[6140]196         CASE(  4  )                   !- 4th order centered
[7753]197            zltu(:,:,jpk) = 0._wp            ! Bottom value : flux set to zero
198            zltv(:,:,jpk) = 0._wp
[6140]199            DO jk = 1, jpkm1                 ! Laplacian
200               DO jj = 1, jpjm1                    ! 1st derivative (gradient)
[5770]201                  DO ji = 1, fs_jpim1   ! vector opt.
202                     ztu(ji,jj,jk) = ( ptn(ji+1,jj  ,jk,jn) - ptn(ji,jj,jk,jn) ) * umask(ji,jj,jk)
203                     ztv(ji,jj,jk) = ( ptn(ji  ,jj+1,jk,jn) - ptn(ji,jj,jk,jn) ) * vmask(ji,jj,jk)
204                  END DO
[503]205               END DO
[6140]206               DO jj = 2, jpjm1                    ! 2nd derivative * 1/ 6
[5770]207                  DO ji = fs_2, fs_jpim1   ! vector opt.
208                     zltu(ji,jj,jk) = (  ztu(ji,jj,jk) + ztu(ji-1,jj,jk)  ) * r1_6
209                     zltv(ji,jj,jk) = (  ztv(ji,jj,jk) + ztv(ji,jj-1,jk)  ) * r1_6
210                  END DO
211               END DO
[503]212            END DO
[9094]213            CALL lbc_lnk_multi( zltu, 'T', 1. , zltv, 'T', 1. )   ! Lateral boundary cond. (unchanged sgn)
[5770]214            !
[6140]215            DO jk = 1, jpkm1                 ! Horizontal advective fluxes
[5770]216               DO jj = 1, jpjm1
217                  DO ji = 1, fs_jpim1   ! vector opt.
218                     zC2t_u = ptn(ji,jj,jk,jn) + ptn(ji+1,jj  ,jk,jn)   ! 2 x C2 interpolation of T at u- & v-points
219                     zC2t_v = ptn(ji,jj,jk,jn) + ptn(ji  ,jj+1,jk,jn)
220                     !                                                  ! C4 minus upstream advective fluxes
221                     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)
222                     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)
223                  END DO
[5120]224               END DO
[5770]225            END DO         
226            !
[6140]227         CASE(  41 )                   !- 4th order centered       ==>>   !!gm coding attempt   need to be tested
[7753]228            ztu(:,:,jpk) = 0._wp             ! Bottom value : flux set to zero
229            ztv(:,:,jpk) = 0._wp
[6140]230            DO jk = 1, jpkm1                 ! 1st derivative (gradient)
231               DO jj = 1, jpjm1
[5770]232                  DO ji = 1, fs_jpim1   ! vector opt.
233                     ztu(ji,jj,jk) = ( ptn(ji+1,jj  ,jk,jn) - ptn(ji,jj,jk,jn) ) * umask(ji,jj,jk)
234                     ztv(ji,jj,jk) = ( ptn(ji  ,jj+1,jk,jn) - ptn(ji,jj,jk,jn) ) * vmask(ji,jj,jk)
235                  END DO
236               END DO
[5120]237            END DO
[9094]238            CALL lbc_lnk_multi( ztu, 'U', -1. , ztv, 'V', -1. )   ! Lateral boundary cond. (unchanged sgn)
[5770]239            !
[6140]240            DO jk = 1, jpkm1                 ! Horizontal advective fluxes
[5770]241               DO jj = 2, jpjm1
242                  DO ji = 2, fs_jpim1   ! vector opt.
243                     zC2t_u = ptn(ji,jj,jk,jn) + ptn(ji+1,jj  ,jk,jn)   ! 2 x C2 interpolation of T at u- & v-points (x2)
244                     zC2t_v = ptn(ji,jj,jk,jn) + ptn(ji  ,jj+1,jk,jn)
245                     !                                                  ! C4 interpolation of T at u- & v-points (x2)
246                     zC4t_u =  zC2t_u + r1_6 * ( ztu(ji-1,jj  ,jk) - ztu(ji+1,jj  ,jk) )
247                     zC4t_v =  zC2t_v + r1_6 * ( ztv(ji  ,jj-1,jk) - ztv(ji  ,jj+1,jk) )
248                     !                                                  ! C4 minus upstream advective fluxes
249                     zwx(ji,jj,jk) =  0.5_wp * pun(ji,jj,jk) * zC4t_u - zwx(ji,jj,jk)
250                     zwy(ji,jj,jk) =  0.5_wp * pvn(ji,jj,jk) * zC4t_v - zwy(ji,jj,jk)
251                  END DO
252               END DO
253            END DO
254            !
255         END SELECT
[6140]256         !                     
257         SELECT CASE( kn_fct_v )    !* vertical anti-diffusive fluxes (w-masked interior values)
[5770]258         !
[6140]259         CASE(  2  )                   !- 2nd order centered
[5770]260            DO jk = 2, jpkm1   
261               DO jj = 2, jpjm1
262                  DO ji = fs_2, fs_jpim1
[6140]263                     zwz(ji,jj,jk) =  (  pwn(ji,jj,jk) * 0.5_wp * ( ptn(ji,jj,jk,jn) + ptn(ji,jj,jk-1,jn) )   &
264                        &              - zwz(ji,jj,jk)  ) * wmask(ji,jj,jk)
[5770]265                  END DO
266               END DO
267            END DO
268            !
[6140]269         CASE(  4  )                   !- 4th order COMPACT
270            CALL interp_4th_cpt( ptn(:,:,:,jn) , ztw )   ! zwt = COMPACT interpolation of T at w-point
[5770]271            DO jk = 2, jpkm1
272               DO jj = 2, jpjm1
273                  DO ji = fs_2, fs_jpim1
274                     zwz(ji,jj,jk) = ( pwn(ji,jj,jk) * ztw(ji,jj,jk) - zwz(ji,jj,jk) ) * wmask(ji,jj,jk)
275                  END DO
276               END DO
277            END DO
278            !
279         END SELECT
[6140]280         IF( ln_linssh ) THEN    ! top ocean value: high order = upstream  ==>>  zwz=0
[7753]281            zwz(:,:,1) = 0._wp   ! only ocean surface as interior zwz values have been w-masked
[6140]282         ENDIF
[5770]283         !
[9094]284         CALL lbc_lnk_multi( zwx, 'U', -1. , zwy, 'V', -1.,  zwz, 'W',  1. )
[6140]285         !
[5770]286         !        !==  monotonicity algorithm  ==!
287         !
[2528]288         CALL nonosc( ptb(:,:,:,jn), zwx, zwy, zwz, zwi, p2dt )
[6140]289         !
[5770]290         !        !==  final trend with corrected fluxes  ==!
291         !
[216]292         DO jk = 1, jpkm1
293            DO jj = 2, jpjm1
[2528]294               DO ji = fs_2, fs_jpim1   ! vector opt. 
[5770]295                  pta(ji,jj,jk,jn) = pta(ji,jj,jk,jn) - (  zwx(ji,jj,jk) - zwx(ji-1,jj  ,jk  )   &
296                     &                                   + zwy(ji,jj,jk) - zwy(ji  ,jj-1,jk  )   &
297                     &                                   + zwz(ji,jj,jk) - zwz(ji  ,jj  ,jk+1) ) &
[6140]298                     &                                * r1_e1e2t(ji,jj) / e3t_n(ji,jj,jk)
[216]299               END DO
300            END DO
301         END DO
[5770]302         !
[9019]303         IF( l_trd .OR. l_hst ) THEN   ! trend diagnostics // heat/salt transport
304            ztrdx(:,:,:) = ztrdx(:,:,:) + zwx(:,:,:)  ! <<< add anti-diffusive fluxes
305            ztrdy(:,:,:) = ztrdy(:,:,:) + zwy(:,:,:)  !     to upstream fluxes
306            ztrdz(:,:,:) = ztrdz(:,:,:) + zwz(:,:,:)  !
[5770]307            !
[9019]308            IF( l_trd ) THEN              ! trend diagnostics
309               CALL trd_tra( kt, cdtype, jn, jptra_xad, ztrdx, pun, ptn(:,:,:,jn) )
310               CALL trd_tra( kt, cdtype, jn, jptra_yad, ztrdy, pvn, ptn(:,:,:,jn) )
311               CALL trd_tra( kt, cdtype, jn, jptra_zad, ztrdz, pwn, ptn(:,:,:,jn) )
312            ENDIF
313            !                             ! heat/salt transport
314            IF( l_hst )   CALL dia_ar5_hst( jn, 'adv', ztrdx(:,:,:), ztrdy(:,:,:) )
[5770]315            !
[9019]316            DEALLOCATE( ztrdx, ztrdy, ztrdz )
317         ENDIF
318         IF( l_ptr ) THEN              ! "Poleward" transports
319            zptry(:,:,:) = zptry(:,:,:) + zwy(:,:,:)  ! <<< add anti-diffusive fluxes
[7646]320            CALL dia_ptr_hst( jn, 'adv', zptry(:,:,:) )
[9019]321            DEALLOCATE( zptry )
[2528]322         ENDIF
[503]323         !
[6140]324      END DO                     ! end of tracer loop
[503]325      !
[9019]326      IF( ln_timing )   CALL timing_stop('tra_adv_fct')
[2528]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      !
[9019]355      IF( ln_timing )   CALL timing_start('nonosc')
[3294]356      !
[2715]357      zbig  = 1.e+40_wp
358      zrtrn = 1.e-15_wp
[7753]359      zbetup(:,:,:) = 0._wp   ;   zbetdo(:,:,:) = 0._wp
[785]360
[3]361      ! Search local extrema
362      ! --------------------
[785]363      ! max/min of pbef & paft with large negative/positive value (-/+zbig) inside land
[4990]364      zbup = MAX( pbef * tmask - zbig * ( 1._wp - tmask ),   &
365         &        paft * tmask - zbig * ( 1._wp - tmask )  )
366      zbdo = MIN( pbef * tmask + zbig * ( 1._wp - tmask ),   &
367         &        paft * tmask + zbig * ( 1._wp - tmask )  )
[785]368
[5120]369      DO jk = 1, jpkm1
370         ikm1 = MAX(jk-1,1)
371         DO jj = 2, jpjm1
372            DO ji = fs_2, fs_jpim1   ! vector opt.
373
[785]374               ! search maximum in neighbourhood
375               zup = MAX(  zbup(ji  ,jj  ,jk  ),   &
376                  &        zbup(ji-1,jj  ,jk  ), zbup(ji+1,jj  ,jk  ),   &
377                  &        zbup(ji  ,jj-1,jk  ), zbup(ji  ,jj+1,jk  ),   &
378                  &        zbup(ji  ,jj  ,ikm1), zbup(ji  ,jj  ,jk+1)  )
[3]379
[785]380               ! search minimum in neighbourhood
381               zdo = MIN(  zbdo(ji  ,jj  ,jk  ),   &
382                  &        zbdo(ji-1,jj  ,jk  ), zbdo(ji+1,jj  ,jk  ),   &
383                  &        zbdo(ji  ,jj-1,jk  ), zbdo(ji  ,jj+1,jk  ),   &
384                  &        zbdo(ji  ,jj  ,ikm1), zbdo(ji  ,jj  ,jk+1)  )
[3]385
[785]386               ! positive part of the flux
[3]387               zpos = MAX( 0., paa(ji-1,jj  ,jk  ) ) - MIN( 0., paa(ji  ,jj  ,jk  ) )   &
388                  & + MAX( 0., pbb(ji  ,jj-1,jk  ) ) - MIN( 0., pbb(ji  ,jj  ,jk  ) )   &
389                  & + MAX( 0., pcc(ji  ,jj  ,jk+1) ) - MIN( 0., pcc(ji  ,jj  ,jk  ) )
[785]390
391               ! negative part of the flux
[3]392               zneg = MAX( 0., paa(ji  ,jj  ,jk  ) ) - MIN( 0., paa(ji-1,jj  ,jk  ) )   &
393                  & + MAX( 0., pbb(ji  ,jj  ,jk  ) ) - MIN( 0., pbb(ji  ,jj-1,jk  ) )   &
394                  & + MAX( 0., pcc(ji  ,jj  ,jk  ) ) - MIN( 0., pcc(ji  ,jj  ,jk+1) )
[785]395
[3]396               ! up & down beta terms
[6140]397               zbt = e1e2t(ji,jj) * e3t_n(ji,jj,jk) / p2dt
[785]398               zbetup(ji,jj,jk) = ( zup            - paft(ji,jj,jk) ) / ( zpos + zrtrn ) * zbt
399               zbetdo(ji,jj,jk) = ( paft(ji,jj,jk) - zdo            ) / ( zneg + zrtrn ) * zbt
[3]400            END DO
401         END DO
402      END DO
[9094]403      CALL lbc_lnk_multi( zbetup, 'T', 1. , zbetdo, 'T', 1. )   ! lateral boundary cond. (unchanged sign)
[3]404
[237]405      ! 3. monotonic flux in the i & j direction (paa & pbb)
406      ! ----------------------------------------
[3]407      DO jk = 1, jpkm1
408         DO jj = 2, jpjm1
409            DO ji = fs_2, fs_jpim1   ! vector opt.
[4990]410               zau = MIN( 1._wp, zbetdo(ji,jj,jk), zbetup(ji+1,jj,jk) )
411               zbu = MIN( 1._wp, zbetup(ji,jj,jk), zbetdo(ji+1,jj,jk) )
[785]412               zcu =       ( 0.5  + SIGN( 0.5 , paa(ji,jj,jk) ) )
[4990]413               paa(ji,jj,jk) = paa(ji,jj,jk) * ( zcu * zau + ( 1._wp - zcu) * zbu )
[3]414
[4990]415               zav = MIN( 1._wp, zbetdo(ji,jj,jk), zbetup(ji,jj+1,jk) )
416               zbv = MIN( 1._wp, zbetup(ji,jj,jk), zbetdo(ji,jj+1,jk) )
[785]417               zcv =       ( 0.5  + SIGN( 0.5 , pbb(ji,jj,jk) ) )
[4990]418               pbb(ji,jj,jk) = pbb(ji,jj,jk) * ( zcv * zav + ( 1._wp - zcv) * zbv )
[3]419
420      ! monotonic flux in the k direction, i.e. pcc
421      ! -------------------------------------------
[785]422               za = MIN( 1., zbetdo(ji,jj,jk+1), zbetup(ji,jj,jk) )
423               zb = MIN( 1., zbetup(ji,jj,jk+1), zbetdo(ji,jj,jk) )
424               zc =       ( 0.5  + SIGN( 0.5 , pcc(ji,jj,jk+1) ) )
[4990]425               pcc(ji,jj,jk+1) = pcc(ji,jj,jk+1) * ( zc * za + ( 1._wp - zc) * zb )
[3]426            END DO
427         END DO
428      END DO
[9094]429      CALL lbc_lnk_multi( paa, 'U', -1. , pbb, 'V', -1. )   ! lateral boundary condition (changed sign)
[503]430      !
[9019]431      IF( ln_timing )   CALL timing_stop('nonosc')
[2715]432      !
[3]433   END SUBROUTINE nonosc
434
[5770]435
[7646]436   SUBROUTINE interp_4th_cpt_org( pt_in, pt_out )
[5770]437      !!----------------------------------------------------------------------
[7646]438      !!                  ***  ROUTINE interp_4th_cpt_org  ***
[5770]439      !!
440      !! **  Purpose :   Compute the interpolation of tracer at w-point
441      !!
442      !! **  Method  :   4th order compact interpolation
443      !!----------------------------------------------------------------------
444      REAL(wp),DIMENSION(jpi,jpj,jpk), INTENT(in   ) ::   pt_in    ! now tracer fields
445      REAL(wp),DIMENSION(jpi,jpj,jpk), INTENT(  out) ::   pt_out   ! now tracer field interpolated at w-pts
446      !
447      INTEGER :: ji, jj, jk   ! dummy loop integers
448      REAL(wp),DIMENSION(jpi,jpj,jpk) :: zwd, zwi, zws, zwrm, zwt
449      !!----------------------------------------------------------------------
450     
451      DO jk = 3, jpkm1        !==  build the three diagonal matrix  ==!
452         DO jj = 1, jpj
453            DO ji = 1, jpi
454               zwd (ji,jj,jk) = 4._wp
455               zwi (ji,jj,jk) = 1._wp
456               zws (ji,jj,jk) = 1._wp
457               zwrm(ji,jj,jk) = 3._wp * ( pt_in(ji,jj,jk-1) + pt_in(ji,jj,jk) )
458               !
459               IF( tmask(ji,jj,jk+1) == 0._wp) THEN   ! Switch to second order centered at bottom
460                  zwd (ji,jj,jk) = 1._wp
461                  zwi (ji,jj,jk) = 0._wp
462                  zws (ji,jj,jk) = 0._wp
463                  zwrm(ji,jj,jk) = 0.5 * ( pt_in(ji,jj,jk-1) + pt_in(ji,jj,jk) )   
464               ENDIF
465            END DO
466         END DO
467      END DO
468      !
[7646]469      jk = 2                                          ! Switch to second order centered at top
470      DO jj = 1, jpj
471         DO ji = 1, jpi
[5770]472            zwd (ji,jj,jk) = 1._wp
473            zwi (ji,jj,jk) = 0._wp
474            zws (ji,jj,jk) = 0._wp
475            zwrm(ji,jj,jk) = 0.5 * ( pt_in(ji,jj,jk-1) + pt_in(ji,jj,jk) )
476         END DO
477      END DO   
478      !
479      !                       !==  tridiagonal solve  ==!
480      DO jj = 1, jpj                ! first recurrence
481         DO ji = 1, jpi
482            zwt(ji,jj,2) = zwd(ji,jj,2)
483         END DO
484      END DO
485      DO jk = 3, jpkm1
486         DO jj = 1, jpj
487            DO ji = 1, jpi
488               zwt(ji,jj,jk) = zwd(ji,jj,jk) - zwi(ji,jj,jk) * zws(ji,jj,jk-1) /zwt(ji,jj,jk-1)
489            END DO
490         END DO
491      END DO
492      !
493      DO jj = 1, jpj                ! second recurrence:    Zk = Yk - Ik / Tk-1  Zk-1
494         DO ji = 1, jpi
495            pt_out(ji,jj,2) = zwrm(ji,jj,2)
496         END DO
497      END DO
498      DO jk = 3, jpkm1
499         DO jj = 1, jpj
500            DO ji = 1, jpi
501               pt_out(ji,jj,jk) = zwrm(ji,jj,jk) - zwi(ji,jj,jk) / zwt(ji,jj,jk-1) *pt_out(ji,jj,jk-1)             
502            END DO
503         END DO
504      END DO
505
506      DO jj = 1, jpj                ! third recurrence: Xk = (Zk - Sk Xk+1 ) / Tk
507         DO ji = 1, jpi
508            pt_out(ji,jj,jpkm1) = pt_out(ji,jj,jpkm1) / zwt(ji,jj,jpkm1)
509         END DO
510      END DO
511      DO jk = jpk-2, 2, -1
512         DO jj = 1, jpj
513            DO ji = 1, jpi
514               pt_out(ji,jj,jk) = ( pt_out(ji,jj,jk) - zws(ji,jj,jk) * pt_out(ji,jj,jk+1) ) / zwt(ji,jj,jk)
515            END DO
516         END DO
517      END DO
518      !   
[7646]519   END SUBROUTINE interp_4th_cpt_org
520   
521
522   SUBROUTINE interp_4th_cpt( pt_in, pt_out )
523      !!----------------------------------------------------------------------
524      !!                  ***  ROUTINE interp_4th_cpt  ***
525      !!
526      !! **  Purpose :   Compute the interpolation of tracer at w-point
527      !!
528      !! **  Method  :   4th order compact interpolation
529      !!----------------------------------------------------------------------
530      REAL(wp),DIMENSION(jpi,jpj,jpk), INTENT(in   ) ::   pt_in    ! field at t-point
531      REAL(wp),DIMENSION(jpi,jpj,jpk), INTENT(  out) ::   pt_out   ! field interpolated at w-point
532      !
533      INTEGER ::   ji, jj, jk   ! dummy loop integers
534      INTEGER ::   ikt, ikb     ! local integers
535      REAL(wp),DIMENSION(jpi,jpj,jpk) :: zwd, zwi, zws, zwrm, zwt
536      !!----------------------------------------------------------------------
537      !
538      !                      !==  build the three diagonal matrix & the RHS  ==!
539      !
540      DO jk = 3, jpkm1                 ! interior (from jk=3 to jpk-1)
541         DO jj = 2, jpjm1
542            DO ji = fs_2, fs_jpim1
543               zwd (ji,jj,jk) = 3._wp * wmask(ji,jj,jk) + 1._wp                 !       diagonal
544               zwi (ji,jj,jk) =         wmask(ji,jj,jk)                         ! lower diagonal
545               zws (ji,jj,jk) =         wmask(ji,jj,jk)                         ! upper diagonal
546               zwrm(ji,jj,jk) = 3._wp * wmask(ji,jj,jk)                     &   ! RHS
547                  &           *       ( pt_in(ji,jj,jk) + pt_in(ji,jj,jk-1) )
548            END DO
549         END DO
550      END DO
551      !
552!!gm
553!      SELECT CASE( kbc )               !* boundary condition
554!      CASE( np_NH   )   ! Neumann homogeneous at top & bottom
555!      CASE( np_CEN2 )   ! 2nd order centered  at top & bottom
556!      END SELECT
557!!gm 
558      !
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
567            zwrm(ji,jj,ikt) = 0.5_wp * ( pt_in(ji,jj,jk-1) + pt_in(ji,jj,jk) )
568            !
569            zwd (ji,jj,ikb) = 1._wp          ! bottom
570            zwi (ji,jj,ikb) = 0._wp
571            zws (ji,jj,ikb) = 0._wp
572            zwrm(ji,jj,ikb) = 0.5_wp * ( pt_in(ji,jj,jk-1) + pt_in(ji,jj,jk) )           
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.