1 |
C $Header: /u/gcmpack/MITgcm/model/src/solve_pentadiagonal.F,v 1.3 2010/03/16 00:08:27 jmc Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "CPP_OPTIONS.h" |
5 |
|
6 |
C o Switch to code that has the k-loop inside the |
7 |
C ij-loops, which matters in adjoint mode. |
8 |
#ifdef ALLOW_AUTODIFF |
9 |
#define ALLOW_SOLVERS_KLOOPINSIDE |
10 |
#endif |
11 |
|
12 |
CBOP |
13 |
C !ROUTINE: SOLVE_PENTADIAGONAL |
14 |
C !INTERFACE: |
15 |
SUBROUTINE SOLVE_PENTADIAGONAL( |
16 |
I iMin,iMax, jMin,jMax, |
17 |
U a5d, b5d, c5d, d5d, e5d, |
18 |
U y5d, |
19 |
O errCode, |
20 |
I bi, bj, myThid ) |
21 |
C !DESCRIPTION: \bv |
22 |
C *==========================================================* |
23 |
C | S/R SOLVE_PENTADIAGONAL |
24 |
C | o Solve a penta-diagonal system A*X=Y (dimension Nr) |
25 |
C *==========================================================* |
26 |
C | o Used to solve implicitly vertical advection & diffusion |
27 |
C *==========================================================* |
28 |
C \ev |
29 |
|
30 |
C !USES: |
31 |
IMPLICIT NONE |
32 |
C == Global data == |
33 |
#include "SIZE.h" |
34 |
#include "EEPARAMS.h" |
35 |
|
36 |
C !INPUT/OUTPUT PARAMETERS: |
37 |
C == Routine Arguments == |
38 |
C INPUT: |
39 |
C iMin,iMax,jMin,jMax :: computational domain |
40 |
C a5d :: 2nd lower diagonal of the pentadiagonal matrix |
41 |
C b5d :: 1rst lower diagonal of the pentadiagonal matrix |
42 |
C c5d :: main diagonal of the pentadiagonal matrix |
43 |
C d5d :: 1rst upper diagonal of the pentadiagonal matrix |
44 |
C e5d :: 2nd upper diagonal of the pentadiagonal matrix |
45 |
C y5d :: Y vector (R.H.S.); |
46 |
C bi,bj :: tile indices |
47 |
C myThid :: thread number |
48 |
C OUTPUT: |
49 |
C y5d :: X = solution of A*X=Y |
50 |
C a5d,b5d,c5d,d5d,e5d :: modified to enable to find Xp solution of |
51 |
C A*Xp=Yp without solving the full system again |
52 |
C errCode :: > 0 if singular matrix |
53 |
INTEGER iMin,iMax,jMin,jMax |
54 |
_RL a5d(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
55 |
_RL b5d(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
56 |
_RL c5d(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
57 |
_RL d5d(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
58 |
_RL e5d(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
59 |
_RL y5d(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr,nSx,nSy) |
60 |
INTEGER errCode |
61 |
INTEGER bi, bj, myThid |
62 |
|
63 |
C !LOCAL VARIABLES: |
64 |
C == Local variables == |
65 |
INTEGER i,j,k |
66 |
#ifdef ALLOW_SOLVERS_KLOOPINSIDE |
67 |
_RL y5d_m1(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
68 |
_RL a5d_prime(Nr), b5d_prime(Nr) |
69 |
_RL c5d_prime(Nr), d5d_prime(Nr), e5d_prime(Nr) |
70 |
_RL y5d_prime(Nr), y5d_update(Nr), tmpval |
71 |
#endif |
72 |
CEOP |
73 |
|
74 |
#ifndef ALLOW_SOLVERS_KLOOPINSIDE |
75 |
|
76 |
errCode = 0 |
77 |
|
78 |
DO k=1,Nr |
79 |
C-- forward sweep (starting from top) |
80 |
IF (k.EQ.1) THEN |
81 |
DO j=jMin,jMax |
82 |
DO i=iMin,iMax |
83 |
IF ( c5d(i,j,k).NE.0. _d 0 ) THEN |
84 |
c5d(i,j,k) = 1. _d 0 / c5d(i,j,k) |
85 |
ELSE |
86 |
c5d(i,j,k) = 0. _d 0 |
87 |
errCode = 1 |
88 |
ENDIF |
89 |
ENDDO |
90 |
ENDDO |
91 |
|
92 |
ELSEIF (k.EQ.2) THEN |
93 |
DO j=jMin,jMax |
94 |
DO i=iMin,iMax |
95 |
C-- [k] <- [k] - b_k/c_k-1 * [k-1] |
96 |
b5d(i,j,k) = b5d(i,j,k)*c5d(i,j,k-1) |
97 |
c5d(i,j,k) = c5d(i,j,k) - b5d(i,j,k)*d5d(i,j,k-1) |
98 |
d5d(i,j,k) = d5d(i,j,k) - b5d(i,j,k)*e5d(i,j,k-1) |
99 |
y5d(i,j,k,bi,bj) = y5d(i,j,k,bi,bj) |
100 |
& - b5d(i,j,k)*y5d(i,j,k-1,bi,bj) |
101 |
IF ( c5d(i,j,k).NE.0. _d 0 ) THEN |
102 |
c5d(i,j,k) = 1. _d 0 / c5d(i,j,k) |
103 |
ELSE |
104 |
c5d(i,j,k) = 0. _d 0 |
105 |
errCode = 1 |
106 |
ENDIF |
107 |
ENDDO |
108 |
ENDDO |
109 |
|
110 |
ELSE |
111 |
C-- Middle of forward sweep |
112 |
DO j=jMin,jMax |
113 |
DO i=iMin,iMax |
114 |
C-- [k] <- [k] - a_k/c_k-2 * [k-2] |
115 |
a5d(i,j,k) = a5d(i,j,k)*c5d(i,j,k-2) |
116 |
b5d(i,j,k) = b5d(i,j,k) - a5d(i,j,k)*d5d(i,j,k-2) |
117 |
c5d(i,j,k) = c5d(i,j,k) - a5d(i,j,k)*e5d(i,j,k-2) |
118 |
C-- [k] <- [k] - b_k/c_k-1 * [k-1] |
119 |
b5d(i,j,k) = b5d(i,j,k)*c5d(i,j,k-1) |
120 |
c5d(i,j,k) = c5d(i,j,k) - b5d(i,j,k)*d5d(i,j,k-1) |
121 |
d5d(i,j,k) = d5d(i,j,k) - b5d(i,j,k)*e5d(i,j,k-1) |
122 |
y5d(i,j,k,bi,bj) = y5d(i,j,k,bi,bj) |
123 |
& - b5d(i,j,k)*y5d(i,j,k-1,bi,bj) |
124 |
& - a5d(i,j,k)*y5d(i,j,k-2,bi,bj) |
125 |
IF ( c5d(i,j,k).NE.0. _d 0 ) THEN |
126 |
c5d(i,j,k) = 1. _d 0 / c5d(i,j,k) |
127 |
ELSE |
128 |
c5d(i,j,k) = 0. _d 0 |
129 |
errCode = 1 |
130 |
ENDIF |
131 |
ENDDO |
132 |
ENDDO |
133 |
C- end if k= .. ; end of k loop |
134 |
ENDIF |
135 |
ENDDO |
136 |
|
137 |
C-- Backward sweep (starting from bottom) |
138 |
DO k=Nr,1,-1 |
139 |
IF (k.EQ.Nr) THEN |
140 |
DO j=jMin,jMax |
141 |
DO i=iMin,iMax |
142 |
y5d(i,j,k,bi,bj) = y5d(i,j,k,bi,bj)*c5d(i,j,k) |
143 |
ENDDO |
144 |
ENDDO |
145 |
ELSEIF (k.EQ.Nr-1) THEN |
146 |
DO j=jMin,jMax |
147 |
DO i=iMin,iMax |
148 |
y5d(i,j,k,bi,bj) = ( y5d(i,j,k,bi,bj) |
149 |
& - d5d(i,j,k)*y5d(i,j,k+1,bi,bj) |
150 |
& )*c5d(i,j,k) |
151 |
ENDDO |
152 |
ENDDO |
153 |
ELSE |
154 |
DO j=jMin,jMax |
155 |
DO i=iMin,iMax |
156 |
y5d(i,j,k,bi,bj) = ( y5d(i,j,k,bi,bj) |
157 |
& - d5d(i,j,k)*y5d(i,j,k+1,bi,bj) |
158 |
& - e5d(i,j,k)*y5d(i,j,k+2,bi,bj) |
159 |
& )*c5d(i,j,k) |
160 |
ENDDO |
161 |
ENDDO |
162 |
C- end if k= .. ; end of k loop |
163 |
ENDIF |
164 |
ENDDO |
165 |
|
166 |
#else /* ALLOW_SOLVERS_KLOOPINSIDE */ |
167 |
|
168 |
errCode = 0 |
169 |
|
170 |
C-- Temporary array |
171 |
DO j=jMin,jMax |
172 |
DO i=iMin,iMax |
173 |
DO k=1,Nr |
174 |
y5d_m1(i,j,k) = y5d(i,j,k,bi,bj) |
175 |
ENDDO |
176 |
ENDDO |
177 |
ENDDO |
178 |
|
179 |
C-- Main loop |
180 |
DO j=jMin,jMax |
181 |
DO i=iMin,iMax |
182 |
|
183 |
DO k=1,Nr |
184 |
a5d_prime(k) = 0. _d 0 |
185 |
b5d_prime(k) = 0. _d 0 |
186 |
c5d_prime(k) = 0. _d 0 |
187 |
d5d_prime(k) = 0. _d 0 |
188 |
e5d_prime(k) = 0. _d 0 |
189 |
y5d_prime(k) = 0. _d 0 |
190 |
y5d_update(k) = 0. _d 0 |
191 |
ENDDO |
192 |
|
193 |
DO k=1,Nr |
194 |
C-- forward sweep (starting from top) |
195 |
|
196 |
IF (k.EQ.1) THEN |
197 |
c just copy terms |
198 |
a5d_prime(k) = 0. _d 0 |
199 |
b5d_prime(k) = 0. _d 0 |
200 |
c5d_prime(k) = c5d(i,j,k) |
201 |
d5d_prime(k) = d5d(i,j,k) |
202 |
e5d_prime(k) = e5d(i,j,k) |
203 |
y5d_prime(k) = y5d_m1(i,j,k) |
204 |
ELSEIF (k.EQ.2) THEN |
205 |
c subtract one term |
206 |
a5d_prime(k) = 0. _d 0 |
207 |
b5d_prime(k) = 0. _d 0 |
208 |
c5d_prime(k) = c5d(i,j,k) |
209 |
& -b5d(i,j,k)*d5d_prime(k-1) |
210 |
d5d_prime(k) = d5d(i,j,k) |
211 |
& -b5d(i,j,k)*e5d_prime(k-1) |
212 |
e5d_prime(k) = e5d(i,j,k) |
213 |
y5d_prime(k) = y5d_m1(i,j,k) |
214 |
& -b5d(i,j,k)*y5d_prime(k-1) |
215 |
ELSE |
216 |
c subtract two terms |
217 |
a5d_prime(k) = 0. _d 0 |
218 |
b5d_prime(k) = 0. _d 0 |
219 |
c5d_prime(k) = c5d(i,j,k) |
220 |
& -a5d(i,j,k)*e5d_prime(k-2) |
221 |
& -(b5d(i,j,k)-a5d(i,j,k)*d5d_prime(k-2))*d5d_prime(k-1) |
222 |
d5d_prime(k) = d5d(i,j,k) |
223 |
& -(b5d(i,j,k)-a5d(i,j,k)*d5d_prime(k-2))*e5d_prime(k-1) |
224 |
e5d_prime(k) = e5d(i,j,k) |
225 |
y5d_prime(k) = y5d_m1(i,j,k) |
226 |
& -a5d(i,j,k)*y5d_prime(k-2) |
227 |
& -(b5d(i,j,k)-a5d(i,j,k)*d5d_prime(k-2))*y5d_prime(k-1) |
228 |
ENDIF |
229 |
|
230 |
c normalization |
231 |
tmpval=c5d_prime(k) |
232 |
IF ( tmpval.NE.0. _d 0 ) THEN |
233 |
a5d_prime(k) = a5d_prime(k) / tmpval |
234 |
b5d_prime(k) = b5d_prime(k) / tmpval |
235 |
c5d_prime(k) = 1. _d 0 |
236 |
d5d_prime(k) = d5d_prime(k) / tmpval |
237 |
e5d_prime(k) = e5d_prime(k) / tmpval |
238 |
y5d_prime(k) = y5d_prime(k) / tmpval |
239 |
ELSE |
240 |
a5d_prime(k) = 0. _d 0 |
241 |
b5d_prime(k) = 0. _d 0 |
242 |
c5d_prime(k) = 0. _d 0 |
243 |
d5d_prime(k) = 0. _d 0 |
244 |
e5d_prime(k) = 0. _d 0 |
245 |
y5d_prime(k) = 0. _d 0 |
246 |
errCode = 1 |
247 |
ENDIF |
248 |
|
249 |
ENDDO |
250 |
|
251 |
C-- Backward sweep (starting from bottom) |
252 |
DO k=Nr,1,-1 |
253 |
IF (k.EQ.Nr) THEN |
254 |
y5d_update(k) = y5d_prime(k) |
255 |
ELSEIF (k.EQ.Nr-1) THEN |
256 |
y5d_update(k) = y5d_prime(k) |
257 |
& - y5d_update(k+1)*d5d_prime(k) |
258 |
ELSE |
259 |
y5d_update(k) = y5d_prime(k) |
260 |
& - y5d_update(k+1)*d5d_prime(k) |
261 |
& - y5d_update(k+2)*e5d_prime(k) |
262 |
ENDIF |
263 |
ENDDO |
264 |
|
265 |
C-- Update array |
266 |
DO k=1,Nr |
267 |
y5d(i,j,k,bi,bj)=y5d_update(k) |
268 |
ENDDO |
269 |
|
270 |
ENDDO |
271 |
ENDDO |
272 |
|
273 |
#endif /* ALLOW_SOLVERS_KLOOPINSIDE */ |
274 |
|
275 |
RETURN |
276 |
END |