/[MITgcm]/MITgcm/model/src/impldiff.F
ViewVC logotype

Contents of /MITgcm/model/src/impldiff.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.36 - (show annotations) (download)
Thu Oct 6 14:22:12 2016 UTC (7 years, 6 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, HEAD
Changes since 1.35: +5 -3 lines
- add parenthesis in calulation of main diagonal of implicit vertical
  diffusision matrix. This affects results at machine truncation level.

1 C $Header: /u/gcmpack/MITgcm/model/src/impldiff.F,v 1.35 2015/06/03 13:39:22 rpa Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6
7 CBOP
8 C !ROUTINE: IMPLDIFF
9 C !INTERFACE:
10 SUBROUTINE IMPLDIFF( bi, bj, iMin, iMax, jMin, jMax,
11 I tracerId, KappaRX, recip_hFac,
12 U gTracer,
13 I myThid )
14 C !DESCRIPTION: \bv
15 C *==========================================================*
16 C | S/R IMPLDIFF
17 C | o Solve implicit diffusion equation for vertical
18 C | diffusivity.
19 C *==========================================================*
20 C | o Recoded from 2d intermediate fields to 3d to reduce
21 C | TAMC storage
22 C | o Fixed missing masks for fields a(), c()
23 C *==========================================================*
24 C \ev
25
26 C !USES:
27 IMPLICIT NONE
28 C == Global data ==
29 #include "SIZE.h"
30 #include "DYNVARS.h"
31 #include "EEPARAMS.h"
32 #include "PARAMS.h"
33 #include "GRID.h"
34 #ifdef ALLOW_GENERIC_ADVDIFF
35 # include "GAD.h"
36 #endif
37 #ifdef ALLOW_LONGSTEP
38 # include "LONGSTEP_PARAMS.h"
39 #endif
40 #ifdef ALLOW_PTRACERS
41 # include "PTRACERS_SIZE.h"
42 # include "PTRACERS_PARAMS.h"
43 #endif
44 c#ifdef ALLOW_AUTODIFF_TAMC
45 c# include "tamc_keys.h"
46 c#endif
47
48 C !INPUT/OUTPUT PARAMETERS:
49 C == Routine Arguments ==
50 C tracerId :: tracer Identificator (if > 0) ; = -1 or -2 when
51 C solving vertical viscosity implicitly for U or V
52 C KappaRk :: vertical diffusion coefficient
53 C recip_hFac :: Inverse of cell open-depth factor
54 C gTracer :: future tracer field
55 INTEGER bi,bj,iMin,iMax,jMin,jMax
56 INTEGER tracerId
57 _RL KappaRX(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
58 _RS recip_hFac(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
59 _RL gTracer(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
60 INTEGER myThid
61
62 C !LOCAL VARIABLES:
63 C == Local variables ==
64 INTEGER i,j,k
65 _RL deltaTX(Nr)
66 _RL locTr(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
67 _RL a(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
68 _RL b(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
69 _RL c(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
70 _RL bet(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
71 _RL gam(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
72 #ifdef ALLOW_DIAGNOSTICS
73 CHARACTER*8 diagName
74 CHARACTER*4 diagSufx
75 #ifdef ALLOW_GENERIC_ADVDIFF
76 CHARACTER*4 GAD_DIAG_SUFX
77 EXTERNAL GAD_DIAG_SUFX
78 #endif
79 LOGICAL DIAGNOSTICS_IS_ON
80 EXTERNAL DIAGNOSTICS_IS_ON
81 _RL df (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
82 #endif /* ALLOW_DIAGNOSTICS */
83 CEOP
84
85 cph(
86 cph Not good for TAF: may create irreducible control flow graph
87 cph IF (Nr.LE.1) RETURN
88 cph)
89
90 #ifdef ALLOW_PTRACERS
91 IF ( tracerId.GE.GAD_TR1) THEN
92 DO k=1,Nr
93 deltaTX(k) = PTRACERS_dTLev(k)
94 ENDDO
95 ELSEIF ( tracerId.GE.1 ) THEN
96 #else
97 IF ( tracerId.GE.1 ) THEN
98 #endif
99 DO k=1,Nr
100 deltaTX(k) = dTtracerLev(k)
101 ENDDO
102 ELSE
103 DO k=1,Nr
104 deltaTX(k) = deltaTMom
105 ENDDO
106 ENDIF
107
108 C-- Initialise
109 DO k=1,Nr
110 #ifdef TARGET_NEC_SX
111 DO j=1-OLy,sNy+OLy
112 DO i=1-OLx,sNx+OLx
113 #else
114 DO j=jMin,jMax
115 DO i=iMin,iMax
116 #endif
117 locTr(i,j,k) = 0. _d 0
118 ENDDO
119 ENDDO
120 ENDDO
121
122 C-- Old aLower
123 #ifdef TARGET_NEC_SX
124 DO j=1-OLy,sNy+OLy
125 DO i=1-OLx,sNx+OLx
126 #else
127 DO j=jMin,jMax
128 DO i=iMin,iMax
129 #endif
130 a(i,j,1) = 0. _d 0
131 ENDDO
132 ENDDO
133 DO k=2,Nr
134 #ifdef TARGET_NEC_SX
135 DO j=1-OLy,sNy+OLy
136 DO i=1-OLx,sNx+OLx
137 #else
138 DO j=jMin,jMax
139 DO i=iMin,iMax
140 #endif
141 a(i,j,k) = -deltaTX(k)*recip_hFac(i,j,k)*recip_drF(k)
142 & *recip_deepFac2C(k)*recip_rhoFacC(k)
143 & *KappaRX(i,j, k )*recip_drC( k )
144 & *deepFac2F(k)*rhoFacF(k)
145 IF (recip_hFac(i,j,k-1).EQ.0.) a(i,j,k)=0.
146 ENDDO
147 ENDDO
148 ENDDO
149
150 C-- Old aUpper
151 DO k=1,Nr-1
152 #ifdef TARGET_NEC_SX
153 DO j=1-OLy,sNy+OLy
154 DO i=1-OLx,sNx+OLx
155 #else
156 DO j=jMin,jMax
157 DO i=iMin,iMax
158 #endif
159 c(i,j,k) = -deltaTX(k)*recip_hFac(i,j,k)*recip_drF(k)
160 & *recip_deepFac2C(k)*recip_rhoFacC(k)
161 & *KappaRX(i,j,k+1)*recip_drC(k+1)
162 & *deepFac2F(k+1)*rhoFacF(k+1)
163 IF (recip_hFac(i,j,k+1).EQ.0.) c(i,j,k)=0.
164 ENDDO
165 ENDDO
166 ENDDO
167 #ifdef TARGET_NEC_SX
168 DO j=1-OLy,sNy+OLy
169 DO i=1-OLx,sNx+OLx
170 #else
171 DO j=jMin,jMax
172 DO i=iMin,iMax
173 #endif
174 c(i,j,Nr) = 0. _d 0
175 ENDDO
176 ENDDO
177
178 C-- Old aCenter
179 DO k=1,Nr
180 #ifdef TARGET_NEC_SX
181 DO j=1-OLy,sNy+OLy
182 DO i=1-OLx,sNx+OLx
183 #else
184 DO j=jMin,jMax
185 DO i=iMin,iMax
186 #endif
187 b(i,j,k) = 1. _d 0 - ( a(i,j,k) + c(i,j,k) )
188 C- to recover older (prior to 2016-10-05) results:
189 c b(i,j,k) = 1. _d 0 - c(i,j,k) - a(i,j,k)
190 ENDDO
191 ENDDO
192 ENDDO
193
194 C-- Old and new gam, bet are the same
195 DO k=1,Nr
196 #ifdef TARGET_NEC_SX
197 DO j=1-OLy,sNy+OLy
198 DO i=1-OLx,sNx+OLx
199 #else
200 DO j=jMin,jMax
201 DO i=iMin,iMax
202 #endif
203 bet(i,j,k) = 1. _d 0
204 gam(i,j,k) = 0. _d 0
205 ENDDO
206 ENDDO
207 ENDDO
208
209 C-- Only need do anything if Nr>1
210 IF (Nr.GT.1) THEN
211
212 k = 1
213 C-- Beginning of forward sweep (top level)
214 #ifdef TARGET_NEC_SX
215 DO j=1-OLy,sNy+OLy
216 DO i=1-OLx,sNx+OLx
217 #else
218 DO j=jMin,jMax
219 DO i=iMin,iMax
220 #endif
221 IF (b(i,j,1).NE.0.) bet(i,j,1) = 1. _d 0 / b(i,j,1)
222 ENDDO
223 ENDDO
224
225 ENDIF
226
227 C-- Middle of forward sweep
228 IF (Nr.GE.2) THEN
229
230 CADJ loop = sequential
231 DO k=2,Nr
232
233 #ifdef TARGET_NEC_SX
234 DO j=1-OLy,sNy+OLy
235 DO i=1-OLx,sNx+OLx
236 #else
237 DO j=jMin,jMax
238 DO i=iMin,iMax
239 #endif
240 gam(i,j,k) = c(i,j,k-1)*bet(i,j,k-1)
241 IF ( ( b(i,j,k) - a(i,j,k)*gam(i,j,k) ) .NE. 0.)
242 & bet(i,j,k) = 1. _d 0 / ( b(i,j,k) - a(i,j,k)*gam(i,j,k) )
243 ENDDO
244 ENDDO
245
246 ENDDO
247
248 ENDIF
249
250 #ifdef TARGET_NEC_SX
251 DO j=1-OLy,sNy+OLy
252 DO i=1-OLx,sNx+OLx
253 #else
254 DO j=jMin,jMax
255 DO i=iMin,iMax
256 #endif
257 locTr(i,j,1) = gTracer(i,j,1)*bet(i,j,1)
258 ENDDO
259 ENDDO
260 DO k=2,Nr
261 #ifdef TARGET_NEC_SX
262 DO j=1-OLy,sNy+OLy
263 DO i=1-OLx,sNx+OLx
264 #else
265 DO j=jMin,jMax
266 DO i=iMin,iMax
267 #endif
268 locTr(i,j,k) = bet(i,j,k)*
269 & (gTracer(i,j,k) - a(i,j,k)*locTr(i,j,k-1))
270 ENDDO
271 ENDDO
272 ENDDO
273
274 C-- Backward sweep
275 CADJ loop = sequential
276 DO k=Nr-1,1,-1
277 #ifdef TARGET_NEC_SX
278 DO j=1-OLy,sNy+OLy
279 DO i=1-OLx,sNx+OLx
280 #else
281 DO j=jMin,jMax
282 DO i=iMin,iMax
283 #endif
284 locTr(i,j,k) = locTr(i,j,k) - gam(i,j,k+1)*locTr(i,j,k+1)
285 ENDDO
286 ENDDO
287 ENDDO
288
289 DO k=1,Nr
290 #ifdef TARGET_NEC_SX
291 DO j=1-OLy,sNy+OLy
292 DO i=1-OLx,sNx+OLx
293 #else
294 DO j=jMin,jMax
295 DO i=iMin,iMax
296 #endif
297 gTracer(i,j,k) = locTr(i,j,k)
298 ENDDO
299 ENDDO
300 ENDDO
301
302 #ifdef ALLOW_DIAGNOSTICS
303 IF ( useDiagnostics .AND.tracerId.NE.0 ) THEN
304 IF ( tracerId.GE. 1 ) THEN
305 C-- Set diagnostic suffix for the current tracer
306 #ifdef ALLOW_GENERIC_ADVDIFF
307 diagSufx = GAD_DIAG_SUFX( tracerId, myThid )
308 #else
309 diagSufx = 'aaaa'
310 #endif
311 diagName = 'DFrI'//diagSufx
312 ELSEIF ( tracerId.EQ. -1 ) THEN
313 diagName = 'VISrI_Um'
314 ELSEIF ( tracerId.EQ. -2 ) THEN
315 diagName = 'VISrI_Vm'
316 ELSE
317 STOP 'IMPLIDIFF: should never reach this point !'
318 ENDIF
319 IF ( DIAGNOSTICS_IS_ON(diagName,myThid) ) THEN
320 DO k= 1,Nr
321 IF ( k.EQ.1 ) THEN
322 C- Note: Needs to call DIAGNOSTICS_FILL at level k=1 even if array == 0
323 C otherwise counter is not incremented !!
324 DO j=1-OLy,sNy+OLy
325 DO i=1-OLx,sNx+OLx
326 df(i,j) = 0. _d 0
327 ENDDO
328 ENDDO
329 ELSEIF ( tracerId.GE.1 ) THEN
330 #ifdef TARGET_NEC_SX
331 DO j=1-OLy,sNy+OLy
332 DO i=1-OLx,sNx+OLx
333 #else
334 DO j=1,sNy
335 DO i=1,sNx
336 #endif
337 df(i,j) =
338 & -rA(i,j,bi,bj)*deepFac2F(k)*rhoFacF(k)
339 & * KappaRX(i,j,k)*recip_drC(k)*rkSign
340 & * (gTracer(i,j,k) - gTracer(i,j,k-1))
341 & * maskC(i,j,k,bi,bj)
342 & * maskC(i,j,k-1,bi,bj)
343 ENDDO
344 ENDDO
345 ELSEIF ( tracerId.EQ.-1 ) THEN
346 #ifdef TARGET_NEC_SX
347 DO j=1-OLy,sNy+OLy
348 DO i=1-OLx,sNx+OLx
349 #else
350 DO j=1,sNy
351 DO i=1,sNx+1
352 #endif
353 df(i,j) =
354 & -rAw(i,j,bi,bj)*deepFac2F(k)*rhoFacF(k)
355 & * KappaRX(i,j,k)*recip_drC(k)*rkSign
356 & * (gTracer(i,j,k) - gTracer(i,j,k-1))
357 & * _maskW(i,j,k,bi,bj)
358 & * _maskW(i,j,k-1,bi,bj)
359 ENDDO
360 ENDDO
361 ELSEIF ( tracerId.EQ.-2 ) THEN
362 #ifdef TARGET_NEC_SX
363 DO j=1-OLy,sNy+OLy
364 DO i=1-OLx,sNx+OLx
365 #else
366 DO j=1,sNy+1
367 DO i=1,sNx
368 #endif
369 df(i,j) =
370 & -rAs(i,j,bi,bj)*deepFac2F(k)*rhoFacF(k)
371 & * KappaRX(i,j,k)*recip_drC(k)*rkSign
372 & * (gTracer(i,j,k) - gTracer(i,j,k-1))
373 & * _maskS(i,j,k,bi,bj)
374 & * _maskS(i,j,k-1,bi,bj)
375 ENDDO
376 ENDDO
377 ENDIF
378 CALL DIAGNOSTICS_FILL(df,diagName, k,1, 2,bi,bj, myThid)
379 #ifdef ALLOW_LAYERS
380 IF ( useLayers ) THEN
381 CALL LAYERS_FILL( df, tracerId, 'DFR',
382 & k, 1, 2,bi,bj, myThid )
383 ENDIF
384 #endif /* ALLOW_LAYERS */
385 ENDDO
386 ENDIF
387 ENDIF
388 #endif /* ALLOW_DIAGNOSTICS */
389
390 RETURN
391 END

  ViewVC Help
Powered by ViewVC 1.1.22