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_crs.F90 in branches/2015/dev_r5003_MERCATOR6_CRS/NEMOGCM/NEMO/OPA_SRC/TRA – NEMO

source: branches/2015/dev_r5003_MERCATOR6_CRS/NEMOGCM/NEMO/OPA_SRC/TRA/trazdf_imp_crs.F90 @ 6772

Last change on this file since 6772 was 6772, checked in by cbricaud, 8 years ago

clean in coarsening branch

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