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.
trazdf_imp.F90 in branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/TRA – NEMO

source: branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/TRA/trazdf_imp.F90 @ 4429

Last change on this file since 4429 was 4429, checked in by trackstand2, 10 years ago

Fixed missing -1 on a z-loop limit

  • Property svn:keywords set to Id
File size: 15.9 KB
RevLine 
[3]1MODULE trazdf_imp
[1438]2   !!======================================================================
[457]3   !!                 ***  MODULE  trazdf_imp  ***
[2528]4   !! Ocean  tracers:  vertical component of the tracer mixing trend
[1438]5   !!======================================================================
6   !! History :  OPA  !  1990-10  (B. Blanke)  Original code
7   !!            7.0  !  1991-11  (G. Madec)
8   !!                 !  1992-06  (M. Imbard) correction on tracer trend loops
9   !!                 !  1996-01  (G. Madec) statement function for e3
10   !!                 !  1997-05  (G. Madec) vertical component of isopycnal
11   !!                 !  1997-07  (G. Madec) geopotential diffusion in s-coord
12   !!                 !  2000-08  (G. Madec) double diffusive mixing
13   !!   NEMO     1.0  !  2002-08  (G. Madec) F90: Free form and module
14   !!            2.0  !  2006-11  (G. Madec) New step reorganisation
15   !!            3.2  !  2009-03  (G. Madec)  heat and salt content trends
[2528]16   !!            3.3  !  2010-06  (C. Ethe, G. Madec) Merge TRA-TRC
[2602]17   !!             -   !  2011-02  (A. Coward, C. Ethe, G. Madec) improvment of surface boundary condition
[3]18   !!----------------------------------------------------------------------
[1438]19 
20   !!----------------------------------------------------------------------
[457]21   !!   tra_zdf_imp : Update the tracer trend with the diagonal vertical 
22   !!                 part of the mixing tensor.
[3]23   !!----------------------------------------------------------------------
[457]24   USE oce             ! ocean dynamics and tracers variables
25   USE dom_oce         ! ocean space and time domain variables
26   USE zdf_oce         ! ocean vertical physics variables
[2528]27   USE trc_oce         ! share passive tracers/ocean variables
28   USE domvvl          ! variable volume
[216]29   USE ldftra_oce      ! ocean active tracers: lateral physics
[2528]30   USE ldftra          ! lateral mixing type
[457]31   USE ldfslp          ! lateral physics: slope of diffusion
32   USE zdfddm          ! ocean vertical physics: double diffusion
[2528]33   USE traldf_iso_grif ! active tracers: Griffies operator
[3]34   USE in_out_manager  ! I/O manager
[457]35   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
[2715]36   USE lib_mpp         ! MPP library
[3]37
38   IMPLICIT NONE
39   PRIVATE
40
[1438]41   PUBLIC   tra_zdf_imp   !  routine called by step.F90
[3]42
[2602]43   REAL(wp) ::  r_vvl     ! variable volume indicator, =1 if lk_vvl=T, =0 otherwise
44
[3211]45   !! * Control permutation of array indices
46#  include "oce_ftrans.h90"
47#  include "dom_oce_ftrans.h90"
48#  include "zdf_oce_ftrans.h90"
49#  include "trc_oce_ftrans.h90"
50#  include "domvvl_ftrans.h90"
51#  include "ldftra_oce_ftrans.h90"
52#  include "ldfslp_ftrans.h90"
53#  include "zdfddm_ftrans.h90"
54#  include "traldf_iso_grif_ftrans.h90"
55
[3]56   !! * Substitutions
57#  include "domzgr_substitute.h90"
[457]58#  include "ldftra_substitute.h90"
[3]59#  include "zdfddm_substitute.h90"
[457]60#  include "vectopt_loop_substitute.h90"
[3]61   !!----------------------------------------------------------------------
[2528]62   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
[1156]63   !! $Id$
[2528]64   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
[457]65   !!----------------------------------------------------------------------
[3]66CONTAINS
[2528]67 
68   SUBROUTINE tra_zdf_imp( kt, cdtype, p2dt, ptb, pta, kjpt ) 
[3]69      !!----------------------------------------------------------------------
70      !!                  ***  ROUTINE tra_zdf_imp  ***
71      !!
[2602]72      !! ** Purpose :   Compute the after tracer through a implicit computation
73      !!     of the vertical tracer diffusion (including the vertical component
74      !!     of lateral mixing (only for 2nd order operator, for fourth order
75      !!     it is already computed and add to the general trend in traldf)
[3]76      !!
[2602]77      !! ** Method  :  The vertical diffusion of the tracer t  is given by:
78      !!                  difft = dz( avt dz(t) ) = 1/e3t dk+1( avt/e3w dk(t) )
[457]79      !!      It is computed using a backward time scheme (t=ta).
[2602]80      !!      If lk_zdfddm=T, use avs for salinity or for passive tracers
[3]81      !!      Surface and bottom boundary conditions: no diffusive flux on
82      !!      both tracers (bottom, applied through the masked field avt).
[2602]83      !!      If iso-neutral mixing, add to avt the contribution due to lateral mixing.
[3]84      !!
[2602]85      !! ** Action  : - pta  becomes the after tracer
[3]86      !!---------------------------------------------------------------------
[3432]87      USE timing,   ONLY:   timing_start, timing_stop
[2715]88      USE wrk_nemo, ONLY:   wrk_in_use, wrk_not_released
89      USE oce     , ONLY:   zwd => ua       , zws => va         ! (ua,va) used as 3D workspace
90      USE wrk_nemo, ONLY:   zwi => wrk_3d_6 , zwt => wrk_3d_7   ! 3D workspace
[3211]91
92      !! DCSE_NEMO: Need additional directives for renamed module variables
93!FTRANS zwd zws :I :I :z
94!FTRANS zwi zwt :I :I :z
[2715]95      !
[2528]96      INTEGER                              , INTENT(in   ) ::   kt       ! ocean time-step index
97      CHARACTER(len=3)                     , INTENT(in   ) ::   cdtype   ! =TRA or TRC (tracer indicator)
98      INTEGER                              , INTENT(in   ) ::   kjpt     ! number of tracers
[4409]99      REAL(wp), DIMENSION(       jpkorig  ), INTENT(in   ) ::   p2dt     ! vertical profile of tracer time-step
[3211]100
101      !! DCSE_NEMO: This style defeats ftrans
102!     REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(in   ) ::   ptb      ! before and now tracer fields
103!     REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::   pta      ! tracer trend
104
105!FTRANS ptb pta :I :I :z :
[4409]106      REAL(wp), INTENT(in   ) ::   ptb(jpi,jpj,jpkorig,kjpt)      ! before and now tracer fields
107      REAL(wp), INTENT(inout) ::   pta(jpi,jpj,jpkorig,kjpt)      ! tracer trend
[2715]108      !
109      INTEGER  ::  ji, jj, jk, jn   ! dummy loop indices
110      REAL(wp) ::  zrhs, ze3tb, ze3tn, ze3ta   ! local scalars
[3]111      !!---------------------------------------------------------------------
112
[3432]113      CALL timing_start('tra_zdf_imp')
114
[2715]115      IF( wrk_in_use(3, 6,7) ) THEN
116         CALL ctl_stop('tra_zdf_imp : requested workspace arrays unavailable.')   ;   RETURN
117      ENDIF
118
[2528]119      IF( kt == nit000 )  THEN
[457]120         IF(lwp)WRITE(numout,*)
[2528]121         IF(lwp)WRITE(numout,*) 'tra_zdf_imp : implicit vertical mixing on ', cdtype
[457]122         IF(lwp)WRITE(numout,*) '~~~~~~~~~~~ '
[2602]123         !
124         IF( lk_vvl ) THEN   ;    r_vvl = 1._wp       ! Variable volume indicator
125         ELSE                ;    r_vvl = 0._wp       
126         ENDIF
[457]127      ENDIF
[2528]128      !
[2602]129      !                                               ! ============= !
130      DO jn = 1, kjpt                                 !  tracer loop  !
131         !                                            ! ============= !
[2528]132         !
133         !  Matrix construction
[2602]134         ! --------------------
135         ! Build matrix if temperature or salinity (only in double diffusion case) or first passive tracer
136         !
[2715]137         IF(  ( cdtype == 'TRA' .AND. ( jn == jp_tem .OR. ( jn == jp_sal .AND. lk_zdfddm ) ) ) .OR.   &
[2602]138            & ( cdtype == 'TRC' .AND. jn == 1 )  )  THEN
139            !
140            ! vertical mixing coef.: avt for temperature, avs for salinity and passive tracers
[3211]141#if defined key_z_first
142            IF( cdtype == 'TRA' .AND. jn == jp_tem ) THEN
143               DO jj = 1, jpj
144                  DO ji = 1, jpi
145                     zwt(ji,jj,1) = 0._wp
[4428]146                     DO jk = 2, mbkmax(ji,jj) ! jpk
[3211]147                        zwt(ji,jj,jk) = avt  (ji,jj,jk)
148                     END DO
149                  END DO
150               END DO
151            ELSE                               
152               DO jj = 1, jpj
153                  DO ji = 1, jpi
154                     zwt(ji,jj,1) = 0._wp
[4428]155                     DO jk = 2, mbkmax(ji,jj) ! jpk
[3211]156                        zwt(ji,jj,jk) = fsavs(ji,jj,jk)
157                     END DO
158                  END DO
159               END DO
160            ENDIF
161#else
[2602]162            IF( cdtype == 'TRA' .AND. jn == jp_tem ) THEN   ;   zwt(:,:,2:jpk) = avt  (:,:,2:jpk)
163            ELSE                                            ;   zwt(:,:,2:jpk) = fsavs(:,:,2:jpk)
[2528]164            ENDIF
[2602]165            zwt(:,:,1) = 0._wp
[3211]166#endif
[2528]167            !
168#if defined key_ldfslp
[2602]169            ! isoneutral diffusion: add the contribution
170            IF( ln_traldf_grif    ) THEN     ! Griffies isoneutral diff
[3211]171#if defined key_z_first
172               DO jj = 2, jpjm1
173                  DO ji = 2, jpim1
[4428]174                     DO jk = 2, mbkmax(ji,jj)-1 ! jpkm1
[3211]175#else
[2528]176               DO jk = 2, jpkm1
177                  DO jj = 2, jpjm1
178                     DO ji = fs_2, fs_jpim1   ! vector opt.
[3211]179#endif
[2602]180                        zwt(ji,jj,jk) = zwt(ji,jj,jk) + ah_wslp2(ji,jj,jk)       
[2528]181                     END DO
182                  END DO
183               END DO
[2602]184            ELSE IF( l_traldf_rot ) THEN     ! standard isoneutral diff
[3211]185#if defined key_z_first
186               DO jj = 2, jpjm1
187                  DO ji = 2, jpim1
[4428]188                     DO jk = 2, mbkmax(ji,jj)-1 ! jpkm1
[3211]189#else
[2528]190               DO jk = 2, jpkm1
191                  DO jj = 2, jpjm1
192                     DO ji = fs_2, fs_jpim1   ! vector opt.
[3211]193#endif
[2602]194                        zwt(ji,jj,jk) = zwt(ji,jj,jk) + fsahtw(ji,jj,jk)                       &
195                           &                          * (  wslpi(ji,jj,jk) * wslpi(ji,jj,jk)   &
196                           &                             + wslpj(ji,jj,jk) * wslpj(ji,jj,jk)  )
[2528]197                     END DO
198                  END DO
199               END DO
200            ENDIF
[592]201#endif
[2602]202            ! Diagonal, lower (i), upper (s)  (including the bottom boundary condition since avt is masked)
[3211]203#if defined key_z_first
204            DO jj = 2, jpjm1
205               DO ji = 2, jpim1
[4429]206                  DO jk = 1, mbkmax(ji,jj)-1 ! jpkm1
[3432]207                     ! after scale factor at T-point
208                     ze3ta =  ( 1. - r_vvl ) +        r_vvl   * fse3t_a(ji,jj,jk) 
209                     ! now   scale factor at T-point
210                     ze3tn =         r_vvl   + ( 1. - r_vvl ) * fse3t_n(ji,jj,jk) 
[3211]211                     zwi(ji,jj,jk) = - p2dt(jk) * zwt(ji,jj,jk  ) / ( ze3tn * fse3w(ji,jj,jk  ) )
212                     zws(ji,jj,jk) = - p2dt(jk) * zwt(ji,jj,jk+1) / ( ze3tn * fse3w(ji,jj,jk+1) )
213                     zwd(ji,jj,jk) = ze3ta - zwi(ji,jj,jk) - zws(ji,jj,jk)
214                 END DO
215#else
[2528]216            DO jk = 1, jpkm1
217               DO jj = 2, jpjm1
218                  DO ji = fs_2, fs_jpim1   ! vector opt.
[2602]219                     ze3ta =  ( 1. - r_vvl ) +        r_vvl   * fse3t_a(ji,jj,jk)   ! after scale factor at T-point
220                     ze3tn =         r_vvl   + ( 1. - r_vvl ) * fse3t_n(ji,jj,jk)   ! now   scale factor at T-point
[2528]221                     zwi(ji,jj,jk) = - p2dt(jk) * zwt(ji,jj,jk  ) / ( ze3tn * fse3w(ji,jj,jk  ) )
222                     zws(ji,jj,jk) = - p2dt(jk) * zwt(ji,jj,jk+1) / ( ze3tn * fse3w(ji,jj,jk+1) )
223                     zwd(ji,jj,jk) = ze3ta - zwi(ji,jj,jk) - zws(ji,jj,jk)
224                 END DO
225               END DO
[3]226            END DO
[3211]227#endif
[2528]228            !
[2602]229            !! Matrix inversion from the first level
230            !!----------------------------------------------------------------------
231            !   solve m.x = y  where m is a tri diagonal matrix ( jpk*jpk )
232            !
233            !        ( zwd1 zws1   0    0    0  )( zwx1 ) ( zwy1 )
234            !        ( zwi2 zwd2 zws2   0    0  )( zwx2 ) ( zwy2 )
235            !        (  0   zwi3 zwd3 zws3   0  )( zwx3 )=( zwy3 )
236            !        (        ...               )( ...  ) ( ...  )
237            !        (  0    0    0   zwik zwdk )( zwxk ) ( zwyk )
238            !
239            !   m is decomposed in the product of an upper and lower triangular matrix.
240            !   The 3 diagonal terms are in 3d arrays: zwd, zws, zwi.
241            !   Suffices i,s and d indicate "inferior" (below diagonal), diagonal
242            !   and "superior" (above diagonal) components of the tridiagonal system.
243            !   The solution will be in the 4d array pta.
244            !   The 3d array zwt is used as a work space array.
245            !   En route to the solution pta is used a to evaluate the rhs and then
246            !   used as a work space array: its value is modified.
247            !
[2528]248            ! first recurrence:   Tk = Dk - Ik Sk-1 / Tk-1   (increasing k)
[2602]249            ! done once for all passive tracers (so included in the IF instruction)
[3211]250#if defined key_z_first
251                  zwt(ji,jj,1) = zwd(ji,jj,1)
[4428]252                  DO jk = 2, mbkmax(ji,jj)-1 ! jpkm1
[3211]253                    zwt(ji,jj,jk) = zwd(ji,jj,jk) - zwi(ji,jj,jk) * zws(ji,jj,jk-1) / zwt(ji,jj,jk-1)
254                  END DO
255               END DO
256            END DO
257#else
[2528]258            DO jj = 2, jpjm1
259               DO ji = fs_2, fs_jpim1
260                  zwt(ji,jj,1) = zwd(ji,jj,1)
261               END DO
[457]262            END DO
[2528]263            DO jk = 2, jpkm1
264               DO jj = 2, jpjm1
265                  DO ji = fs_2, fs_jpim1
266                    zwt(ji,jj,jk) = zwd(ji,jj,jk) - zwi(ji,jj,jk) * zws(ji,jj,jk-1) / zwt(ji,jj,jk-1)
267                  END DO
268               END DO
269            END DO
[3211]270#endif
[2528]271            !
272         END IF 
[2602]273         !         
[2528]274         ! second recurrence:    Zk = Yk - Ik / Tk-1  Zk-1
[3211]275#if defined key_z_first
[457]276         DO jj = 2, jpjm1
[3211]277            DO ji = 2, jpim1
278               ze3tb = ( 1. - r_vvl ) + r_vvl * fse3t_b(ji,jj,1)
279               ze3tn = ( 1. - r_vvl ) + r_vvl * fse3t(ji,jj,1)
280               pta(ji,jj,1,jn) = ze3tb * ptb(ji,jj,1,jn) + p2dt(1) * ze3tn * pta(ji,jj,1,jn)
[4428]281               DO jk = 2, mbkmax(ji,jj)-1 ! jpkm1
[3211]282                  ze3tb = ( 1. - r_vvl ) + r_vvl * fse3t_b(ji,jj,jk)
283                  ze3tn = ( 1. - r_vvl ) + r_vvl * fse3t  (ji,jj,jk)
284                  zrhs = ze3tb * ptb(ji,jj,jk,jn) + p2dt(jk) * ze3tn * pta(ji,jj,jk,jn)   ! zrhs=right hand side
285                  pta(ji,jj,jk,jn) = zrhs - zwi(ji,jj,jk) / zwt(ji,jj,jk-1) * pta(ji,jj,jk-1,jn)
286               END DO
287#else
288         DO jj = 2, jpjm1
[457]289            DO ji = fs_2, fs_jpim1
[2602]290               ze3tb = ( 1. - r_vvl ) + r_vvl * fse3t_b(ji,jj,1)
291               ze3tn = ( 1. - r_vvl ) + r_vvl * fse3t(ji,jj,1)
[2528]292               pta(ji,jj,1,jn) = ze3tb * ptb(ji,jj,1,jn) + p2dt(1) * ze3tn * pta(ji,jj,1,jn)
[457]293            END DO
294         END DO
[2528]295         DO jk = 2, jpkm1
296            DO jj = 2, jpjm1
297               DO ji = fs_2, fs_jpim1
[2602]298                  ze3tb = ( 1. - r_vvl ) + r_vvl * fse3t_b(ji,jj,jk)
299                  ze3tn = ( 1. - r_vvl ) + r_vvl * fse3t  (ji,jj,jk)
[2528]300                  zrhs = ze3tb * ptb(ji,jj,jk,jn) + p2dt(jk) * ze3tn * pta(ji,jj,jk,jn)   ! zrhs=right hand side
301                  pta(ji,jj,jk,jn) = zrhs - zwi(ji,jj,jk) / zwt(ji,jj,jk-1) * pta(ji,jj,jk-1,jn)
302               END DO
[3]303            END DO
304         END DO
[3211]305#endif
[457]306
[2602]307         ! third recurrence:    Xk = (Zk - Sk Xk+1 ) / Tk   (result is the after tracer)
[3211]308#if defined key_z_first
309               pta(ji,jj,jpkm1,jn) = pta(ji,jj,jpkm1,jn) / zwt(ji,jj,jpkm1) * tmask(ji,jj,jpkm1)
[4428]310               DO jk = mbkmax(ji,jj)-2, 1, -1 ! jpk-2
[3211]311                  pta(ji,jj,jk,jn) = ( pta(ji,jj,jk,jn) - zws(ji,jj,jk) * pta(ji,jj,jk+1,jn) )   &
312                     &             / zwt(ji,jj,jk) * tmask(ji,jj,jk)
313               END DO
314            END DO
315         END DO
316#else
[457]317         DO jj = 2, jpjm1
318            DO ji = fs_2, fs_jpim1
[2528]319               pta(ji,jj,jpkm1,jn) = pta(ji,jj,jpkm1,jn) / zwt(ji,jj,jpkm1) * tmask(ji,jj,jpkm1)
[3]320            END DO
321         END DO
[2528]322         DO jk = jpk-2, 1, -1
323            DO jj = 2, jpjm1
324               DO ji = fs_2, fs_jpim1
325                  pta(ji,jj,jk,jn) = ( pta(ji,jj,jk,jn) - zws(ji,jj,jk) * pta(ji,jj,jk+1,jn) )   &
326                     &             / zwt(ji,jj,jk) * tmask(ji,jj,jk)
327               END DO
[457]328            END DO
329         END DO
[3211]330#endif
[2602]331         !                                            ! ================= !
332      END DO                                          !  end tracer loop  !
333      !                                               ! ================= !
[1438]334      !
[2715]335      IF( wrk_not_released(3, 6,7) )   CALL ctl_stop('tra_zdf_imp: failed to release workspace arrays')
336      !
[3432]337      CALL timing_stop('tra_zdf_imp','section')
338      !
[3]339   END SUBROUTINE tra_zdf_imp
340
341   !!==============================================================================
342END MODULE trazdf_imp
Note: See TracBrowser for help on using the repository browser.