1 |
C $Header: /u/gcmpack/MITgcm/model/src/solve_tridiagonal.F,v 1.4 2011/10/13 14:48:44 mlosch Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "CPP_OPTIONS.h" |
5 |
|
6 |
CMLC o Switch to code that has the k-loop inside the |
7 |
CMLC ij-loops, which matters in adjoint mode. |
8 |
CML#ifdef ALLOW_AUTODIFF |
9 |
CML#define ALLOW_SOLVERS_KLOOPINSIDE |
10 |
CML#endif |
11 |
|
12 |
CBOP |
13 |
C !ROUTINE: SOLVE_TRIDIAGONAL |
14 |
C !INTERFACE: |
15 |
SUBROUTINE SOLVE_TRIDIAGONAL( |
16 |
I iMin,iMax, jMin,jMax, |
17 |
I a3d, b3d, c3d, |
18 |
U y3d, |
19 |
O errCode, |
20 |
I bi, bj, myThid ) |
21 |
C !DESCRIPTION: \bv |
22 |
C *==========================================================* |
23 |
C | S/R SOLVE_TRIDIAGONAL |
24 |
C | o Solve a tri-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 a3d :: matrix lower diagnonal |
39 |
C b3d :: matrix main diagnonal |
40 |
C c3d :: matrix upper diagnonal |
41 |
C y3d :: Input = Y vector ; Output = X = solution of A*X=Y |
42 |
C errCode :: > 0 if singular matrix |
43 |
INTEGER iMin,iMax,jMin,jMax |
44 |
_RL a3d(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
45 |
_RL b3d(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
46 |
_RL c3d(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
47 |
_RL y3d(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr,nSx,nSy) |
48 |
INTEGER errCode |
49 |
INTEGER bi, bj, myThid |
50 |
|
51 |
C !LOCAL VARIABLES: |
52 |
C == Local variables == |
53 |
INTEGER i,j,k |
54 |
_RL tmpvar |
55 |
#ifndef ALLOW_SOLVERS_KLOOPINSIDE |
56 |
_RL bet(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
57 |
# if (defined ALLOW_AUTODIFF_TAMC && defined TARGET_NEC_SX) |
58 |
C This field is required for efficient adjoint vector code |
59 |
C The extra storage that is involved is avoided in the case |
60 |
C of forward simulations at the cost of even uglier code. |
61 |
_RL y3dloc(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
62 |
# endif |
63 |
#else |
64 |
_RL c3d_prime(Nr), y3d_prime(Nr),y3d_update(Nr) |
65 |
_RL c3d_m1(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
66 |
_RL y3d_m1(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
67 |
#endif |
68 |
CEOP |
69 |
|
70 |
#ifndef ALLOW_SOLVERS_KLOOPINSIDE |
71 |
|
72 |
errCode = 0 |
73 |
|
74 |
# if (defined ALLOW_AUTODIFF_TAMC && defined TARGET_NEC_SX) |
75 |
DO k=1,Nr |
76 |
DO j=1-Oly,sNy+Oly |
77 |
DO i=1-Olx,sNx+Olx |
78 |
y3dloc(i,j,k) = y3d(i,j,k,bi,bj) |
79 |
y3d (i,j,k,bi,bj) = 0. _d 0 |
80 |
ENDDO |
81 |
ENDDO |
82 |
ENDDO |
83 |
#endif /* ALLOW_AUTODIFF_TAMC and TARGET_NEC_SX */ |
84 |
C-- Beginning of forward sweep (top level) |
85 |
# ifdef TARGET_NEC_SX |
86 |
DO j=1-OLy,sNy+OLy |
87 |
DO i=1-OLx,sNx+OLx |
88 |
# else |
89 |
DO j=jMin,jMax |
90 |
DO i=iMin,iMax |
91 |
# endif |
92 |
IF ( b3d(i,j,1).NE.0. _d 0 ) THEN |
93 |
bet(i,j,1) = 1. _d 0 / b3d(i,j,1) |
94 |
ELSE |
95 |
bet(i,j,1) = 0. _d 0 |
96 |
# ifndef TARGET_NEC_SX |
97 |
C keep your fingers crossed and enjoy the better performance |
98 |
errCode = 1 |
99 |
# endif |
100 |
ENDIF |
101 |
ENDDO |
102 |
ENDDO |
103 |
|
104 |
C-- Middle of forward sweep |
105 |
CADJ loop = sequential |
106 |
DO k=2,Nr |
107 |
# ifdef TARGET_NEC_SX |
108 |
DO j=1-OLy,sNy+OLy |
109 |
DO i=1-OLx,sNx+OLx |
110 |
# else |
111 |
DO j=jMin,jMax |
112 |
DO i=iMin,iMax |
113 |
# endif |
114 |
tmpvar = b3d(i,j,k) - a3d(i,j,k)*c3d(i,j,k-1)*bet(i,j,k-1) |
115 |
IF ( tmpvar .NE. 0. _d 0 ) THEN |
116 |
bet(i,j,k) = 1. _d 0 / tmpvar |
117 |
ELSE |
118 |
bet(i,j,k) = 0. _d 0 |
119 |
# ifndef TARGET_NEC_SX |
120 |
C keep your fingers crossed and enjoy the better performance |
121 |
errCode = 1 |
122 |
# endif |
123 |
ENDIF |
124 |
ENDDO |
125 |
ENDDO |
126 |
ENDDO |
127 |
|
128 |
# ifdef TARGET_NEC_SX |
129 |
DO j=1-OLy,sNy+OLy |
130 |
DO i=1-OLx,sNx+OLx |
131 |
# else |
132 |
DO j=jMin,jMax |
133 |
DO i=iMin,iMax |
134 |
# endif |
135 |
# if (defined ALLOW_AUTODIFF_TAMC && defined TARGET_NEC_SX) |
136 |
y3d(i,j,1,bi,bj) = y3dloc(i,j,1)*bet(i,j,1) |
137 |
# else |
138 |
y3d(i,j,1,bi,bj) = y3d(i,j,1,bi,bj)*bet(i,j,1) |
139 |
# endif |
140 |
ENDDO |
141 |
ENDDO |
142 |
DO k=2,Nr |
143 |
# ifdef TARGET_NEC_SX |
144 |
DO j=1-OLy,sNy+OLy |
145 |
DO i=1-OLx,sNx+OLx |
146 |
# else |
147 |
DO j=jMin,jMax |
148 |
DO i=iMin,iMax |
149 |
# endif |
150 |
# if (defined ALLOW_AUTODIFF_TAMC && defined TARGET_NEC_SX) |
151 |
y3d(i,j,k,bi,bj) = ( y3dloc(i,j,k) |
152 |
# else |
153 |
y3d(i,j,k,bi,bj) = ( y3d(i,j,k,bi,bj) |
154 |
# endif |
155 |
& - a3d(i,j,k)*y3d(i,j,k-1,bi,bj) |
156 |
& )*bet(i,j,k) |
157 |
ENDDO |
158 |
ENDDO |
159 |
ENDDO |
160 |
|
161 |
C-- Backward sweep |
162 |
CADJ loop = sequential |
163 |
DO k=Nr-1,1,-1 |
164 |
# ifdef TARGET_NEC_SX |
165 |
DO j=1-OLy,sNy+OLy |
166 |
DO i=1-OLx,sNx+OLx |
167 |
# else |
168 |
DO j=jMin,jMax |
169 |
DO i=iMin,iMax |
170 |
# endif |
171 |
y3d(i,j,k,bi,bj) = y3d(i,j,k,bi,bj) |
172 |
& - c3d(i,j,k)*bet(i,j,k)*y3d(i,j,k+1,bi,bj) |
173 |
ENDDO |
174 |
ENDDO |
175 |
ENDDO |
176 |
|
177 |
#else /* ALLOW_SOLVERS_KLOOPINSIDE */ |
178 |
|
179 |
errCode = 0 |
180 |
|
181 |
C-- Temporary array |
182 |
DO j=jMin,jMax |
183 |
DO i=iMin,iMax |
184 |
DO k=1,Nr |
185 |
c3d_m1(i,j,k) = c3d(i,j,k) |
186 |
y3d_m1(i,j,k) = y3d(i,j,k,bi,bj) |
187 |
ENDDO |
188 |
ENDDO |
189 |
ENDDO |
190 |
|
191 |
C-- Main loop |
192 |
DO j=jMin,jMax |
193 |
DO i=iMin,iMax |
194 |
|
195 |
DO k=1,Nr |
196 |
c3d_prime(k) = 0. _d 0 |
197 |
y3d_prime(k) = 0. _d 0 |
198 |
y3d_update(k) = 0. _d 0 |
199 |
ENDDO |
200 |
|
201 |
C-- Forward sweep |
202 |
DO k=1,Nr |
203 |
IF ( k.EQ.1 ) THEN |
204 |
IF ( b3d(i,j,1).NE.0. _d 0 ) THEN |
205 |
c3d_prime(1) = c3d_m1(i,j,1) / b3d(i,j,1) |
206 |
y3d_prime(1) = y3d_m1(i,j,1) / b3d(i,j,1) |
207 |
ELSE |
208 |
c3d_prime(1) = 0. _d 0 |
209 |
y3d_prime(1) = 0. _d 0 |
210 |
errCode = 1 |
211 |
ENDIF |
212 |
ELSE |
213 |
tmpvar = b3d(i,j,k) - a3d(i,j,k)*c3d_prime(k-1) |
214 |
IF ( tmpvar .NE. 0. _d 0 ) THEN |
215 |
c3d_prime(k) = c3d_m1(i,j,k) / tmpvar |
216 |
y3d_prime(k) = (y3d_m1(i,j,k) - y3d_prime(k-1)*a3d(i,j,k)) |
217 |
& / tmpvar |
218 |
ELSE |
219 |
c3d_prime(k) = 0. _d 0 |
220 |
y3d_prime(k) = 0. _d 0 |
221 |
errCode = 1 |
222 |
ENDIF |
223 |
ENDIF |
224 |
ENDDO |
225 |
|
226 |
C-- Backward sweep |
227 |
DO k=Nr,1,-1 |
228 |
IF ( k.EQ.Nr ) THEN |
229 |
y3d_update(k)=y3d_prime(k) |
230 |
ELSE |
231 |
y3d_update(k)=y3d_prime(k)-c3d_prime(k)*y3d_update(k+1) |
232 |
ENDIF |
233 |
ENDDO |
234 |
|
235 |
C-- Update array |
236 |
DO k=1,Nr |
237 |
y3d(i,j,k,bi,bj)=y3d_update(k) |
238 |
ENDDO |
239 |
|
240 |
ENDDO |
241 |
ENDDO |
242 |
|
243 |
#endif /* ALLOW_SOLVERS_KLOOPINSIDE */ |
244 |
|
245 |
RETURN |
246 |
END |