-
Notifications
You must be signed in to change notification settings - Fork 237
/
solve_pentadiagonal.F
432 lines (397 loc) · 12.6 KB
/
solve_pentadiagonal.F
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
#include "PACKAGES_CONFIG.h"
#include "CPP_OPTIONS.h"
#ifdef ALLOW_AUTODIFF
# include "AUTODIFF_OPTIONS.h"
#endif
CBOP
C !ROUTINE: SOLVE_PENTADIAGONAL
C !INTERFACE:
SUBROUTINE SOLVE_PENTADIAGONAL(
I iMin,iMax, jMin,jMax,
U a5d, b5d, c5d, d5d, e5d,
U y5d,
O errCode,
I bi, bj, myThid )
C !DESCRIPTION: \bv
C *==========================================================*
C | S/R SOLVE_PENTADIAGONAL
C | o Solve a penta-diagonal system A*X=Y (dimension Nr)
C *==========================================================*
C | o Used to solve implicitly vertical advection & diffusion
#ifdef SOLVE_DIAGONAL_LOWMEMORY
C | o Allows 2 types of call:
C | 1) with INPUT errCode < 0 (first call):
C | Solve system A*X=Y by LU decomposition and return
C | inverse matrix as output (in a5d,b5d,c5d,d5d,e5d)
C | 2) with INPUT errCode = 0 (subsequent calls):
C | Solve system A*Xp=Yp by using inverse matrix coeff
C | from first call output.
#endif /* SOLVE_DIAGONAL_LOWMEMORY */
C *==========================================================*
C \ev
C !USES:
IMPLICIT NONE
C == Global data ==
#include "SIZE.h"
#include "EEPARAMS.h"
C !INPUT/OUTPUT PARAMETERS:
C -- INPUT: --
C iMin,iMax,jMin,jMax :: computational domain
C a5d :: 2nd lower diagonal of the pentadiagonal matrix
C b5d :: 1rst lower diagonal of the pentadiagonal matrix
C c5d :: main diagonal of the pentadiagonal matrix
C d5d :: 1rst upper diagonal of the pentadiagonal matrix
C e5d :: 2nd upper diagonal of the pentadiagonal matrix
C y5d :: Y vector (R.H.S.);
C bi,bj :: tile indices
C myThid :: thread number
C -- OUTPUT: --
C y5d :: X = solution of A*X=Y
C errCode :: > 0 if singular matrix
#ifdef SOLVE_DIAGONAL_LOWMEMORY
C a5d,b5d,c5d,d5d,e5d :: inverse matrix coeff to enable to find Xp solution
C of A*Xp=Yp with subsequent call to this routine
#endif /* SOLVE_DIAGONAL_LOWMEMORY */
INTEGER iMin,iMax,jMin,jMax
_RL a5d(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
_RL b5d(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
_RL c5d(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
_RL d5d(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
_RL e5d(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
_RL y5d(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
INTEGER errCode
INTEGER bi, bj, myThid
C/* This flag combination is only false if INCLUDE_IMPLVERTADV_CODE is
C not defined while ALLOW_AUTODIFF is. */
#if ( defined INCLUDE_IMPLVERTADV_CODE || !defined ALLOW_AUTODIFF )
C !LOCAL VARIABLES:
C == Local variables ==
INTEGER i,j,k
#ifndef SOLVE_DIAGONAL_LOWMEMORY
_RL tmpVar, recVar
_RL y5d_m1(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
# ifdef SOLVE_DIAGONAL_KINNER
_RL c5d_prime(Nr)
_RL d5d_prime(Nr)
_RL e5d_prime(Nr)
_RL y5d_prime(Nr)
_RL y5d_update(Nr)
# ifdef ALLOW_AUTODIFF_TAMC
C kkey :: tape key, depends on i,j,k
INTEGER kkey
# endif
# else
_RL c5d_prime(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
_RL d5d_prime(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
_RL e5d_prime(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
_RL y5d_prime(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
# endif
#endif /* SOLVE_DIAGONAL_LOWMEMORY */
CEOP
#ifdef SOLVE_DIAGONAL_LOWMEMORY
IF ( errCode .LT. 0 ) THEN
errCode = 0
DO k=1,Nr
C-- forward sweep (starting from top) part-1: matrix L.U decomposition
IF (k.EQ.1) THEN
DO j=jMin,jMax
DO i=iMin,iMax
IF ( c5d(i,j,k).NE.0. _d 0 ) THEN
c5d(i,j,k) = 1. _d 0 / c5d(i,j,k)
ELSE
c5d(i,j,k) = 0. _d 0
errCode = 1
ENDIF
ENDDO
ENDDO
ELSEIF (k.EQ.2) THEN
DO j=jMin,jMax
DO i=iMin,iMax
C-- [k] <- [k] - b_k/c_k-1 * [k-1]
b5d(i,j,k) = b5d(i,j,k)*c5d(i,j,k-1)
c5d(i,j,k) = c5d(i,j,k) - b5d(i,j,k)*d5d(i,j,k-1)
d5d(i,j,k) = d5d(i,j,k) - b5d(i,j,k)*e5d(i,j,k-1)
IF ( c5d(i,j,k).NE.0. _d 0 ) THEN
c5d(i,j,k) = 1. _d 0 / c5d(i,j,k)
ELSE
c5d(i,j,k) = 0. _d 0
errCode = 1
ENDIF
ENDDO
ENDDO
ELSE
C-- Middle of forward sweep
DO j=jMin,jMax
DO i=iMin,iMax
C-- [k] <- [k] - a_k/c_k-2 * [k-2]
a5d(i,j,k) = a5d(i,j,k)*c5d(i,j,k-2)
b5d(i,j,k) = b5d(i,j,k) - a5d(i,j,k)*d5d(i,j,k-2)
c5d(i,j,k) = c5d(i,j,k) - a5d(i,j,k)*e5d(i,j,k-2)
C-- [k] <- [k] - b_k/c_k-1 * [k-1]
b5d(i,j,k) = b5d(i,j,k)*c5d(i,j,k-1)
c5d(i,j,k) = c5d(i,j,k) - b5d(i,j,k)*d5d(i,j,k-1)
d5d(i,j,k) = d5d(i,j,k) - b5d(i,j,k)*e5d(i,j,k-1)
IF ( c5d(i,j,k).NE.0. _d 0 ) THEN
c5d(i,j,k) = 1. _d 0 / c5d(i,j,k)
ELSE
c5d(i,j,k) = 0. _d 0
errCode = 1
ENDIF
ENDDO
ENDDO
C- end if k= .. ; end of k loop
ENDIF
ENDDO
C-- end if errCode < 0
ENDIF
DO k=2,Nr
C-- forward sweep (starting from top) part-2: process RHS
IF (k.EQ.2) THEN
DO j=jMin,jMax
DO i=iMin,iMax
y5d(i,j,k) = y5d(i,j,k) - b5d(i,j,k)*y5d(i,j,k-1)
ENDDO
ENDDO
ELSE
DO j=jMin,jMax
DO i=iMin,iMax
y5d(i,j,k) = y5d(i,j,k) - b5d(i,j,k)*y5d(i,j,k-1)
& - a5d(i,j,k)*y5d(i,j,k-2)
ENDDO
ENDDO
C- end if k= .. ; end of k loop
ENDIF
ENDDO
C-- Backward sweep (starting from bottom)
DO k=Nr,1,-1
IF (k.EQ.Nr) THEN
DO j=jMin,jMax
DO i=iMin,iMax
y5d(i,j,k) = y5d(i,j,k)*c5d(i,j,k)
ENDDO
ENDDO
ELSEIF (k.EQ.Nr-1) THEN
DO j=jMin,jMax
DO i=iMin,iMax
y5d(i,j,k) = ( y5d(i,j,k)
& - d5d(i,j,k)*y5d(i,j,k+1)
& )*c5d(i,j,k)
ENDDO
ENDDO
ELSE
DO j=jMin,jMax
DO i=iMin,iMax
y5d(i,j,k) = ( y5d(i,j,k)
& - d5d(i,j,k)*y5d(i,j,k+1)
& - e5d(i,j,k)*y5d(i,j,k+2)
& )*c5d(i,j,k)
ENDDO
ENDDO
C- end if k= .. ; end of k loop
ENDIF
ENDDO
#else /* ndef SOLVE_DIAGONAL_LOWMEMORY */
errCode = 0
#ifdef SOLVE_DIAGONAL_KINNER
#ifdef ALLOW_AUTODIFF_TAMC
CADJ INIT loctape_solvepenta = COMMON, Nr*(sNy+2*OLy)*(sNx+2*OLx)
#endif
C-- Temporary array
DO k=1,Nr
DO j=1-OLy,sNy+OLy
DO i=1-OLx,sNx+OLx
y5d_m1(i,j,k) = y5d(i,j,k)
ENDDO
ENDDO
ENDDO
C-- Main loop
DO j=1-OLy,sNy+OLy
DO i=1-OLx,sNx+OLx
DO k=1,Nr
c5d_prime(k) = 0. _d 0
d5d_prime(k) = 0. _d 0
e5d_prime(k) = 0. _d 0
y5d_prime(k) = 0. _d 0
y5d_update(k) = 0. _d 0
ENDDO
DO k=1,Nr
C-- forward sweep (starting from top)
#ifdef ALLOW_AUTODIFF_TAMC
kkey = k + Nr*(i-1+OLx) + Nr*(sNx+2*OLx)*(j-1+OLy)
CADJ STORE d5d_prime(k) = loctape_solvepenta, key = kkey
CADJ STORE e5d_prime(k) = loctape_solvepenta, key = kkey
CADJ STORE y5d_prime(k) = loctape_solvepenta, key = kkey
#endif
IF (k.EQ.1) THEN
c just copy terms
c5d_prime(k) = c5d(i,j,k)
d5d_prime(k) = d5d(i,j,k)
e5d_prime(k) = e5d(i,j,k)
y5d_prime(k) = y5d_m1(i,j,k)
ELSEIF (k.EQ.2) THEN
c subtract one term
c5d_prime(k) = c5d(i,j,k)
& -b5d(i,j,k)*d5d_prime(k-1)
d5d_prime(k) = d5d(i,j,k)
& -b5d(i,j,k)*e5d_prime(k-1)
e5d_prime(k) = e5d(i,j,k)
y5d_prime(k) = y5d_m1(i,j,k)
& -b5d(i,j,k)*y5d_prime(k-1)
ELSE
c subtract two terms
#ifdef ALLOW_AUTODIFF_TAMC
CADJ STORE d5d_prime(k) = loctape_solvepenta, key = kkey
#endif
c5d_prime(k) = c5d(i,j,k)
& -a5d(i,j,k)*e5d_prime(k-2)
& -(b5d(i,j,k)-a5d(i,j,k)*d5d_prime(k-2))*d5d_prime(k-1)
d5d_prime(k) = d5d(i,j,k)
& -(b5d(i,j,k)-a5d(i,j,k)*d5d_prime(k-2))*e5d_prime(k-1)
e5d_prime(k) = e5d(i,j,k)
y5d_prime(k) = y5d_m1(i,j,k)
& -a5d(i,j,k)*y5d_prime(k-2)
& -(b5d(i,j,k)-a5d(i,j,k)*d5d_prime(k-2))*y5d_prime(k-1)
ENDIF
c normalization
tmpVar = c5d_prime(k)
IF ( tmpVar.NE.0. _d 0 ) THEN
recVar = 1. _d 0 / tmpVar
d5d_prime(k) = d5d_prime(k)*recVar
e5d_prime(k) = e5d_prime(k)*recVar
y5d_prime(k) = y5d_prime(k)*recVar
ELSE
d5d_prime(k) = 0. _d 0
e5d_prime(k) = 0. _d 0
y5d_prime(k) = 0. _d 0
errCode = 1
ENDIF
C-- k loop
ENDDO
C-- Backward sweep (starting from bottom)
DO k=Nr,1,-1
IF (k.EQ.Nr) THEN
y5d_update(k) = y5d_prime(k)
ELSEIF (k.EQ.Nr-1) THEN
y5d_update(k) = y5d_prime(k)
& - y5d_update(k+1)*d5d_prime(k)
ELSE
y5d_update(k) = y5d_prime(k)
& - y5d_update(k+1)*d5d_prime(k)
& - y5d_update(k+2)*e5d_prime(k)
ENDIF
ENDDO
C-- Update array
DO k=1,Nr
y5d(i,j,k) = y5d_update(k)
ENDDO
C- end i,j loop
ENDDO
ENDDO
#else /* ndef SOLVE_DIAGONAL_KINNER */
#ifdef ALLOW_AUTODIFF_TAMC
CADJ INIT loctape_solvepenta = COMMON, Nr
#endif
C-- Init. + copy to temp. array
DO k=1,Nr
DO j=1-OLy,sNy+OLy
DO i=1-OLx,sNx+OLx
c5d_prime(i,j,k) = 0. _d 0
d5d_prime(i,j,k) = 0. _d 0
e5d_prime(i,j,k) = 0. _d 0
y5d_prime(i,j,k) = 0. _d 0
y5d_m1(i,j,k) = y5d(i,j,k)
ENDDO
ENDDO
ENDDO
DO k=1,Nr
C-- Forward sweep (starting from top)
#ifdef ALLOW_AUTODIFF_TAMC
CADJ STORE d5d_prime(:,:,k) = loctape_solvepenta, key = k
CADJ STORE e5d_prime(:,:,k) = loctape_solvepenta, key = k
CADJ STORE y5d_prime(:,:,k) = loctape_solvepenta, key = k
#endif
IF ( k .EQ. 1 ) THEN
C- just copy terms
DO j=1-OLy,sNy+OLy
DO i=1-OLx,sNx+OLx
c5d_prime(i,j,k) = c5d(i,j,k)
d5d_prime(i,j,k) = d5d(i,j,k)
e5d_prime(i,j,k) = e5d(i,j,k)
y5d_prime(i,j,k) = y5d_m1(i,j,k)
ENDDO
ENDDO
ELSEIF ( k .EQ. 2 ) THEN
C- subtract one term
DO j=1-OLy,sNy+OLy
DO i=1-OLx,sNx+OLx
tmpVar = b5d(i,j,k)
c5d_prime(i,j,k) = c5d(i,j,k) - tmpVar*d5d_prime(i,j,k-1)
d5d_prime(i,j,k) = d5d(i,j,k) - tmpVar*e5d_prime(i,j,k-1)
e5d_prime(i,j,k) = e5d(i,j,k)
y5d_prime(i,j,k) = y5d_m1(i,j,k) - tmpVar*y5d_prime(i,j,k-1)
ENDDO
ENDDO
ELSE
C- subtract two terms
DO j=1-OLy,sNy+OLy
DO i=1-OLx,sNx+OLx
tmpVar = b5d(i,j,k)-a5d(i,j,k)*d5d_prime(i,j,k-2)
c5d_prime(i,j,k) = c5d(i,j,k) - tmpVar*d5d_prime(i,j,k-1)
& -a5d(i,j,k)*e5d_prime(i,j,k-2)
d5d_prime(i,j,k) = d5d(i,j,k) - tmpVar*e5d_prime(i,j,k-1)
e5d_prime(i,j,k) = e5d(i,j,k)
y5d_prime(i,j,k) = y5d_m1(i,j,k) - tmpVar*y5d_prime(i,j,k-1)
& -a5d(i,j,k)*y5d_prime(i,j,k-2)
ENDDO
ENDDO
ENDIF
C- normalization
DO j=1-OLy,sNy+OLy
DO i=1-OLx,sNx+OLx
tmpVar = c5d_prime(i,j,k)
IF ( tmpVar.NE.0. _d 0 ) THEN
recVar = 1. _d 0 / tmpVar
d5d_prime(i,j,k) = d5d_prime(i,j,k)*recVar
e5d_prime(i,j,k) = e5d_prime(i,j,k)*recVar
y5d_prime(i,j,k) = y5d_prime(i,j,k)*recVar
ELSE
d5d_prime(i,j,k) = 0. _d 0
e5d_prime(i,j,k) = 0. _d 0
y5d_prime(i,j,k) = 0. _d 0
errCode = 1
ENDIF
ENDDO
ENDDO
C- end k-loop
ENDDO
C-- Backward sweep (starting from bottom)
DO k=Nr,1,-1
IF ( k .EQ. Nr ) THEN
DO j=1-OLy,sNy+OLy
DO i=1-OLx,sNx+OLx
y5d(i,j,k) = y5d_prime(i,j,k)
ENDDO
ENDDO
ELSEIF ( k .EQ. Nr-1 ) THEN
DO j=1-OLy,sNy+OLy
DO i=1-OLx,sNx+OLx
y5d(i,j,k) = y5d_prime(i,j,k)
& - y5d(i,j,k+1)*d5d_prime(i,j,k)
ENDDO
ENDDO
ELSE
DO j=1-OLy,sNy+OLy
DO i=1-OLx,sNx+OLx
y5d(i,j,k) = y5d_prime(i,j,k)
& - y5d(i,j,k+1)*d5d_prime(i,j,k)
& - y5d(i,j,k+2)*e5d_prime(i,j,k)
ENDDO
ENDDO
ENDIF
C- k-loop
ENDDO
#endif /* SOLVE_DIAGONAL_KINNER */
#endif /* SOLVE_DIAGONAL_LOWMEMORY */
#endif /* INCLUDE_IMPLVERTADV_CODE */
RETURN
END