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

source: branches/UKMO/dev_r5107_iceshelf_fw_input_coupled_model/NEMOGCM/NEMO/OPA_SRC/TRA/trazdf_imp.F90 @ 5511

Last change on this file since 5511 was 5511, checked in by davestorkey, 9 years ago

UKMO/dev_r5107_iceshelf_fw_input_coupled_model branch: clear SVN keywords

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