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

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

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


Revision 1.16 - (show annotations) (download)
Tue Nov 8 06:18:10 2005 UTC (18 years, 6 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57y_post, checkpoint57y_pre, checkpoint57x_post
Changes since 1.15: +12 -53 lines
- remove 2 unnecessary exchanges (cg3d_b & cg3d_s).
- fix exchange calls for CS-grid (using the new EXCH_S3D_RL)

1 C $Header: /u/gcmpack/MITgcm/model/src/cg3d.F,v 1.15 2005/02/04 19:30:33 jmc Exp $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5
6 CBOP
7 C !ROUTINE: CG3D
8 C !INTERFACE:
9 SUBROUTINE CG3D(
10 I cg3d_b,
11 U cg3d_x,
12 O firstResidual,
13 O lastResidual,
14 U numIters,
15 I myThid )
16 C !DESCRIPTION: \bv
17 C *==========================================================*
18 C | SUBROUTINE CG3D
19 C | o Three-dimensional grid problem conjugate-gradient
20 C | inverter (with preconditioner).
21 C *==========================================================*
22 C | Con. grad is an iterative procedure for solving Ax = b.
23 C | It requires the A be symmetric.
24 C | This implementation assumes A is a seven-diagonal
25 C | matrix of the form that arises in the discrete
26 C | representation of the del^2 operator in a
27 C | three-dimensional space.
28 C | Notes:
29 C | ======
30 C | This implementation can support shared-memory
31 C | multi-threaded execution. In order to do this COMMON
32 C | blocks are used for many of the arrays - even ones that
33 C | are only used for intermedaite results. This design is
34 C | OK if you want to all the threads to collaborate on
35 C | solving the same problem. On the other hand if you want
36 C | the threads to solve several different problems
37 C | concurrently this implementation will not work.
38 C *==========================================================*
39 C \ev
40
41 C !USES:
42 IMPLICIT NONE
43 C === Global data ===
44 #include "SIZE.h"
45 #include "EEPARAMS.h"
46 #include "PARAMS.h"
47 #include "GRID.h"
48 #include "CG3D.h"
49
50 C !INPUT/OUTPUT PARAMETERS:
51 C === Routine arguments ===
52 C myThid - Thread on which I am working.
53 C cg2d_b - The source term or "right hand side"
54 C cg2d_x - The solution
55 C firstResidual - the initial residual before any iterations
56 C lastResidual - the actual residual reached
57 C numIters - Entry: the maximum number of iterations allowed
58 C Exit: the actual number of iterations used
59 _RL cg3d_b(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
60 _RL cg3d_x(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
61 _RL firstResidual
62 _RL lastResidual
63 INTEGER numIters
64 INTEGER myThid
65
66
67 #ifdef ALLOW_NONHYDROSTATIC
68
69 C !LOCAL VARIABLES:
70 C === Local variables ====
71 C actualIts - Number of iterations taken
72 C actualResidual - residual
73 C bi - Block index in X and Y.
74 C bj
75 C eta_qrN - Used in computing search directions
76 C eta_qrNM1 suffix N and NM1 denote current and
77 C cgBeta previous iterations respectively.
78 C alpha
79 C sumRHS - Sum of right-hand-side. Sometimes this is a
80 C useful debuggin/trouble shooting diagnostic.
81 C For neumann problems sumRHS needs to be ~0.
82 C or they converge at a non-zero residual.
83 C err - Measure of residual of Ax - b, usually the norm.
84 C I, J, N - Loop counters ( N counts CG iterations )
85 INTEGER actualIts
86 _RL actualResidual
87 INTEGER bi, bj
88 INTEGER I, J, K, it3d
89 INTEGER KM1, KP1
90 _RL err, errTile
91 _RL eta_qrN, eta_qrNtile
92 _RL eta_qrNM1
93 _RL cgBeta
94 _RL alpha , alphaTile
95 _RL sumRHS, sumRHStile
96 _RL rhsMax
97 _RL rhsNorm
98 _RL topLevTerm
99 CEOP
100
101
102 C-- Initialise inverter
103 eta_qrNM1 = 1. D0
104
105 C-- Normalise RHS
106 rhsMax = 0. _d 0
107 DO bj=myByLo(myThid),myByHi(myThid)
108 DO bi=myBxLo(myThid),myBxHi(myThid)
109 DO K=1,Nr
110 DO J=1,sNy
111 DO I=1,sNx
112 cg3d_b(I,J,K,bi,bj) = cg3d_b(I,J,K,bi,bj)*cg3dNorm
113 rhsMax = MAX(ABS(cg3d_b(I,J,K,bi,bj)),rhsMax)
114 ENDDO
115 ENDDO
116 ENDDO
117 ENDDO
118 ENDDO
119 _GLOBAL_MAX_R8( rhsMax, myThid )
120 rhsNorm = 1. _d 0
121 IF ( rhsMax .NE. 0. ) rhsNorm = 1. _d 0 / rhsMax
122 DO bj=myByLo(myThid),myByHi(myThid)
123 DO bi=myBxLo(myThid),myBxHi(myThid)
124 DO K=1,Nr
125 DO J=1,sNy
126 DO I=1,sNx
127 cg3d_b(I,J,K,bi,bj) = cg3d_b(I,J,K,bi,bj)*rhsNorm
128 cg3d_x(I,J,K,bi,bj) = cg3d_x(I,J,K,bi,bj)*rhsNorm
129 ENDDO
130 ENDDO
131 ENDDO
132 ENDDO
133 ENDDO
134
135 C-- Update overlaps
136 c _EXCH_XYZ_R8( cg3d_b, myThid )
137 _EXCH_XYZ_R8( cg3d_x, myThid )
138
139 C-- Initial residual calculation (with free-Surface term)
140 err = 0. _d 0
141 sumRHS = 0. _d 0
142 DO bj=myByLo(myThid),myByHi(myThid)
143 DO bi=myBxLo(myThid),myBxHi(myThid)
144 errTile = 0. _d 0
145 sumRHStile = 0. _d 0
146 DO K=1,Nr
147 KM1 = K-1
148 IF ( K .EQ. 1 ) KM1 = 1
149 KP1 = K+1
150 IF ( K .EQ. Nr ) KP1 = 1
151 topLevTerm = 0.
152 IF ( K .EQ. 1) topLevTerm = freeSurfFac*cg3dNorm*
153 & (horiVertRatio/gravity)/deltaTMom/deltaTMom
154 DO J=1,sNy
155 DO I=1,sNx
156 cg3d_r(I,J,K,bi,bj) = cg3d_b(I,J,K,bi,bj) -( 0.
157 & +aW3d(I ,J ,K ,bi,bj)*cg3d_x(I-1,J ,K ,bi,bj)
158 & +aW3d(I+1,J ,K ,bi,bj)*cg3d_x(I+1,J ,K ,bi,bj)
159 & +aS3d(I ,J ,K ,bi,bj)*cg3d_x(I ,J-1,K ,bi,bj)
160 & +aS3d(I ,J+1,K ,bi,bj)*cg3d_x(I ,J+1,K ,bi,bj)
161 & +aV3d(I ,J ,K ,bi,bj)*cg3d_x(I ,J ,KM1,bi,bj)
162 & +aV3d(I ,J ,KP1,bi,bj)*cg3d_x(I ,J ,KP1,bi,bj)
163 & -aW3d(I ,J ,K ,bi,bj)*cg3d_x(I ,J ,K ,bi,bj)
164 & -aW3d(I+1,J ,K ,bi,bj)*cg3d_x(I ,J ,K ,bi,bj)
165 & -aS3d(I ,J ,K ,bi,bj)*cg3d_x(I ,J ,K ,bi,bj)
166 & -aS3d(I ,J+1,K ,bi,bj)*cg3d_x(I ,J ,K ,bi,bj)
167 & -aV3d(I ,J ,K ,bi,bj)*cg3d_x(I ,J ,K ,bi,bj)
168 & -aV3d(I ,J ,KP1,bi,bj)*cg3d_x(I ,J ,K ,bi,bj)
169 & -topLevTerm*_rA(I,J,bi,bj)*cg3d_x(I,J,K,bi,bj)
170 & )
171 errTile = errTile
172 & +cg3d_r(I,J,K,bi,bj)*cg3d_r(I,J,K,bi,bj)
173 sumRHStile = sumRHStile
174 & +cg3d_b(I,J,K,bi,bj)
175 ENDDO
176 ENDDO
177 DO J=1-1,sNy+1
178 DO I=1-1,sNx+1
179 cg3d_s(I,J,K,bi,bj) = 0.
180 ENDDO
181 ENDDO
182 ENDDO
183 err = err + errTile
184 sumRHS = sumRHS + sumRHStile
185 ENDDO
186 ENDDO
187 C _EXCH_XYZ_R8( cg3d_r, myThid )
188 CALL EXCH_S3D_RL( cg3d_r, myThid )
189 C _EXCH_XYZ_R8( cg3d_s, myThid )
190 c CALL EXCH_S3D_RL( cg3d_s, myThid )
191 _GLOBAL_SUM_R8( sumRHS, myThid )
192 _GLOBAL_SUM_R8( err , myThid )
193
194 IF ( debugLevel .GE. debLevZero ) THEN
195 _BEGIN_MASTER( myThid )
196 write(standardmessageunit,'(A,1P2E22.14)')
197 & ' cg3d: Sum(rhs),rhsMax = ',sumRHS,rhsMax
198 _END_MASTER( myThid )
199 ENDIF
200
201 actualIts = 0
202 actualResidual = SQRT(err)
203 C _BARRIER
204 c _BEGIN_MASTER( myThid )
205 c WRITE(*,'(A,I6,1PE30.14)') ' CG3D iters, err = ',
206 c & actualIts, actualResidual
207 c _END_MASTER( myThid )
208 firstResidual=actualResidual
209
210 C >>>>>>>>>>>>>>> BEGIN SOLVER <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
211 DO 10 it3d=1, cg3dMaxIters
212
213 CcnhDebugStarts
214 c IF ( mod(it3d-1,10).EQ.0)
215 c & WRITE(*,*) ' CG3D: Iteration ',it3d-1,
216 c & ' residual = ',actualResidual
217 CcnhDebugEnds
218 IF ( actualResidual .LT. cg3dTargetResidual ) GOTO 11
219 C-- Solve preconditioning equation and update
220 C-- conjugate direction vector "s".
221 C Note. On the next to loops over all tiles the inner loop ranges
222 C in sNx and sNy are expanded by 1 to avoid a communication
223 C step. However this entails a bit of gynamastics because we only
224 C want eta_qrN for the interior points.
225 eta_qrN = 0. _d 0
226 DO bj=myByLo(myThid),myByHi(myThid)
227 DO bi=myBxLo(myThid),myBxHi(myThid)
228 eta_qrNtile = 0. _d 0
229 DO K=1,1
230 DO J=1-1,sNy+1
231 DO I=1-1,sNx+1
232 cg3d_q(I,J,K,bi,bj) =
233 & zMC(I ,J ,K,bi,bj)*cg3d_r(I ,J ,K,bi,bj)
234 ENDDO
235 ENDDO
236 ENDDO
237 DO K=2,Nr
238 DO J=1-1,sNy+1
239 DO I=1-1,sNx+1
240 cg3d_q(I,J,K,bi,bj) =
241 & zMC(I,J,K,bi,bj)*(cg3d_r(I,J,K ,bi,bj)
242 & -zML(I,J,K,bi,bj)*cg3d_q(I,J,K-1,bi,bj))
243 ENDDO
244 ENDDO
245 ENDDO
246 DO K=Nr,Nr
247 caja IF (Nr .GT. 1) THEN
248 caja DO J=1-1,sNy+1
249 caja DO I=1-1,sNx+1
250 caja cg3d_q(I,J,K,bi,bj) =
251 caja & zMC(i,j,k,bi,bj)*(cg3d_r(i,j,k ,bi,bj)
252 caja & -zML(i,j,k,bi,bj)*cg3d_q(i,j,k-1,bi,bj))
253 caja ENDDO
254 caja ENDDO
255 caja ENDIF
256 DO J=1,sNy
257 DO I=1,sNx
258 eta_qrNtile = eta_qrNtile
259 & +cg3d_q(I,J,K,bi,bj)*cg3d_r(I,J,K,bi,bj)
260 ENDDO
261 ENDDO
262 ENDDO
263 DO K=Nr-1,1,-1
264 DO J=1-1,sNy+1
265 DO I=1-1,sNx+1
266 cg3d_q(I,J,K,bi,bj) =
267 & cg3d_q(I,J,K,bi,bj)
268 & -zMU(I,J,K,bi,bj)*cg3d_q(I,J,K+1,bi,bj)
269 ENDDO
270 ENDDO
271 DO J=1,sNy
272 DO I=1,sNx
273 eta_qrNtile = eta_qrNtile
274 & +cg3d_q(I,J,K,bi,bj)*cg3d_r(I,J,K,bi,bj)
275 ENDDO
276 ENDDO
277 ENDDO
278 eta_qrN = eta_qrN + eta_qrNtile
279 ENDDO
280 ENDDO
281 caja
282 caja eta_qrN=0.
283 caja DO bj=myByLo(myThid),myByHi(myThid)
284 caja DO bi=myBxLo(myThid),myBxHi(myThid)
285 caja DO K=1,Nr
286 caja DO J=1,sNy
287 caja DO I=1,sNx
288 caja eta_qrN = eta_qrN
289 caja & +cg3d_q(I,J,K,bi,bj)*cg3d_r(I,J,K,bi,bj)
290 caja ENDDO
291 caja ENDDO
292 caja ENDDO
293 caja ENDDO
294 caja ENDDO
295 caja
296
297 _GLOBAL_SUM_R8(eta_qrN, myThid)
298 CcnhDebugStarts
299 C WRITE(*,*) ' CG3D: Iteration ',it3d-1,' eta_qrN = ',eta_qrN
300 CcnhDebugEnds
301 cgBeta = eta_qrN/eta_qrNM1
302 CcnhDebugStarts
303 C WRITE(*,*) ' CG3D: Iteration ',it3d-1,' beta = ',cgBeta
304 CcnhDebugEnds
305 eta_qrNM1 = eta_qrN
306
307 DO bj=myByLo(myThid),myByHi(myThid)
308 DO bi=myBxLo(myThid),myBxHi(myThid)
309 DO K=1,Nr
310 DO J=1-1,sNy+1
311 DO I=1-1,sNx+1
312 cg3d_s(I,J,K,bi,bj) = cg3d_q(I,J,K,bi,bj)
313 & + cgBeta*cg3d_s(I,J,K,bi,bj)
314 ENDDO
315 ENDDO
316 ENDDO
317 ENDDO
318 ENDDO
319
320 C== Evaluate laplace operator on conjugate gradient vector
321 C== q = A.s
322 alpha = 0. _d 0
323 topLevTerm = freeSurfFac*cg3dNorm*
324 & (horiVertRatio/gravity)/deltaTMom/deltaTMom
325 DO bj=myByLo(myThid),myByHi(myThid)
326 DO bi=myBxLo(myThid),myBxHi(myThid)
327 alphaTile = 0. _d 0
328 IF ( Nr .GT. 1 ) THEN
329 DO K=1,1
330 DO J=1,sNy
331 DO I=1,sNx
332 cg3d_q(I,J,K,bi,bj) =
333 & aW3d(I ,J ,K ,bi,bj)*cg3d_s(I-1,J ,K ,bi,bj)
334 & +aW3d(I+1,J ,K ,bi,bj)*cg3d_s(I+1,J ,K ,bi,bj)
335 & +aS3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J-1,K ,bi,bj)
336 & +aS3d(I ,J+1,K ,bi,bj)*cg3d_s(I ,J+1,K ,bi,bj)
337 & +aV3d(I ,J ,K+1,bi,bj)*cg3d_s(I ,J ,K+1,bi,bj)
338 & -aW3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
339 & -aW3d(I+1,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
340 & -aS3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
341 & -aS3d(I ,J+1,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
342 & -aV3d(I ,J ,K+1,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
343 & -topLevTerm*_rA(I,J,bi,bj)*cg3d_s(I,J,K,bi,bj)
344 alphaTile = alphaTile
345 & +cg3d_s(I,J,K,bi,bj)*cg3d_q(I,J,K,bi,bj)
346 ENDDO
347 ENDDO
348 ENDDO
349 ELSE
350 DO K=1,1
351 DO J=1,sNy
352 DO I=1,sNx
353 cg3d_q(I,J,K,bi,bj) =
354 & aW3d(I ,J ,K ,bi,bj)*cg3d_s(I-1,J ,K ,bi,bj)
355 & +aW3d(I+1,J ,K ,bi,bj)*cg3d_s(I+1,J ,K ,bi,bj)
356 & +aS3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J-1,K ,bi,bj)
357 & +aS3d(I ,J+1,K ,bi,bj)*cg3d_s(I ,J+1,K ,bi,bj)
358 & -aW3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
359 & -aW3d(I+1,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
360 & -aS3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
361 & -aS3d(I ,J+1,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
362 & -topLevTerm*_rA(I,J,bi,bj)*cg3d_s(I,J,K,bi,bj)
363 alphaTile = alphaTile
364 & +cg3d_s(I,J,K,bi,bj)*cg3d_q(I,J,K,bi,bj)
365 ENDDO
366 ENDDO
367 ENDDO
368 ENDIF
369 DO K=2,Nr-1
370 DO J=1,sNy
371 DO I=1,sNx
372 cg3d_q(I,J,K,bi,bj) =
373 & aW3d(I ,J ,K ,bi,bj)*cg3d_s(I-1,J ,K ,bi,bj)
374 & +aW3d(I+1,J ,K ,bi,bj)*cg3d_s(I+1,J ,K ,bi,bj)
375 & +aS3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J-1,K ,bi,bj)
376 & +aS3d(I ,J+1,K ,bi,bj)*cg3d_s(I ,J+1,K ,bi,bj)
377 & +aV3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K-1,bi,bj)
378 & +aV3d(I ,J ,K+1,bi,bj)*cg3d_s(I ,J ,K+1,bi,bj)
379 & -aW3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
380 & -aW3d(I+1,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
381 & -aS3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
382 & -aS3d(I ,J+1,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
383 & -aV3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
384 & -aV3d(I ,J ,K+1,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
385 alphaTile = alphaTile
386 & +cg3d_s(I,J,K,bi,bj)*cg3d_q(I,J,K,bi,bj)
387 ENDDO
388 ENDDO
389 ENDDO
390 IF ( Nr .GT. 1 ) THEN
391 DO K=Nr,Nr
392 DO J=1,sNy
393 DO I=1,sNx
394 cg3d_q(I,J,K,bi,bj) =
395 & aW3d(I ,J ,K ,bi,bj)*cg3d_s(I-1,J ,K ,bi,bj)
396 & +aW3d(I+1,J ,K ,bi,bj)*cg3d_s(I+1,J ,K ,bi,bj)
397 & +aS3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J-1,K ,bi,bj)
398 & +aS3d(I ,J+1,K ,bi,bj)*cg3d_s(I ,J+1,K ,bi,bj)
399 & +aV3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K-1,bi,bj)
400 & -aW3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
401 & -aW3d(I+1,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
402 & -aS3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
403 & -aS3d(I ,J+1,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
404 & -aV3d(I ,J ,K ,bi,bj)*cg3d_s(I ,J ,K ,bi,bj)
405 alphaTile = alphaTile
406 & +cg3d_s(I,J,K,bi,bj)*cg3d_q(I,J,K,bi,bj)
407 ENDDO
408 ENDDO
409 ENDDO
410 ENDIF
411 alpha = alpha + alphaTile
412 ENDDO
413 ENDDO
414 _GLOBAL_SUM_R8(alpha,myThid)
415 CcnhDebugStarts
416 C WRITE(*,*) ' CG3D: Iteration ',it3d-1,' SUM(s*q)= ',alpha
417 CcnhDebugEnds
418 alpha = eta_qrN/alpha
419 CcnhDebugStarts
420 C WRITE(*,*) ' CG3D: Iteration ',it3d-1,' alpha= ',alpha
421 CcnhDebugEnds
422
423 C== Update solution and residual vectors
424 C Now compute "interior" points.
425 err = 0. _d 0
426 DO bj=myByLo(myThid),myByHi(myThid)
427 DO bi=myBxLo(myThid),myBxHi(myThid)
428 errTile = 0. _d 0
429 DO K=1,Nr
430 DO J=1,sNy
431 DO I=1,sNx
432 cg3d_x(I,J,K,bi,bj)=cg3d_x(I,J,K,bi,bj)
433 & +alpha*cg3d_s(I,J,K,bi,bj)
434 cg3d_r(I,J,K,bi,bj)=cg3d_r(I,J,K,bi,bj)
435 & -alpha*cg3d_q(I,J,K,bi,bj)
436 errTile = errTile
437 & +cg3d_r(I,J,K,bi,bj)*cg3d_r(I,J,K,bi,bj)
438 ENDDO
439 ENDDO
440 ENDDO
441 err = err + errTile
442 ENDDO
443 ENDDO
444
445 _GLOBAL_SUM_R8( err , myThid )
446 err = SQRT(err)
447 actualIts = it3d
448 actualResidual = err
449 IF ( actualResidual .LT. cg3dTargetResidual ) GOTO 11
450 C _EXCH_XYZ_R8(cg3d_r, myThid )
451 CALL EXCH_S3D_RL( cg3d_r, myThid )
452
453 10 CONTINUE
454 11 CONTINUE
455
456 C-- Un-normalise the answer
457 DO bj=myByLo(myThid),myByHi(myThid)
458 DO bi=myBxLo(myThid),myBxHi(myThid)
459 DO K=1,Nr
460 DO J=1,sNy
461 DO I=1,sNx
462 cg3d_x(I,J,K,bi,bj) = cg3d_x(I,J,K,bi,bj)/rhsNorm
463 ENDDO
464 ENDDO
465 ENDDO
466 ENDDO
467 ENDDO
468
469 Cadj _EXCH_XYZ_R8(cg3d_x, myThid )
470 c _BEGIN_MASTER( myThid )
471 c WRITE(*,'(A,I6,1PE30.14)') ' CG3D iters, err = ',
472 c & actualIts, actualResidual
473 c _END_MASTER( myThid )
474 lastResidual=actualResidual
475 numIters=actualIts
476
477 #endif /* ALLOW_NONHYDROSTATIC */
478
479 RETURN
480 END

  ViewVC Help
Powered by ViewVC 1.1.22