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-2012, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
9:
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
24: #include <slepc-private/svdimpl.h> /*I "slepcsvd.h" I*/
28: /*@
29: SVDSolve - Solves the singular value problem.
31: Collective on SVD
33: Input Parameter:
34: . svd - singular value solver context obtained from SVDCreate()
36: Options Database:
37: . -svd_view - print information about the solver used
39: Level: beginner
41: .seealso: SVDCreate(), SVDSetUp(), SVDDestroy()
42: @*/
43: PetscErrorCode SVDSolve(SVD svd)
44: {
46: PetscBool flg;
47: PetscInt i,*workperm;
48: char filename[PETSC_MAX_PATH_LEN];
49: PetscViewer viewer;
50: PetscErrorCode (*which_func)(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
54: if (!svd->setupcalled) { SVDSetUp(svd); }
55: svd->its = 0;
56: svd->nconv = 0;
57: svd->reason = SVD_CONVERGED_ITERATING;
58: for (i=0;i<svd->ncv;i++) svd->sigma[i]=svd->errest[i]=0.0;
59: SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,svd->ncv);
61: which_func = (svd->which==SVD_LARGEST)? SlepcCompareLargestReal: SlepcCompareSmallestReal;
62: DSSetEigenvalueComparison(svd->ds,which_func,PETSC_NULL);
64: PetscLogEventBegin(SVD_Solve,svd,0,0,0);
65: (*svd->ops->solve)(svd);
66: PetscLogEventEnd(SVD_Solve,svd,0,0,0);
68: /* sort singular triplets */
69: if (svd->which == SVD_SMALLEST) {
70: for (i=0;i<svd->nconv;i++) svd->perm[i] = i;
71: PetscSortRealWithPermutation(svd->nconv,svd->sigma,svd->perm);
72: } else {
73: PetscMalloc(sizeof(PetscInt)*svd->nconv,&workperm);
74: for (i=0;i<svd->nconv;i++) workperm[i] = i;
75: PetscSortRealWithPermutation(svd->nconv,svd->sigma,workperm);
76: for (i=0;i<svd->nconv;i++) svd->perm[i] = workperm[svd->nconv-i-1];
77: PetscFree(workperm);
78: }
80: PetscOptionsGetString(((PetscObject)svd)->prefix,"-svd_view",filename,PETSC_MAX_PATH_LEN,&flg);
81: if (flg && !PetscPreLoadingOn) {
82: PetscViewerASCIIOpen(((PetscObject)svd)->comm,filename,&viewer);
83: SVDView(svd,viewer);
84: PetscViewerDestroy(&viewer);
85: }
87: /* Remove the initial subspace */
88: svd->nini = 0;
89: return(0);
90: }
94: /*@
95: SVDGetIterationNumber - Gets the current iteration number. If the
96: call to SVDSolve() is complete, then it returns the number of iterations
97: carried out by the solution method.
98:
99: Not Collective
101: Input Parameter:
102: . svd - the singular value solver context
104: Output Parameter:
105: . its - number of iterations
107: Level: intermediate
109: Notes:
110: During the i-th iteration this call returns i-1. If SVDSolve() is
111: complete, then parameter "its" contains either the iteration number at
112: which convergence was successfully reached, or failure was detected.
113: Call SVDGetConvergedReason() to determine if the solver converged or
114: failed and why.
116: @*/
117: PetscErrorCode SVDGetIterationNumber(SVD svd,PetscInt *its)
118: {
122: *its = svd->its;
123: return(0);
124: }
128: /*@C
129: SVDGetConvergedReason - Gets the reason why the SVDSolve() iteration was
130: stopped.
132: Not Collective
134: Input Parameter:
135: . svd - the singular value solver context
137: Output Parameter:
138: . reason - negative value indicates diverged, positive value converged
139: (see SVDConvergedReason)
141: Possible values for reason:
142: + SVD_CONVERGED_TOL - converged up to tolerance
143: . SVD_DIVERGED_ITS - required more than its to reach convergence
144: - SVD_DIVERGED_BREAKDOWN - generic breakdown in method
146: Level: intermediate
148: Notes: Can only be called after the call to SVDSolve() is complete.
150: .seealso: SVDSetTolerances(), SVDSolve(), SVDConvergedReason
151: @*/
152: PetscErrorCode SVDGetConvergedReason(SVD svd,SVDConvergedReason *reason)
153: {
157: *reason = svd->reason;
158: return(0);
159: }
163: /*@
164: SVDGetConverged - Gets the number of converged singular values.
166: Not Collective
168: Input Parameter:
169: . svd - the singular value solver context
170:
171: Output Parameter:
172: . nconv - number of converged singular values
174: Note:
175: This function should be called after SVDSolve() has finished.
177: Level: beginner
179: @*/
180: PetscErrorCode SVDGetConverged(SVD svd,PetscInt *nconv)
181: {
185: *nconv = svd->nconv;
186: return(0);
187: }
191: /*@
192: SVDGetSingularTriplet - Gets the i-th triplet of the singular value decomposition
193: as computed by SVDSolve(). The solution consists in the singular value and its left
194: and right singular vectors.
196: Not Collective, but vectors are shared by all processors that share the SVD
198: Input Parameters:
199: + svd - singular value solver context
200: - i - index of the solution
202: Output Parameters:
203: + sigma - singular value
204: . u - left singular vector
205: - v - right singular vector
207: Note:
208: The index i should be a value between 0 and nconv-1 (see SVDGetConverged()).
209: Both U or V can be PETSC_NULL if singular vectors are not required.
211: Level: beginner
213: .seealso: SVDSolve(), SVDGetConverged()
214: @*/
215: PetscErrorCode SVDGetSingularTriplet(SVD svd,PetscInt i,PetscReal *sigma,Vec u,Vec v)
216: {
218: PetscReal norm;
219: PetscInt j,M,N;
220: Vec w;
226: if (svd->reason == SVD_CONVERGED_ITERATING) SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_WRONGSTATE,"SVDSolve must be called first");
227: if (i<0 || i>=svd->nconv) SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
228: *sigma = svd->sigma[svd->perm[i]];
229: MatGetSize(svd->OP,&M,&N);
230: if (M<N) { w = u; u = v; v = w; }
231: if (u) {
232: if (!svd->U) {
233: VecDuplicateVecs(svd->tl,svd->ncv,&svd->U);
234: for (j=0;j<svd->nconv;j++) {
235: SVDMatMult(svd,PETSC_FALSE,svd->V[j],svd->U[j]);
236: IPOrthogonalize(svd->ip,0,PETSC_NULL,j,PETSC_NULL,svd->U,svd->U[j],PETSC_NULL,&norm,PETSC_NULL);
237: VecScale(svd->U[j],1.0/norm);
238: }
239: }
240: VecCopy(svd->U[svd->perm[i]],u);
241: }
242: if (v) {
243: VecCopy(svd->V[svd->perm[i]],v);
244: }
245: return(0);
246: }
250: /*@
251: SVDComputeResidualNorms - Computes the norms of the residual vectors associated with
252: the i-th computed singular triplet.
254: Collective on SVD
256: Input Parameters:
257: + svd - the singular value solver context
258: - i - the solution index
260: Output Parameters:
261: + norm1 - the norm ||A*v-sigma*u||_2 where sigma is the
262: singular value, u and v are the left and right singular vectors.
263: - norm2 - the norm ||A^T*u-sigma*v||_2 with the same sigma, u and v
265: Note:
266: The index i should be a value between 0 and nconv-1 (see SVDGetConverged()).
267: Both output parameters can be PETSC_NULL on input if not needed.
269: Level: beginner
271: .seealso: SVDSolve(), SVDGetConverged(), SVDComputeRelativeError()
272: @*/
273: PetscErrorCode SVDComputeResidualNorms(SVD svd,PetscInt i,PetscReal *norm1,PetscReal *norm2)
274: {
276: Vec u,v,x = PETSC_NULL,y = PETSC_NULL;
277: PetscReal sigma;
278: PetscInt M,N;
283: if (svd->reason == SVD_CONVERGED_ITERATING) SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_WRONGSTATE,"SVDSolve must be called first");
284: if (i<0 || i>=svd->nconv) SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
285:
286: MatGetVecs(svd->OP,&v,&u);
287: SVDGetSingularTriplet(svd,i,&sigma,u,v);
288: if (norm1) {
289: VecDuplicate(u,&x);
290: MatMult(svd->OP,v,x);
291: VecAXPY(x,-sigma,u);
292: VecNorm(x,NORM_2,norm1);
293: }
294: if (norm2) {
295: VecDuplicate(v,&y);
296: if (svd->A && svd->AT) {
297: MatGetSize(svd->OP,&M,&N);
298: if (M<N) {
299: MatMult(svd->A,u,y);
300: } else {
301: MatMult(svd->AT,u,y);
302: }
303: } else {
304: #if defined(PETSC_USE_COMPLEX)
305: MatMultHermitianTranspose(svd->OP,u,y);
306: #else
307: MatMultTranspose(svd->OP,u,y);
308: #endif
309: }
310: VecAXPY(y,-sigma,v);
311: VecNorm(y,NORM_2,norm2);
312: }
314: VecDestroy(&v);
315: VecDestroy(&u);
316: VecDestroy(&x);
317: VecDestroy(&y);
318: return(0);
319: }
323: /*@
324: SVDComputeRelativeError - Computes the relative error bound associated
325: with the i-th singular triplet.
327: Collective on SVD
329: Input Parameter:
330: + svd - the singular value solver context
331: - i - the solution index
333: Output Parameter:
334: . error - the relative error bound, computed as sqrt(n1^2+n2^2)/sigma
335: where n1 = ||A*v-sigma*u||_2 , n2 = ||A^T*u-sigma*v||_2 , sigma is the singular value,
336: u and v are the left and right singular vectors.
337: If sigma is too small the relative error is computed as sqrt(n1^2+n2^2).
339: Level: beginner
341: .seealso: SVDSolve(), SVDComputeResidualNorms()
342: @*/
343: PetscErrorCode SVDComputeRelativeError(SVD svd,PetscInt i,PetscReal *error)
344: {
346: PetscReal sigma,norm1,norm2;
352: SVDGetSingularTriplet(svd,i,&sigma,PETSC_NULL,PETSC_NULL);
353: SVDComputeResidualNorms(svd,i,&norm1,&norm2);
354: *error = PetscSqrtReal(norm1*norm1+norm2*norm2);
355: if (sigma>*error) *error /= sigma;
356: return(0);
357: }
361: /*@
362: SVDGetOperationCounters - Gets the total number of matrix vector and dot
363: products used by the SVD object during the last SVDSolve() call.
365: Not Collective
367: Input Parameter:
368: . svd - SVD context
370: Output Parameter:
371: + matvecs - number of matrix vector product operations
372: - dots - number of dot product operations
374: Notes:
375: These counters are reset to zero at each successive call to SVDSolve().
377: Level: intermediate
379: @*/
380: PetscErrorCode SVDGetOperationCounters(SVD svd,PetscInt* matvecs,PetscInt* dots)
381: {
383:
386: if (matvecs) *matvecs = svd->matvecs;
387: if (dots) {
388: if (!svd->ip) { SVDGetIP(svd,&svd->ip); }
389: IPGetOperationCounters(svd->ip,dots);
390: }
391: return(0);
392: }