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

source: branches/2017/dev_r7881_HPC09_ZDF/NEMOGCM/NEMO/OPA_SRC/TRA/trazdf_imp.F90 @ 7953

Last change on this file since 7953 was 7931, checked in by gm, 7 years ago

#1880 (HPC-09): remove key_zdfddm + phasing with last changes of HPC08 branch

  • Property svn:keywords set to Id
File size: 11.1 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
[6140]18   !!            3.7  !  2015-11  (G. Madec, A. Coward)  non linear free surface by default
[3]19   !!----------------------------------------------------------------------
[1438]20 
21   !!----------------------------------------------------------------------
[6140]22   !!   tra_zdf_imp   : Update the tracer trend with vertical mixing, nad compute the after tracer field
[3]23   !!----------------------------------------------------------------------
[5836]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
27   USE trc_oce        ! share passive tracers/ocean variables
28   USE domvvl         ! variable volume
29   USE ldftra         ! lateral mixing type
30   USE ldfslp         ! lateral physics: slope of diffusion
31   USE zdfddm         ! ocean vertical physics: double diffusion
32   USE traldf_triad   ! active tracers: Method of Stabilizing Correction
33   !
34   USE in_out_manager ! I/O manager
35   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
36   USE lib_mpp        ! MPP library
37   USE wrk_nemo       ! Memory Allocation
38   USE timing         ! Timing
[3]39
40   IMPLICIT NONE
41   PRIVATE
42
[1438]43   PUBLIC   tra_zdf_imp   !  routine called by step.F90
[3]44
45   !! * Substitutions
[457]46#  include "vectopt_loop_substitute.h90"
[3]47   !!----------------------------------------------------------------------
[5836]48   !! NEMO/OPA 3.7 , NEMO Consortium (2015)
[1156]49   !! $Id$
[2528]50   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
[457]51   !!----------------------------------------------------------------------
[3]52CONTAINS
[2528]53 
[3294]54   SUBROUTINE tra_zdf_imp( kt, kit000, cdtype, p2dt, ptb, pta, kjpt ) 
[3]55      !!----------------------------------------------------------------------
56      !!                  ***  ROUTINE tra_zdf_imp  ***
57      !!
[2602]58      !! ** Purpose :   Compute the after tracer through a implicit computation
59      !!     of the vertical tracer diffusion (including the vertical component
60      !!     of lateral mixing (only for 2nd order operator, for fourth order
61      !!     it is already computed and add to the general trend in traldf)
[3]62      !!
[6140]63      !! ** Method  :  The vertical diffusion of a tracer ,t , is given by:
64      !!          difft = dz( avt dz(t) ) = 1/e3t dk+1( avt/e3w dk(t) )
65      !!      It is computed using a backward time scheme (t=after field)
66      !!      which provide directly the after tracer field.
[7931]67      !!      If ln_zdfddm=T, use avs for salinity or for passive tracers
[3]68      !!      Surface and bottom boundary conditions: no diffusive flux on
69      !!      both tracers (bottom, applied through the masked field avt).
[2602]70      !!      If iso-neutral mixing, add to avt the contribution due to lateral mixing.
[3]71      !!
[2602]72      !! ** Action  : - pta  becomes the after tracer
[3]73      !!---------------------------------------------------------------------
[2528]74      INTEGER                              , INTENT(in   ) ::   kt       ! ocean time-step index
[6140]75      INTEGER                              , INTENT(in   ) ::   kit000   ! first time step index
[2528]76      CHARACTER(len=3)                     , INTENT(in   ) ::   cdtype   ! =TRA or TRC (tracer indicator)
77      INTEGER                              , INTENT(in   ) ::   kjpt     ! number of tracers
[6140]78      REAL(wp)                             , INTENT(in   ) ::   p2dt     ! tracer time-step
[2528]79      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(in   ) ::   ptb      ! before and now tracer fields
[6140]80      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::   pta      ! in: tracer trend ; out: after tracer field
[2715]81      !
82      INTEGER  ::  ji, jj, jk, jn   ! dummy loop indices
[6140]83      REAL(wp) ::  zrhs             ! local scalars
[5836]84      REAL(wp), POINTER, DIMENSION(:,:,:) ::  zwi, zwt, zwd, zws
[3]85      !!---------------------------------------------------------------------
[3294]86      !
87      IF( nn_timing == 1 )  CALL timing_start('tra_zdf_imp')
88      !
[5836]89      CALL wrk_alloc( jpi,jpj,jpk,   zwi, zwt, zwd, zws ) 
[3294]90      !
91      IF( kt == kit000 )  THEN
[457]92         IF(lwp)WRITE(numout,*)
[2528]93         IF(lwp)WRITE(numout,*) 'tra_zdf_imp : implicit vertical mixing on ', cdtype
[457]94         IF(lwp)WRITE(numout,*) '~~~~~~~~~~~ '
95      ENDIF
[2602]96      !                                               ! ============= !
97      DO jn = 1, kjpt                                 !  tracer loop  !
98         !                                            ! ============= !
[2528]99         !  Matrix construction
[2602]100         ! --------------------
101         ! Build matrix if temperature or salinity (only in double diffusion case) or first passive tracer
102         !
[7931]103         IF(  ( cdtype == 'TRA' .AND. ( jn == jp_tem .OR. ( jn == jp_sal .AND. ln_zdfddm ) ) ) .OR.   &
[2602]104            & ( cdtype == 'TRC' .AND. jn == 1 )  )  THEN
105            !
106            ! vertical mixing coef.: avt for temperature, avs for salinity and passive tracers
[7931]107            IF( cdtype == 'TRA' .AND. jn == jp_tem ) THEN   ;   zwt(:,:,2:jpk) = avt(:,:,2:jpk)
108            ELSE                                            ;   zwt(:,:,2:jpk) = avs(:,:,2:jpk)
[2528]109            ENDIF
[7753]110            zwt(:,:,1) = 0._wp
[5836]111            !
112            IF( l_ldfslp ) THEN            ! isoneutral diffusion: add the contribution
113               IF( ln_traldf_msc  ) THEN     ! MSC iso-neutral operator
114                  DO jk = 2, jpkm1
115                     DO jj = 2, jpjm1
116                        DO ji = fs_2, fs_jpim1   ! vector opt.
117                           zwt(ji,jj,jk) = zwt(ji,jj,jk) + akz(ji,jj,jk) 
118                        END DO
[2528]119                     END DO
120                  END DO
[5836]121               ELSE                          ! standard or triad iso-neutral operator
122                  DO jk = 2, jpkm1
123                     DO jj = 2, jpjm1
124                        DO ji = fs_2, fs_jpim1   ! vector opt.
125                           zwt(ji,jj,jk) = zwt(ji,jj,jk) + ah_wslp2(ji,jj,jk)
126                        END DO
[2528]127                     END DO
128                  END DO
[5836]129               ENDIF
[2528]130            ENDIF
[5836]131            !
[2602]132            ! Diagonal, lower (i), upper (s)  (including the bottom boundary condition since avt is masked)
[2528]133            DO jk = 1, jpkm1
134               DO jj = 2, jpjm1
135                  DO ji = fs_2, fs_jpim1   ! vector opt.
[6140]136!!gm BUG  I think, use e3w_a instead of e3w_n
137                     zwi(ji,jj,jk) = - p2dt * zwt(ji,jj,jk  ) / e3w_n(ji,jj,jk  )
138                     zws(ji,jj,jk) = - p2dt * zwt(ji,jj,jk+1) / e3w_n(ji,jj,jk+1)
139                     zwd(ji,jj,jk) = e3t_a(ji,jj,jk) - zwi(ji,jj,jk) - zws(ji,jj,jk)
[2528]140                 END DO
141               END DO
[3]142            END DO
[2528]143            !
[2602]144            !! Matrix inversion from the first level
145            !!----------------------------------------------------------------------
146            !   solve m.x = y  where m is a tri diagonal matrix ( jpk*jpk )
147            !
148            !        ( zwd1 zws1   0    0    0  )( zwx1 ) ( zwy1 )
149            !        ( zwi2 zwd2 zws2   0    0  )( zwx2 ) ( zwy2 )
150            !        (  0   zwi3 zwd3 zws3   0  )( zwx3 )=( zwy3 )
151            !        (        ...               )( ...  ) ( ...  )
152            !        (  0    0    0   zwik zwdk )( zwxk ) ( zwyk )
153            !
154            !   m is decomposed in the product of an upper and lower triangular matrix.
155            !   The 3 diagonal terms are in 3d arrays: zwd, zws, zwi.
156            !   Suffices i,s and d indicate "inferior" (below diagonal), diagonal
157            !   and "superior" (above diagonal) components of the tridiagonal system.
158            !   The solution will be in the 4d array pta.
159            !   The 3d array zwt is used as a work space array.
160            !   En route to the solution pta is used a to evaluate the rhs and then
161            !   used as a work space array: its value is modified.
162            !
[6140]163            DO jj = 2, jpjm1        !* 1st recurrence:   Tk = Dk - Ik Sk-1 / Tk-1   (increasing k)
164               DO ji = fs_2, fs_jpim1            ! done one for all passive tracers (so included in the IF instruction)
[5120]165                  zwt(ji,jj,1) = zwd(ji,jj,1)
166               END DO
167            END DO
168            DO jk = 2, jpkm1
169               DO jj = 2, jpjm1
170                  DO ji = fs_2, fs_jpim1
[4990]171                     zwt(ji,jj,jk) = zwd(ji,jj,jk) - zwi(ji,jj,jk) * zws(ji,jj,jk-1) / zwt(ji,jj,jk-1)
[2528]172                  END DO
173               END DO
174            END DO
175            !
[6140]176         ENDIF 
[2602]177         !         
[6140]178         DO jj = 2, jpjm1           !* 2nd recurrence:    Zk = Yk - Ik / Tk-1  Zk-1
[457]179            DO ji = fs_2, fs_jpim1
[6140]180               pta(ji,jj,1,jn) = e3t_b(ji,jj,1) * ptb(ji,jj,1,jn) + p2dt * e3t_n(ji,jj,1) * pta(ji,jj,1,jn)
[5120]181            END DO
182         END DO
183         DO jk = 2, jpkm1
184            DO jj = 2, jpjm1
185               DO ji = fs_2, fs_jpim1
[6140]186                  zrhs = e3t_b(ji,jj,jk) * ptb(ji,jj,jk,jn) + p2dt * e3t_n(ji,jj,jk) * pta(ji,jj,jk,jn)   ! zrhs=right hand side
[2528]187                  pta(ji,jj,jk,jn) = zrhs - zwi(ji,jj,jk) / zwt(ji,jj,jk-1) * pta(ji,jj,jk-1,jn)
188               END DO
[3]189            END DO
190         END DO
[6140]191         !
192         DO jj = 2, jpjm1           !* 3d recurrence:    Xk = (Zk - Sk Xk+1 ) / Tk   (result is the after tracer)
[457]193            DO ji = fs_2, fs_jpim1
[2528]194               pta(ji,jj,jpkm1,jn) = pta(ji,jj,jpkm1,jn) / zwt(ji,jj,jpkm1) * tmask(ji,jj,jpkm1)
[5120]195            END DO
196         END DO
197         DO jk = jpk-2, 1, -1
198            DO jj = 2, jpjm1
199               DO ji = fs_2, fs_jpim1
[2528]200                  pta(ji,jj,jk,jn) = ( pta(ji,jj,jk,jn) - zws(ji,jj,jk) * pta(ji,jj,jk+1,jn) )   &
201                     &             / zwt(ji,jj,jk) * tmask(ji,jj,jk)
202               END DO
[457]203            END DO
204         END DO
[2602]205         !                                            ! ================= !
206      END DO                                          !  end tracer loop  !
207      !                                               ! ================= !
[1438]208      !
[5836]209      CALL wrk_dealloc( jpi,jpj,jpk,   zwi, zwt, zwd, zws ) 
[2715]210      !
[3294]211      IF( nn_timing == 1 )  CALL timing_stop('tra_zdf_imp')
212      !
[3]213   END SUBROUTINE tra_zdf_imp
214
215   !!==============================================================================
216END MODULE trazdf_imp
Note: See TracBrowser for help on using the repository browser.