1 | MODULE 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 | !! 3.7 ! 2015-11 (G. Madec, A. Coward) non linear free surface by default |
---|
19 | !!---------------------------------------------------------------------- |
---|
20 | |
---|
21 | !!---------------------------------------------------------------------- |
---|
22 | !! tra_zdf_imp : Update the tracer trend with vertical mixing, nad compute the after tracer field |
---|
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 ! 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 |
---|
39 | |
---|
40 | IMPLICIT NONE |
---|
41 | PRIVATE |
---|
42 | |
---|
43 | PUBLIC tra_zdf_imp ! routine called by step.F90 |
---|
44 | |
---|
45 | !! * Substitutions |
---|
46 | # include "zdfddm_substitute.h90" |
---|
47 | # include "vectopt_loop_substitute.h90" |
---|
48 | !!---------------------------------------------------------------------- |
---|
49 | !! NEMO/OPA 3.7 , NEMO Consortium (2015) |
---|
50 | !! $Id$ |
---|
51 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
52 | !!---------------------------------------------------------------------- |
---|
53 | CONTAINS |
---|
54 | |
---|
55 | SUBROUTINE tra_zdf_imp( kt, kit000, cdtype, p2dt, ptb, pta, kjpt ) |
---|
56 | !!---------------------------------------------------------------------- |
---|
57 | !! *** ROUTINE tra_zdf_imp *** |
---|
58 | !! |
---|
59 | !! ** Purpose : Compute the after tracer through a implicit computation |
---|
60 | !! of the vertical tracer diffusion (including the vertical component |
---|
61 | !! of lateral mixing (only for 2nd order operator, for fourth order |
---|
62 | !! it is already computed and add to the general trend in traldf) |
---|
63 | !! |
---|
64 | !! ** Method : The vertical diffusion of a tracer ,t , is given by: |
---|
65 | !! difft = dz( avt dz(t) ) = 1/e3t dk+1( avt/e3w dk(t) ) |
---|
66 | !! It is computed using a backward time scheme (t=after field) |
---|
67 | !! which provide directly the after tracer field. |
---|
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) , INTENT(in ) :: p2dt ! 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 ! in: tracer trend ; out: after tracer field |
---|
82 | ! |
---|
83 | INTEGER :: ji, jj, jk, jn ! dummy loop indices |
---|
84 | REAL(wp) :: zrhs ! 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 | ENDIF |
---|
97 | ! ! ============= ! |
---|
98 | DO jn = 1, kjpt ! tracer loop ! |
---|
99 | ! ! ============= ! |
---|
100 | ! Matrix construction |
---|
101 | ! -------------------- |
---|
102 | ! Build matrix if temperature or salinity (only in double diffusion case) or first passive tracer |
---|
103 | ! |
---|
104 | IF( ( cdtype == 'TRA' .AND. ( jn == jp_tem .OR. ( jn == jp_sal .AND. lk_zdfddm ) ) ) .OR. & |
---|
105 | & ( cdtype == 'TRC' .AND. jn == 1 ) ) THEN |
---|
106 | ! |
---|
107 | ! vertical mixing coef.: avt for temperature, avs for salinity and passive tracers |
---|
108 | IF( cdtype == 'TRA' .AND. jn == jp_tem ) THEN ; zwt(:,:,2:jpk) = avt (:,:,2:jpk) |
---|
109 | ELSE ; zwt(:,:,2:jpk) = fsavs(:,:,2:jpk) |
---|
110 | ENDIF |
---|
111 | zwt(:,:,1) = 0._wp |
---|
112 | ! |
---|
113 | IF( l_ldfslp ) THEN ! isoneutral diffusion: add the contribution |
---|
114 | IF( ln_traldf_msc ) THEN ! MSC iso-neutral operator |
---|
115 | DO jk = 2, jpkm1 |
---|
116 | DO jj = 2, jpjm1 |
---|
117 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
118 | zwt(ji,jj,jk) = zwt(ji,jj,jk) + akz(ji,jj,jk) |
---|
119 | END DO |
---|
120 | END DO |
---|
121 | END DO |
---|
122 | ELSE ! standard or triad iso-neutral operator |
---|
123 | DO jk = 2, jpkm1 |
---|
124 | DO jj = 2, jpjm1 |
---|
125 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
126 | zwt(ji,jj,jk) = zwt(ji,jj,jk) + ah_wslp2(ji,jj,jk) |
---|
127 | END DO |
---|
128 | END DO |
---|
129 | END DO |
---|
130 | ENDIF |
---|
131 | ENDIF |
---|
132 | ! |
---|
133 | ! Diagonal, lower (i), upper (s) (including the bottom boundary condition since avt is masked) |
---|
134 | DO jk = 1, jpkm1 |
---|
135 | DO jj = 2, jpjm1 |
---|
136 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
137 | !!gm BUG I think, use e3w_a instead of e3w_n |
---|
138 | zwi(ji,jj,jk) = - p2dt * zwt(ji,jj,jk ) / e3w_n(ji,jj,jk ) |
---|
139 | zws(ji,jj,jk) = - p2dt * zwt(ji,jj,jk+1) / e3w_n(ji,jj,jk+1) |
---|
140 | zwd(ji,jj,jk) = e3t_a(ji,jj,jk) - zwi(ji,jj,jk) - zws(ji,jj,jk) |
---|
141 | END DO |
---|
142 | END DO |
---|
143 | END DO |
---|
144 | ! |
---|
145 | !! Matrix inversion from the first level |
---|
146 | !!---------------------------------------------------------------------- |
---|
147 | ! solve m.x = y where m is a tri diagonal matrix ( jpk*jpk ) |
---|
148 | ! |
---|
149 | ! ( zwd1 zws1 0 0 0 )( zwx1 ) ( zwy1 ) |
---|
150 | ! ( zwi2 zwd2 zws2 0 0 )( zwx2 ) ( zwy2 ) |
---|
151 | ! ( 0 zwi3 zwd3 zws3 0 )( zwx3 )=( zwy3 ) |
---|
152 | ! ( ... )( ... ) ( ... ) |
---|
153 | ! ( 0 0 0 zwik zwdk )( zwxk ) ( zwyk ) |
---|
154 | ! |
---|
155 | ! m is decomposed in the product of an upper and lower triangular matrix. |
---|
156 | ! The 3 diagonal terms are in 3d arrays: zwd, zws, zwi. |
---|
157 | ! Suffices i,s and d indicate "inferior" (below diagonal), diagonal |
---|
158 | ! and "superior" (above diagonal) components of the tridiagonal system. |
---|
159 | ! The solution will be in the 4d array pta. |
---|
160 | ! The 3d array zwt is used as a work space array. |
---|
161 | ! En route to the solution pta is used a to evaluate the rhs and then |
---|
162 | ! used as a work space array: its value is modified. |
---|
163 | ! |
---|
164 | DO jj = 2, jpjm1 !* 1st recurrence: Tk = Dk - Ik Sk-1 / Tk-1 (increasing k) |
---|
165 | DO ji = fs_2, fs_jpim1 ! done one for all passive tracers (so included in the IF instruction) |
---|
166 | zwt(ji,jj,1) = zwd(ji,jj,1) |
---|
167 | END DO |
---|
168 | END DO |
---|
169 | DO jk = 2, jpkm1 |
---|
170 | DO jj = 2, jpjm1 |
---|
171 | DO ji = fs_2, fs_jpim1 |
---|
172 | zwt(ji,jj,jk) = zwd(ji,jj,jk) - zwi(ji,jj,jk) * zws(ji,jj,jk-1) / zwt(ji,jj,jk-1) |
---|
173 | END DO |
---|
174 | END DO |
---|
175 | END DO |
---|
176 | ! |
---|
177 | ENDIF |
---|
178 | ! |
---|
179 | DO jj = 2, jpjm1 !* 2nd recurrence: Zk = Yk - Ik / Tk-1 Zk-1 |
---|
180 | DO ji = fs_2, fs_jpim1 |
---|
181 | 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) |
---|
182 | END DO |
---|
183 | END DO |
---|
184 | DO jk = 2, jpkm1 |
---|
185 | DO jj = 2, jpjm1 |
---|
186 | DO ji = fs_2, fs_jpim1 |
---|
187 | 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 |
---|
188 | pta(ji,jj,jk,jn) = zrhs - zwi(ji,jj,jk) / zwt(ji,jj,jk-1) * pta(ji,jj,jk-1,jn) |
---|
189 | END DO |
---|
190 | END DO |
---|
191 | END DO |
---|
192 | ! |
---|
193 | DO jj = 2, jpjm1 !* 3d recurrence: Xk = (Zk - Sk Xk+1 ) / Tk (result is the after tracer) |
---|
194 | DO ji = fs_2, fs_jpim1 |
---|
195 | pta(ji,jj,jpkm1,jn) = pta(ji,jj,jpkm1,jn) / zwt(ji,jj,jpkm1) * tmask(ji,jj,jpkm1) |
---|
196 | END DO |
---|
197 | END DO |
---|
198 | DO jk = jpk-2, 1, -1 |
---|
199 | DO jj = 2, jpjm1 |
---|
200 | DO ji = fs_2, fs_jpim1 |
---|
201 | pta(ji,jj,jk,jn) = ( pta(ji,jj,jk,jn) - zws(ji,jj,jk) * pta(ji,jj,jk+1,jn) ) & |
---|
202 | & / zwt(ji,jj,jk) * tmask(ji,jj,jk) |
---|
203 | END DO |
---|
204 | END DO |
---|
205 | END DO |
---|
206 | ! ! ================= ! |
---|
207 | END DO ! end tracer loop ! |
---|
208 | ! ! ================= ! |
---|
209 | ! |
---|
210 | CALL wrk_dealloc( jpi,jpj,jpk, zwi, zwt, zwd, zws ) |
---|
211 | ! |
---|
212 | IF( nn_timing == 1 ) CALL timing_stop('tra_zdf_imp') |
---|
213 | ! |
---|
214 | END SUBROUTINE tra_zdf_imp |
---|
215 | |
---|
216 | !!============================================================================== |
---|
217 | END MODULE trazdf_imp |
---|