Actual source code: svdsolve.c
1: /*
2: SVD routines related to the solution process.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain
8: This file is part of SLEPc. See the README file for conditions of use
9: and additional information.
10: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
11: */
13: #include src/svd/svdimpl.h
17: /*@
18: SVDSolve - Solves the singular value problem.
20: Collective on SVD
22: Input Parameter:
23: . svd - singular value solver context obtained from SVDCreate()
25: Options Database:
26: . -svd_view - print information about the solver used
28: Level: beginner
30: .seealso: SVDCreate(), SVDSetUp(), SVDDestroy()
31: @*/
32: PetscErrorCode SVDSolve(SVD svd)
33: {
35: PetscTruth flg;
36: int i;
41: if (!svd->setupcalled) { SVDSetUp(svd); }
42: svd->its = 0;
43: svd->matvecs = 0;
44: svd->nconv = 0;
45: svd->reason = SVD_CONVERGED_ITERATING;
46: IPResetOperationCounters(svd->ip);
47: for (i=0;i<svd->ncv;i++) svd->sigma[i]=svd->errest[i]=0.0;
48: SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,svd->ncv);
50: PetscLogEventBegin(SVD_Solve,svd,0,0,0);
51: (*svd->ops->solve)(svd);
52: PetscLogEventEnd(SVD_Solve,svd,0,0,0);
54: PetscOptionsHasName(svd->prefix,"-svd_view",&flg);
55: if (flg && !PetscPreLoadingOn) { SVDView(svd,PETSC_VIEWER_STDOUT_WORLD); }
57: return(0);
58: }
62: /*@
63: SVDGetIterationNumber - Gets the current iteration number. If the
64: call to SVDSolve() is complete, then it returns the number of iterations
65: carried out by the solution method.
66:
67: Not Collective
69: Input Parameter:
70: . svd - the singular value solver context
72: Output Parameter:
73: . its - number of iterations
75: Level: intermediate
77: Notes:
78: During the i-th iteration this call returns i-1. If SVDSolve() is
79: complete, then parameter "its" contains either the iteration number at
80: which convergence was successfully reached, or failure was detected.
81: Call SVDGetConvergedReason() to determine if the solver converged or
82: failed and why.
84: @*/
85: PetscErrorCode SVDGetIterationNumber(SVD svd,int *its)
86: {
90: *its = svd->its;
91: return(0);
92: }
96: /*@C
97: SVDGetConvergedReason - Gets the reason why the SVDSolve() iteration was
98: stopped.
100: Not Collective
102: Input Parameter:
103: . svd - the singular value solver context
105: Output Parameter:
106: . reason - negative value indicates diverged, positive value converged
107: (see SVDConvergedReason)
109: Possible values for reason:
110: + SVD_CONVERGED_TOL - converged up to tolerance
111: . SVD_DIVERGED_ITS - required more than its to reach convergence
112: - SVD_DIVERGED_BREAKDOWN - generic breakdown in method
114: Level: intermediate
116: Notes: Can only be called after the call to SVDSolve() is complete.
118: .seealso: SVDSetTolerances(), SVDSolve(), SVDConvergedReason
119: @*/
120: PetscErrorCode SVDGetConvergedReason(SVD svd,SVDConvergedReason *reason)
121: {
125: *reason = svd->reason;
126: return(0);
127: }
131: /*@
132: SVDGetConverged - Gets the number of converged singular values.
134: Not Collective
136: Input Parameter:
137: . svd - the singular value solver context
138:
139: Output Parameter:
140: . nconv - number of converged singular values
142: Note:
143: This function should be called after SVDSolve() has finished.
145: Level: beginner
147: @*/
148: PetscErrorCode SVDGetConverged(SVD svd,int *nconv)
149: {
153: if (svd->reason == SVD_CONVERGED_ITERATING) {
154: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "SVDSolve must be called first");
155: }
156: *nconv = svd->nconv;
157: return(0);
158: }
162: /*@
163: SVDGetSingularTriplet - Gets the i-th triplet of the singular value decomposition
164: as computed by SVDSolve(). The solution consists in the singular value and its left
165: and right singular vectors.
167: Not Collective
169: Input Parameters:
170: + svd - singular value solver context
171: - i - index of the solution
173: Output Parameters:
174: + sigma - singular value
175: . u - left singular vector
176: - v - right singular vector
178: Note:
179: The index i should be a value between 0 and nconv-1 (see SVDGetConverged()).
180: Both U or V can be PETSC_NULL if singular vectors are not required.
182: Level: beginner
184: .seealso: SVDSolve(), SVDGetConverged()
185: @*/
186: PetscErrorCode SVDGetSingularTriplet(SVD svd, int i, PetscReal *sigma, Vec u, Vec v)
187: {
193: if (svd->reason == SVD_CONVERGED_ITERATING) {
194: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "SVDSolve must be called first");
195: }
196: if (i<0 || i>=svd->nconv) {
197: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
198: }
199: *sigma = svd->sigma[i];
200: if (u) {
202: if (svd->U) {
203: VecCopy(svd->U[i],u);
204: } else {
205: SVDMatMult(svd,PETSC_FALSE,svd->V[i],u);
206: VecScale(u,1.0/svd->sigma[i]);
207: }
208: }
209: if (v) {
211: VecCopy(svd->V[i],v);
212: }
213: return(0);
214: }
218: /*@
219: SVDComputeResidualNorms - Computes the norms of the residual vectors associated with
220: the i-th computed singular triplet.
222: Collective on SVD
224: Input Parameters:
225: + svd - the singular value solver context
226: - i - the solution index
228: Output Parameters:
229: + norm1 - the norm ||A*v-sigma*u||_2 where sigma is the
230: singular value, u and v are the left and right singular vectors.
231: - norm2 - the norm ||A^T*u-sigma*v||_2 with the same sigma, u and v
233: Note:
234: The index i should be a value between 0 and nconv-1 (see SVDGetConverged()).
235: Both output parameters can be PETSC_NULL on input if not needed.
237: Level: beginner
239: .seealso: SVDSolve(), SVDGetConverged(), SVDComputeRelativeError()
240: @*/
241: PetscErrorCode SVDComputeResidualNorms(SVD svd, int i, PetscReal *norm1, PetscReal *norm2)
242: {
244: Vec u,v,x = PETSC_NULL,y = PETSC_NULL;
245: PetscReal sigma;
249: if (svd->reason == SVD_CONVERGED_ITERATING) {
250: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "SVDSolve must be called first");
251: }
252: if (i<0 || i>=svd->nconv) {
253: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
254: }
255:
256: SVDMatGetVecs(svd,&v,&u);
257: SVDGetSingularTriplet(svd,i,&sigma,u,v);
258: if (norm1) {
259: VecDuplicate(u,&x);
260: SVDMatMult(svd,PETSC_FALSE,v,x);
261: VecAXPY(x,-sigma,u);
262: VecNorm(x,NORM_2,norm1);
263: }
264: if (norm2) {
265: VecDuplicate(v,&y);
266: SVDMatMult(svd,PETSC_TRUE,u,y);
267: VecAXPY(y,-sigma,v);
268: VecNorm(y,NORM_2,norm2);
269: }
271: VecDestroy(v);
272: VecDestroy(u);
273: if (x) { VecDestroy(x); }
274: if (y) { VecDestroy(y); }
275: return(0);
276: }
280: /*@
281: SVDComputeRelativeError - Computes the relative error bound associated
282: with the i-th singular triplet.
284: Collective on SVD
286: Input Parameter:
287: + svd - the singular value solver context
288: - i - the solution index
290: Output Parameter:
291: . error - the relative error bound, computed as sqrt(n1^2+n2^2)/(sigma*sqrt(||u||_2^2+||v||_2^2))
292: where n1 = ||A*v-sigma*u||_2 , n2 = ||A^T*u-sigma*v||_2 , sigma is the singular value,
293: u and v are the left and right singular vectors.
294: If sigma is too small the relative error is computed as sqrt(n1^2+n2^2)/sqrt(||u||_2^2+||v||_2^2).
296: Level: beginner
298: .seealso: SVDSolve(), SVDComputeResidualNorms()
299: @*/
300: PetscErrorCode SVDComputeRelativeError(SVD svd, int i, PetscReal *error)
301: {
303: PetscReal sigma,norm1,norm2,norm3,norm4;
304: Vec u,v;
309: SVDMatGetVecs(svd,&v,&u);
310: SVDGetSingularTriplet(svd,i,&sigma,u,v);
311: SVDComputeResidualNorms(svd,i,&norm1,&norm2);
312: VecNorm(u,NORM_2,&norm3);
313: VecNorm(v,NORM_2,&norm4);
314: *error = sqrt(norm1*norm1+norm2*norm2) / sqrt(norm3*norm3+norm4*norm4);
315: if (sigma>*error) *error /= sigma;
316: VecDestroy(v);
317: VecDestroy(u);
318: return(0);
319: }
323: /*@
324: SVDGetOperationCounters - Gets the total number of matrix vector and dot
325: products used by the SVD object during the last SVDSolve() call.
327: Not Collective
329: Input Parameter:
330: . svd - SVD context
332: Output Parameter:
333: + matvecs - number of matrix vector product operations
334: - dots - number of dot product operations
336: Notes:
337: These counters are reset to zero at each successive call to SVDSolve().
339: Level: intermediate
341: @*/
342: PetscErrorCode SVDGetOperationCounters(SVD svd,int* matvecs,int* dots)
343: {
345:
348: if (matvecs) *matvecs = svd->matvecs;
349: if (dots) {
350: IPGetOperationCounters(svd->ip,dots);
351: }
352: return(0);
353: }