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