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/2015/dev_r5836_NOC3_vvl_by_default/NEMOGCM/NEMO/OPA_SRC/TRA – NEMO

source: branches/2015/dev_r5836_NOC3_vvl_by_default/NEMOGCM/NEMO/OPA_SRC/TRA/trazdf_imp.F90 @ 5866

Last change on this file since 5866 was 5866, checked in by gm, 8 years ago

#1613: vvl by default: add ln_linssh and remove key_vvl

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