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

source: branches/dev_r2586_dynamic_mem/NEMOGCM/NEMO/OPA_SRC/TRA/trazdf_imp.F90 @ 2696

Last change on this file since 2696 was 2696, checked in by rblod, 13 years ago

Compatibility with key_traldf_ano for dynamic memory

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