1: /*
2: SVD routines related to the solution process.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
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: PetscErrorCode SVDComputeVectors(SVD svd) 29: {
31: Vec tl,uj,vj;
32: PetscInt j,oldsize;
33: PetscReal norm;
36: SVDCheckSolved(svd,1);
37: switch (svd->state) {
38: case SVD_STATE_SOLVED:
39: /* generate left singular vectors on U */
40: if (!svd->U) { SVDGetBV(svd,NULL,&svd->U); }
41: BVGetSizes(svd->U,NULL,NULL,&oldsize);
42: if (!oldsize) {
43: if (!((PetscObject)(svd->U))->type_name) {
44: BVSetType(svd->U,BVSVEC);
45: }
46: SVDMatCreateVecs(svd,NULL,&tl);
47: BVSetSizesFromVec(svd->U,tl,svd->ncv);
48: VecDestroy(&tl);
49: }
50: for (j=0;j<svd->nconv;j++) {
51: BVGetColumn(svd->V,j,&vj);
52: BVGetColumn(svd->U,j,&uj);
53: SVDMatMult(svd,PETSC_FALSE,vj,uj);
54: BVRestoreColumn(svd->V,j,&vj);
55: BVRestoreColumn(svd->U,j,&uj);
56: BVOrthogonalizeColumn(svd->U,j,NULL,&norm,NULL);
57: BVScaleColumn(svd->U,j,1.0/norm);
58: }
59: break;
60: default: 61: break;
62: }
63: svd->state = SVD_STATE_VECTORS;
64: return(0);
65: }
69: /*@
70: SVDSolve - Solves the singular value problem.
72: Collective on SVD 74: Input Parameter:
75: . svd - singular value solver context obtained from SVDCreate()
77: Options Database Keys:
78: + -svd_view - print information about the solver used
79: . -svd_view_mat binary - save the matrix to the default binary viewer
80: . -svd_view_vectors binary - save the computed singular vectors to the default binary viewer
81: . -svd_view_values - print computed singular values
82: . -svd_converged_reason - print reason for convergence, and number of iterations
83: . -svd_error_absolute - print absolute errors of each singular triplet
84: - -svd_error_relative - print relative errors of each singular triplet
86: Level: beginner
88: .seealso: SVDCreate(), SVDSetUp(), SVDDestroy()
89: @*/
90: PetscErrorCode SVDSolve(SVD svd) 91: {
93: PetscInt i,*workperm;
97: if (svd->state>=SVD_STATE_SOLVED) return(0);
98: PetscLogEventBegin(SVD_Solve,svd,0,0,0);
100: /* call setup */
101: SVDSetUp(svd);
102: svd->its = 0;
103: svd->nconv = 0;
104: for (i=0;i<svd->ncv;i++) {
105: svd->sigma[i] = 0.0;
106: svd->errest[i] = 0.0;
107: svd->perm[i] = i;
108: }
109: SVDViewFromOptions(svd,NULL,"-svd_view_pre");
111: (*svd->ops->solve)(svd);
112: svd->state = (svd->leftbasis)? SVD_STATE_VECTORS: SVD_STATE_SOLVED;
114: /* sort singular triplets */
115: if (svd->which == SVD_SMALLEST) {
116: PetscSortRealWithPermutation(svd->nconv,svd->sigma,svd->perm);
117: } else {
118: PetscMalloc(sizeof(PetscInt)*svd->nconv,&workperm);
119: for (i=0;i<svd->nconv;i++) workperm[i] = i;
120: PetscSortRealWithPermutation(svd->nconv,svd->sigma,workperm);
121: for (i=0;i<svd->nconv;i++) svd->perm[i] = workperm[svd->nconv-i-1];
122: PetscFree(workperm);
123: }
124: PetscLogEventEnd(SVD_Solve,svd,0,0,0);
126: /* various viewers */
127: SVDViewFromOptions(svd,NULL,"-svd_view");
128: SVDReasonViewFromOptions(svd);
129: SVDErrorViewFromOptions(svd);
130: SVDValuesViewFromOptions(svd);
131: SVDVectorsViewFromOptions(svd);
132: MatViewFromOptions(svd->OP,(PetscObject)svd,"-svd_view_mat");
134: /* Remove the initial subspaces */
135: svd->nini = 0;
136: svd->ninil = 0;
137: return(0);
138: }
142: /*@
143: SVDGetIterationNumber - Gets the current iteration number. If the
144: call to SVDSolve() is complete, then it returns the number of iterations
145: carried out by the solution method.
147: Not Collective
149: Input Parameter:
150: . svd - the singular value solver context
152: Output Parameter:
153: . its - number of iterations
155: Level: intermediate
157: Note:
158: During the i-th iteration this call returns i-1. If SVDSolve() is
159: complete, then parameter "its" contains either the iteration number at
160: which convergence was successfully reached, or failure was detected.
161: Call SVDGetConvergedReason() to determine if the solver converged or
162: failed and why.
164: .seealso: SVDGetConvergedReason(), SVDSetTolerances()
165: @*/
166: PetscErrorCode SVDGetIterationNumber(SVD svd,PetscInt *its)167: {
171: *its = svd->its;
172: return(0);
173: }
177: /*@
178: SVDGetConvergedReason - Gets the reason why the SVDSolve() iteration was
179: stopped.
181: Not Collective
183: Input Parameter:
184: . svd - the singular value solver context
186: Output Parameter:
187: . reason - negative value indicates diverged, positive value converged
188: (see SVDConvergedReason)
190: Notes:
192: Possible values for reason are
193: + SVD_CONVERGED_TOL - converged up to tolerance
194: . SVD_CONVERGED_USER - converged due to a user-defined condition
195: . SVD_DIVERGED_ITS - required more than max_it iterations to reach convergence
196: - SVD_DIVERGED_BREAKDOWN - generic breakdown in method
198: Can only be called after the call to SVDSolve() is complete.
200: Level: intermediate
202: .seealso: SVDSetTolerances(), SVDSolve(), SVDConvergedReason203: @*/
204: PetscErrorCode SVDGetConvergedReason(SVD svd,SVDConvergedReason *reason)205: {
209: SVDCheckSolved(svd,1);
210: *reason = svd->reason;
211: return(0);
212: }
216: /*@
217: SVDGetConverged - Gets the number of converged singular values.
219: Not Collective
221: Input Parameter:
222: . svd - the singular value solver context
224: Output Parameter:
225: . nconv - number of converged singular values
227: Note:
228: This function should be called after SVDSolve() has finished.
230: Level: beginner
232: @*/
233: PetscErrorCode SVDGetConverged(SVD svd,PetscInt *nconv)234: {
238: SVDCheckSolved(svd,1);
239: *nconv = svd->nconv;
240: return(0);
241: }
245: /*@
246: SVDGetSingularTriplet - Gets the i-th triplet of the singular value decomposition
247: as computed by SVDSolve(). The solution consists in the singular value and its left
248: and right singular vectors.
250: Not Collective, but vectors are shared by all processors that share the SVD252: Input Parameters:
253: + svd - singular value solver context
254: - i - index of the solution
256: Output Parameters:
257: + sigma - singular value
258: . u - left singular vector
259: - v - right singular vector
261: Note:
262: Both U or V can be NULL if singular vectors are not required.
263: Otherwise, the caller must provide valid Vec objects, i.e.,
264: they must be created by the calling program with e.g. MatCreateVecs().
266: The index i should be a value between 0 and nconv-1 (see SVDGetConverged()).
267: Singular triplets are indexed according to the ordering criterion established
268: with SVDSetWhichSingularTriplets().
270: Level: beginner
272: .seealso: SVDSolve(), SVDGetConverged(), SVDSetWhichSingularTriplets()
273: @*/
274: PetscErrorCode SVDGetSingularTriplet(SVD svd,PetscInt i,PetscReal *sigma,Vec u,Vec v)275: {
277: PetscInt M,N;
278: Vec w;
283: SVDCheckSolved(svd,1);
286: if (i<0 || i>=svd->nconv) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
287: *sigma = svd->sigma[svd->perm[i]];
288: MatGetSize(svd->OP,&M,&N);
289: if (M<N) { w = u; u = v; v = w; }
290: if (u) {
291: SVDComputeVectors(svd);
292: BVCopyVec(svd->U,svd->perm[i],u);
293: }
294: if (v) {
295: BVCopyVec(svd->V,svd->perm[i],v);
296: }
297: return(0);
298: }
302: /*
303: SVDComputeResidualNorms_Private - Computes the norms of the left and
304: right residuals associated with the i-th computed singular triplet.
305: @*/
306: static PetscErrorCode SVDComputeResidualNorms_Private(SVD svd,PetscInt i,PetscReal *norm1,PetscReal *norm2)307: {
309: Vec u,v,x = NULL,y = NULL;
310: PetscReal sigma;
311: PetscInt M,N;
314: MatCreateVecs(svd->OP,&v,&u);
315: SVDGetSingularTriplet(svd,i,&sigma,u,v);
316: /* norm1 = ||A*v-sigma*u||_2 */
317: if (norm1) {
318: VecDuplicate(u,&x);
319: MatMult(svd->OP,v,x);
320: VecAXPY(x,-sigma,u);
321: VecNorm(x,NORM_2,norm1);
322: }
323: /* norm2 = ||A^T*u-sigma*v||_2 */
324: if (norm2) {
325: VecDuplicate(v,&y);
326: if (svd->A && svd->AT) {
327: MatGetSize(svd->OP,&M,&N);
328: if (M<N) {
329: MatMult(svd->A,u,y);
330: } else {
331: MatMult(svd->AT,u,y);
332: }
333: } else {
334: #if defined(PETSC_USE_COMPLEX)
335: MatMultHermitianTranspose(svd->OP,u,y);
336: #else
337: MatMultTranspose(svd->OP,u,y);
338: #endif
339: }
340: VecAXPY(y,-sigma,v);
341: VecNorm(y,NORM_2,norm2);
342: }
344: VecDestroy(&v);
345: VecDestroy(&u);
346: VecDestroy(&x);
347: VecDestroy(&y);
348: return(0);
349: }
353: /*@
354: SVDComputeError - Computes the error (based on the residual norm) associated
355: with the i-th singular triplet.
357: Collective on SVD359: Input Parameter:
360: + svd - the singular value solver context
361: . i - the solution index
362: - type - the type of error to compute
364: Output Parameter:
365: . error - the error
367: Notes:
368: The error can be computed in various ways, all of them based on the residual
369: norm obtained as sqrt(n1^2+n2^2) with n1 = ||A*v-sigma*u||_2 and
370: n2 = ||A^T*u-sigma*v||_2, where sigma is the singular value, u is the left
371: singular vector and v is the right singular vector.
373: Level: beginner
375: .seealso: SVDErrorType, SVDSolve()
376: @*/
377: PetscErrorCode SVDComputeError(SVD svd,PetscInt i,SVDErrorType type,PetscReal *error)378: {
380: PetscReal sigma,norm1,norm2;
387: SVDCheckSolved(svd,1);
388: SVDGetSingularTriplet(svd,i,&sigma,NULL,NULL);
389: SVDComputeResidualNorms_Private(svd,i,&norm1,&norm2);
390: *error = PetscSqrtReal(norm1*norm1+norm2*norm2);
391: switch (type) {
392: case SVD_ERROR_ABSOLUTE:
393: break;
394: case SVD_ERROR_RELATIVE:
395: *error /= sigma;
396: break;
397: default:398: SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Invalid error type");
399: }
400: return(0);
401: }