/[MITgcm]/MITgcm/pkg/admtlm/admtlm_dsvd.F
ViewVC logotype

Contents of /MITgcm/pkg/admtlm/admtlm_dsvd.F

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


Revision 1.3 - (show annotations) (download)
Wed Nov 17 03:06:12 2004 UTC (19 years, 7 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint57d_post, checkpoint57g_post, checkpoint57b_post, checkpoint57c_pre, checkpoint56b_post, checkpoint57e_post, checkpoint57g_pre, checkpoint56c_post, checkpoint57f_pre, checkpoint57a_post, checkpoint57a_pre, checkpoint57, checkpoint56, eckpoint57e_pre, checkpoint57f_post, checkpoint57c_post, checkpoint56a_post
Changes since 1.2: +72 -50 lines
o finishing interface btw. dsvd and MITgcm
o SVD seems complete w.r.t. L2 norm of SST, SSS

1 C $Header: $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5
6 subroutine admtlm_dsvd ( mythid )
7 c
8 c This example program is intended to illustrate the
9 c the use of ARPACK to compute the Singular Value Decomposition.
10 c
11 c This code shows how to use ARPACK to find a few of the
12 c largest singular values(sigma) and corresponding right singular
13 c vectors (v) for the the matrix A by solving the symmetric problem:
14 c
15 c (A'*A)*v = sigma*v
16 c
17 c where A is an m by n real matrix.
18 c
19 c This code may be easily modified to estimate the 2-norm
20 c condition number largest(sigma)/smallest(sigma) by setting
21 c which = 'BE' below. This will ask for a few of the smallest
22 c and a few of the largest singular values simultaneously.
23 c The condition number could then be estimated by taking
24 c the ratio of the largest and smallest singular values.
25 c
26 c This formulation is appropriate when m .ge. n.
27 c Reverse the roles of A and A' in the case that m .le. n.
28 c
29 c The main points illustrated here are
30 c
31 c 1) How to declare sufficient memory to find NEV
32 c largest singular values of A .
33 c
34 c 2) Illustration of the reverse communication interface
35 c needed to utilize the top level ARPACK routine DSAUPD
36 c that computes the quantities needed to construct
37 c the desired singular values and vectors(if requested).
38 c
39 c 3) How to extract the desired singular values and vectors
40 c using the ARPACK routine DSEUPD.
41 c
42 c 4) How to construct the left singular vectors U from the
43 c right singular vectors V to obtain the decomposition
44 c
45 c A*V = U*S
46 c
47 c where S = diag(sigma_1, sigma_2, ..., sigma_k).
48 c
49 c The only thing that must be supplied in order to use this
50 c routine on your problem is to change the array dimensions
51 c appropriately, to specify WHICH singular values you want to
52 c compute and to supply a the matrix-vector products
53 c
54 c w <- Ax
55 c y <- A'w
56 c
57 c in place of the calls to AV( ) and ATV( ) respectively below.
58 c
59 c Further documentation is available in the header of DSAUPD
60 c which may be found in the SRC directory.
61 c
62 c This codes implements
63 c
64 c\Example-1
65 c ... Suppose we want to solve A'A*v = sigma*v in regular mode,
66 c where A is derived from the simplest finite difference
67 c discretization of the 2-dimensional kernel K(s,t)dt where
68 c
69 c K(s,t) = s(t-1) if 0 .le. s .le. t .le. 1,
70 c t(s-1) if 0 .le. t .lt. s .le. 1.
71 c
72 c See subroutines AV and ATV for details.
73 c ... OP = A'*A and B = I.
74 c ... Assume "call av (n,x,y)" computes y = A*x
75 c ... Assume "call atv (n,y,w)" computes w = A'*y
76 c ... Assume exact shifts are used
77 c ...
78 c
79 c\BeginLib
80 c
81 c\Routines called:
82 c dsaupd ARPACK reverse communication interface routine.
83 c dseupd ARPACK routine that returns Ritz values and (optionally)
84 c Ritz vectors.
85 c dnrm2 Level 1 BLAS that computes the norm of a vector.
86 c daxpy Level 1 BLAS that computes y <- alpha*x+y.
87 c dscal Level 1 BLAS thst computes x <- x*alpha.
88 c dcopy Level 1 BLAS thst computes y <- x.
89 c
90 c\Author
91 c Richard Lehoucq
92 c Danny Sorensen
93 c Chao Yang
94 c Dept. of Computational &
95 c Applied Mathematics
96 c Rice University
97 c Houston, Texas
98 c
99 c\SCCS Information: @(#)
100 c FILE: svd.F SID: 2.4 DATE OF SID: 10/17/00 RELEASE: 2
101 c
102 c\Remarks
103 c 1. None
104 c
105 c\EndLib
106 c
107 c-----------------------------------------------------------------------
108 c
109 c %------------------------------------------------------%
110 c | Storage Declarations: |
111 c | |
112 c | It is assumed that A is M by N with M .ge. N. |
113 c | |
114 c | The maximum dimensions for all arrays are |
115 c | set here to accommodate a problem size of |
116 c | M .le. MAXM and N .le. MAXN |
117 c | |
118 c | The NEV right singular vectors will be computed in |
119 c | the N by NCV array V. |
120 c | |
121 c | The NEV left singular vectors will be computed in |
122 c | the M by NEV array U. |
123 c | |
124 c | NEV is the number of singular values requested. |
125 c | See specifications for ARPACK usage below. |
126 c | |
127 c | NCV is the largest number of basis vectors that will |
128 c | be used in the Implicitly Restarted Arnoldi |
129 c | Process. Work per major iteration is |
130 c | proportional to N*NCV*NCV. |
131 c | |
132 c | You must set: |
133 c | |
134 c | MAXM: Maximum number of rows of the A allowed. |
135 c | MAXN: Maximum number of columns of the A allowed. |
136 c | MAXNEV: Maximum NEV allowed |
137 c | MAXNCV: Maximum NCV allowed |
138 c %------------------------------------------------------%
139 c
140 C !USES:
141 C == Global variables ===
142 #include "SIZE.h"
143 #include "EEPARAMS.h"
144 #include "PARAMS.h"
145 #ifdef ALLOW_ADMTLM
146 # include "tamc.h"
147 # include "ctrl.h"
148 # include "optim.h"
149 #endif
150
151 C !INPUT/OUTPUT PARAMETERS:
152 C == Routine arguments ==
153 INTEGER mythid
154
155 C == PARAMETERS ==
156 integer maxm, maxn, maxnev, maxncv, ldv, ldu
157 cph(
158 parameter (maxm = 2*admtlmrec, maxn=2*admtlmrec,
159 & maxnev=15, maxncv=30,
160 & ldu = maxm, ldv=maxn )
161 character*(80) fnameGlobal
162 cph)
163 c
164 c %--------------%
165 c | Local Arrays |
166 c %--------------%
167 c
168 Double precision
169 & v(ldv,maxncv), u(ldu, maxnev),
170 & workl(maxncv*(maxncv+8)), workd(3*maxn),
171 & s(maxncv,2), resid(maxn), ax(maxm)
172 logical select(maxncv)
173 integer iparam(11), ipntr(11)
174 c
175 c %---------------%
176 c | Local Scalars |
177 c %---------------%
178 c
179 character bmat*1, which*2
180 integer ido, m, n, nev, ncv, lworkl, info, ierr,
181 & j, ishfts, maxitr, mode, nconv
182 logical rvec
183 Double precision
184 & tol, sigma, temp
185 c
186 c %------------%
187 c | Parameters |
188 c %------------%
189 c
190 Double precision
191 & one, zero
192 parameter (one = 1.0D+0, zero = 0.0D+0)
193 c
194 c %-----------------------------%
195 c | BLAS & LAPACK routines used |
196 c %-----------------------------%
197 c
198 Double precision
199 & dnrm2
200 external dnrm2, daxpy, dcopy, dscal
201
202 cph(
203 integer l, ilinsysinfo, ii, jj
204 integer ipiv(maxn)
205 integer phiwork(maxn)
206
207 double precision ferr, berr
208 double precision phtmpin(maxn), phtmpout(maxn)
209 double precision phxout(maxn)
210 double precision tmpvec1(maxn), tmpvec2(maxn)
211 double precision phrwork(3*maxn)
212 cph double precision metricLoc(maxn,maxn)
213 cph double precision metricTriag(maxn,maxn)
214 cph double precision metricInv(maxn,maxn)
215
216 cph DATA (metricInv(ii,1),ii=1,maxn)
217 cph & / 9.9896D+07, 0.0519D+07, -0.0415D+07,
218 cph & 0.0486D+07, -0.2432D+07, 0.1945D+07 /
219 cph DATA (metricInv(ii,2),ii=1,maxn)
220 cph & / 0.0519D+07, 9.7403D+07, 0.2077D+07,
221 cph & -0.2432D+07, 1.2158D+07, -0.9726D+07 /
222 cph DATA (metricInv(ii,3),ii=1,maxn)
223 cph & / -0.0415D+07, 0.2077D+07, 9.8338D+07,
224 cph & 0.1945D+07, -0.9726D+07, 0.7781D+07 /
225 cph DATA (metricInv(ii,4),ii=1,maxn)
226 cph & / 0.0486D+07, -0.2432D+07, 0.1945D+07,
227 cph & 9.7723D+07, 1.1385D+07, -0.9108D+07 /
228 cph DATA (metricInv(ii,5),ii=1,maxn)
229 cph & / -0.2432D+07, 1.2158D+07, -0.9726D+07,
230 cph & 1.1385D+07, 4.3073D+07, 4.5542D+07 /
231 cph DATA (metricInv(ii,6),ii=1,maxn)
232 cph & / 0.1945D+07, -0.9726D+07, 0.7781D+07,
233 cph & -0.9108D+07, 4.5542D+07, 6.3567D+07 /
234 cph)
235 c
236 c %-----------------------%
237 c | Executable Statements |
238 c %-----------------------%
239 c
240 c %-------------------------------------------------%
241 c | The following include statement and assignments |
242 c | initiate trace output from the internal |
243 c | actions of ARPACK. See debug.doc in the |
244 c | DOCUMENTS directory for usage. Initially, the |
245 c | most useful information will be a breakdown of |
246 c | time spent in the various stages of computation |
247 c | given by setting msaupd = 1. |
248 c %-------------------------------------------------%
249 c
250 include 'arpack_debug.h'
251 ndigit = -3
252 logfil = 6
253 msgets = 0
254 msaitr = 0
255 msapps = 0
256 msaupd = 1
257 msaup2 = 0
258 mseigt = 0
259 mseupd = 0
260 c
261 c %-------------------------------------------------%
262 c | The following sets dimensions for this problem. |
263 c %-------------------------------------------------%
264 c
265 cph(
266 m = 2*admtlmrec
267 n = 2*admtlmrec
268 cph)
269 c
270 c %------------------------------------------------%
271 c | Specifications for ARPACK usage are set |
272 c | below: |
273 c | |
274 c | 1) NEV = 4 asks for 4 singular values to be |
275 c | computed. |
276 c | |
277 c | 2) NCV = 20 sets the length of the Arnoldi |
278 c | factorization |
279 c | |
280 c | 3) This is a standard problem |
281 c | (indicated by bmat = 'I') |
282 c | |
283 c | 4) Ask for the NEV singular values of |
284 c | largest magnitude |
285 c | (indicated by which = 'LM') |
286 c | See documentation in DSAUPD for the |
287 c | other options SM, BE. |
288 c | |
289 c | Note: NEV and NCV must satisfy the following |
290 c | conditions: |
291 c | NEV <= MAXNEV, |
292 c | NEV + 1 <= NCV <= MAXNCV |
293 c %------------------------------------------------%
294 c
295 cph(
296 nev = 15
297 ncv = 30
298 bmat = 'I'
299 cph)
300 which = 'LM'
301 c
302 if ( n .gt. maxn ) then
303 print *, ' ERROR with _SVD: N is greater than MAXN '
304 go to 9000
305 else if ( m .gt. maxm ) then
306 print *, ' ERROR with _SVD: M is greater than MAXM '
307 go to 9000
308 else if ( nev .gt. maxnev ) then
309 print *, ' ERROR with _SVD: NEV is greater than MAXNEV '
310 go to 9000
311 else if ( ncv .gt. maxncv ) then
312 print *, ' ERROR with _SVD: NCV is greater than MAXNCV '
313 go to 9000
314 end if
315 c
316 c %-----------------------------------------------------%
317 c | Specification of stopping rules and initial |
318 c | conditions before calling DSAUPD |
319 c | |
320 c | abs(sigmaC - sigmaT) < TOL*abs(sigmaC) |
321 c | computed true |
322 c | |
323 c | If TOL .le. 0, then TOL <- macheps |
324 c | (machine precision) is used. |
325 c | |
326 c | IDO is the REVERSE COMMUNICATION parameter |
327 c | used to specify actions to be taken on return |
328 c | from DSAUPD. (See usage below.) |
329 c | |
330 c | It MUST initially be set to 0 before the first |
331 c | call to DSAUPD. |
332 c | |
333 c | INFO on entry specifies starting vector information |
334 c | and on return indicates error codes |
335 c | |
336 c | Initially, setting INFO=0 indicates that a |
337 c | random starting vector is requested to |
338 c | start the ARNOLDI iteration. Setting INFO to |
339 c | a nonzero value on the initial call is used |
340 c | if you want to specify your own starting |
341 c | vector (This vector must be placed in RESID.) |
342 c | |
343 c | The work array WORKL is used in DSAUPD as |
344 c | workspace. Its dimension LWORKL is set as |
345 c | illustrated below. |
346 c %-----------------------------------------------------%
347 c
348 lworkl = ncv*(ncv+8)
349 cph(
350 tol = zero
351 cph tol = 1.D-10
352 cph)
353 info = 0
354 ido = 0
355 c
356 c %---------------------------------------------------%
357 c | Specification of Algorithm Mode: |
358 c | |
359 c | This program uses the exact shift strategy |
360 c | (indicated by setting IPARAM(1) = 1.) |
361 c | IPARAM(3) specifies the maximum number of Arnoldi |
362 c | iterations allowed. Mode 1 of DSAUPD is used |
363 c | (IPARAM(7) = 1). All these options can be changed |
364 c | by the user. For details see the documentation in |
365 c | DSAUPD. |
366 c %---------------------------------------------------%
367 c
368 ishfts = 1
369 cph(
370 maxitr = 10
371 c maxitr = 5
372 c
373 mode = 1
374 cph)
375 iparam(1) = ishfts
376 c
377 iparam(3) = maxitr
378 c
379 iparam(7) = mode
380 c
381 c %------------------------------------------------%
382 c | M A I N L O O P (Reverse communication loop) |
383 c %------------------------------------------------%
384 c
385 C-- Set model configuration (fixed arrays)
386 CALL INITIALISE_FIXED( myThid )
387 c
388 10 continue
389 c
390 print *, 'ph----------------------------------------------------'
391 print *, 'ph----------------------------------------------------'
392 print *, 'ph----------------------------------------------------'
393 c
394 c %---------------------------------------------%
395 c | Repeatedly call the routine DSAUPD and take |
396 c | actions indicated by parameter IDO until |
397 c | either convergence is indicated or maxitr |
398 c | has been exceeded. |
399 c %---------------------------------------------%
400 c
401 cph(
402 cph if ( iphcount .GT. 3 ) then
403 cph endif
404 cph)
405 CALL DSAUPD ( ido, bmat, n, which, nev, tol, resid,
406 & ncv, v, ldv, iparam, ipntr, workd, workl,
407 & lworkl, info )
408 c
409 cph(
410 print *, 'ph-count: optimcycle, ido, info, ipntr(1) ',
411 & optimcycle, ido, info, ipntr(1)
412 cph)
413
414 if (ido .eq. -1 .or. ido .eq. 1) then
415 c
416 c %---------------------------------------%
417 c | Perform matrix vector multiplications |
418 c | w <--- A*x (av()) |
419 c | y <--- A'*w (atv()) |
420 c | The user should supply his/her own |
421 c | matrix vector multiplication routines |
422 c | here that takes workd(ipntr(1)) as |
423 c | the input, and returns the result in |
424 c | workd(ipntr(2)). |
425 c %---------------------------------------%
426 c
427 c call av (m, n, workd(ipntr(1)), ax)
428 c call atv (m, n, ax, workd(ipntr(2)))
429 c
430 cph(
431 c do l = 1, n
432 c print *, 'ph-test A ', l,
433 c & workd(ipntr(1)+l), workd(ipntr(2)+l), workd(l)
434 c enddo
435 cph)
436 do l = 1, n
437 phtmpin(l) = workd(ipntr(1)+l)
438 enddo
439
440 if (optimcycle .GT. 0 ) then
441 write(fnameGlobal(1:80),'(a)') ' '
442 write(fnameGlobal,'(a,i4.4)')
443 & 'admtlm_vector.it', optimcycle
444 call mdswritevector( fnameGlobal, 64, .false., 'RL',
445 & 2*admtlmrec, phtmpin, 1, 1, 1, 0, mythid )
446 endif
447
448 c-- MWMWMWMWMWMWMW
449 CALL ADMTLM_DRIVER( mythid )
450 cph call box_main( phtmpin, phtmpout, metricLoc, ldoadjoint )
451 c-- MWMWMWMWMWMWMW
452
453 write(fnameGlobal(1:80),'(a)') ' '
454 write(fnameGlobal,'(a,i4.4)')
455 & 'admtlm_vector.it', optimcycle
456 call mdsreadvector( fnameGlobal, 64, 'RL',
457 & 2*admtlmrec, phtmpout, 1, 1, 1, mythid )
458
459 do l = 1, n
460 workd(ipntr(2)+l) = phtmpout(l)
461 enddo
462 c
463
464 if (optimcycle .EQ. 8)
465 & STOP 'in ADMTLM_DSVD after the_model_main'
466
467 cph(
468 cph Since we solve for a generalized EV problem, we have to solve
469 cph M*y=A*x for y with known matrix M and vector A*x
470 c do ii = 1, n
471 c phxout(ii) = 0.
472 c do jj = 1, n
473 c phxout(ii) = phxout(ii) +
474 c & metricInv(ii,jj)*phtmpout(jj)
475 c enddo
476 c print *, 'ph-test C ', ii, phxout(ii)
477 c enddo
478 c
479 c call dposvx (
480 c & 'E', 'U', maxn, 1, metricLoc, 6, metricTriag, 6,
481 c & 'N', tmpvec1, phtmpout, 6, phxout, 6, tmpvec2,
482 c & ferr, berr, phrwork, phiwork, ilinsysinfo)
483 cph
484 c call dgesv ( maxn, 1, metricLoc, maxn,
485 c & ipiv, phtmpout, maxn, ilinsysinfo )
486 cph
487 c
488 cph Finally schift result y -> workd(ipntr(2))
489 cph call dcopy ( maxn, tmpvec1, 1, workd(ipntr(2)), 1 )
490 cph We have restored orig. standard EV problem
491 cph OP*x = lambda*x for OP = INV(M)*A
492 cph)
493 c
494 c %-----------------------------------------%
495 c | L O O P B A C K to call DSAUPD again. |
496 c %-----------------------------------------%
497 c
498 optimcycle = optimcycle + 1
499
500 go to 10
501 c
502 end if
503 c
504 cph(
505 1001 continue
506 print *, 'ph-continue ', info, ierr
507 cph)
508 c
509 c %----------------------------------------%
510 c | Either we have convergence or there is |
511 c | an error. |
512 c %----------------------------------------%
513 c
514 if ( info .lt. 0 ) then
515 c
516 c %--------------------------%
517 c | Error message. Check the |
518 c | documentation in DSAUPD. |
519 c %--------------------------%
520 c
521 print *, ' '
522 print *, ' Error with _saupd, info = ', info
523 print *, ' Check documentation in _saupd '
524 print *, ' '
525 c
526 else
527 c
528 c %--------------------------------------------%
529 c | No fatal errors occurred. |
530 c | Post-Process using DSEUPD. |
531 c | |
532 c | Computed singular values may be extracted. |
533 c | |
534 c | Singular vectors may also be computed now |
535 c | if desired. (indicated by rvec = .true.) |
536 c | |
537 c | The routine DSEUPD now called to do this |
538 c | post processing |
539 c %--------------------------------------------%
540 c
541 rvec = .true.
542 c
543 call dseupd ( rvec, 'All', select, s, v, ldv, sigma,
544 & bmat, n, which, nev, tol, resid, ncv, v, ldv,
545 & iparam, ipntr, workd, workl, lworkl, ierr )
546 c
547 c %-----------------------------------------------%
548 c | Singular values are returned in the first |
549 c | column of the two dimensional array S |
550 c | and the corresponding right singular vectors |
551 c | are returned in the first NEV columns of the |
552 c | two dimensional array V as requested here. |
553 c %-----------------------------------------------%
554 c
555 if ( ierr .ne. 0) then
556 c
557 c %------------------------------------%
558 c | Error condition: |
559 c | Check the documentation of DSEUPD. |
560 c %------------------------------------%
561 c
562 print *, ' '
563 print *, ' Error with _seupd, info = ', ierr
564 print *, ' Check the documentation of _seupd. '
565 print *, ' '
566 c
567 else
568 c
569 nconv = iparam(5)
570 cph(
571 do j=1, nconv
572 print *, 'ph-ev ', j, s(j,1), s(j,2)
573 enddo
574
575 STOP 'TEST AFTER dneupd'
576 cph)
577 do 20 j=1, nconv
578 c
579 s(j,1) = sqrt(s(j,1))
580 c
581 c %-----------------------------%
582 c | Compute the left singular |
583 c | vectors from the formula |
584 c | |
585 c | u = Av/sigma |
586 c | |
587 c | u should have norm 1 so |
588 c | divide by norm(Av) instead. |
589 c %-----------------------------%
590 c
591 do l = 1, n
592 tmpvec1(l) = v(l,j)
593 enddo
594
595 c
596 c-- MWMWMWMWMWMWMW
597 cph call box_main( tmpvec1, tmpvec2, metricLoc, ldoadjoint )
598 c-- MWMWMWMWMWMWMW
599 c
600 call dcopy( m, tmpvec2, 1, u(l,j), 1 )
601 temp = one / dnrm2( m, u(l,j), 1 )
602 cph(
603 print *, 'ph-print ', j, dnrm2( m, u(l,j), 1 )
604 cph)
605 call dscal(m , temp, u(l,j), 1 )
606
607 cph do l = 1, n
608 cph u(l,j) = tmpvec2(l)
609 cph enddo
610 cph temp = one/dnrm2(m, tmpvec2, 1)
611 cph call dscal(m, temp, tmpvec2, 1)
612 c
613 c %---------------------------%
614 c | |
615 c | Compute the residual norm |
616 c | |
617 c | || A*v - sigma*u || |
618 c | |
619 c | for the NCONV accurately |
620 c | computed singular values |
621 c | and vectors. (iparam(5) |
622 c | indicates how many are |
623 c | accurate to the requested |
624 c | tolerance). |
625 c | Store the result in 2nd |
626 c | column of array S. |
627 c %---------------------------%
628 c
629 call daxpy(m, -s(j,1), u(1,j), 1, tmpvec2, 1)
630 s(j,2) = dnrm2(m, tmpvec2, 1)
631 c
632 20 continue
633 c
634 c %-------------------------------%
635 c | Display computed residuals |
636 c %-------------------------------%
637 c
638 call dmout(6, nconv, 2, s, maxncv, -6,
639 & 'Singular values and direct residuals')
640 end if
641 c
642 c %------------------------------------------%
643 c | Print additional convergence information |
644 c %------------------------------------------%
645 c
646 if ( info .eq. 1) then
647 print *, ' '
648 print *, ' Maximum number of iterations reached.'
649 print *, ' '
650 else if ( info .eq. 3) then
651 print *, ' '
652 print *, ' No shifts could be applied during implicit',
653 & ' Arnoldi update, try increasing NCV.'
654 print *, ' '
655 end if
656 c
657 print *, ' '
658 print *, ' _SVD '
659 print *, ' ==== '
660 print *, ' '
661 print *, ' Size of the matrix is ', n
662 print *, ' The number of Ritz values requested is ', nev
663 print *, ' The number of Arnoldi vectors generated',
664 & ' (NCV) is ', ncv
665 print *, ' What portion of the spectrum: ', which
666 print *, ' The number of converged Ritz values is ',
667 & nconv
668 print *, ' The number of Implicit Arnoldi update',
669 & ' iterations taken is ', iparam(3)
670 print *, ' The number of OP*x is ', iparam(9)
671 print *, ' The convergence criterion is ', tol
672 print *, ' '
673 c
674 end if
675 c
676 c %-------------------------%
677 c | Done with program dsvd. |
678 c %-------------------------%
679 c
680 9000 continue
681 c
682 end

  ViewVC Help
Powered by ViewVC 1.1.22