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: SVDMatCreateVecs(svd,NULL,&tl);
44: BVSetSizesFromVec(svd->U,tl,svd->ncv);
45: VecDestroy(&tl);
46: }
47: for (j=0;j<svd->nconv;j++) {
48: BVGetColumn(svd->V,j,&vj);
49: BVGetColumn(svd->U,j,&uj);
50: SVDMatMult(svd,PETSC_FALSE,vj,uj);
51: BVRestoreColumn(svd->V,j,&vj);
52: BVRestoreColumn(svd->U,j,&uj);
53: BVOrthogonalizeColumn(svd->U,j,NULL,&norm,NULL);
54: BVScaleColumn(svd->U,j,1.0/norm);
55: }
56: break;
57: default: 58: break;
59: }
60: svd->state = SVD_STATE_VECTORS;
61: return(0);
62: }
66: /*@
67: SVDSolve - Solves the singular value problem.
69: Collective on SVD 71: Input Parameter:
72: . svd - singular value solver context obtained from SVDCreate()
74: Options Database Keys:
75: + -svd_view - print information about the solver used
76: . -svd_view_mat binary - save the matrix to the default binary viewer
77: . -svd_view_vectors binary - save the computed singular vectors to the default binary viewer
78: . -svd_view_values - print computed singular values
79: . -svd_converged_reason - print reason for convergence, and number of iterations
80: . -svd_error_absolute - print absolute errors of each singular triplet
81: - -svd_error_relative - print relative errors of each singular triplet
83: Level: beginner
85: .seealso: SVDCreate(), SVDSetUp(), SVDDestroy()
86: @*/
87: PetscErrorCode SVDSolve(SVD svd) 88: {
90: PetscInt i,*workperm;
94: if (svd->state>=SVD_STATE_SOLVED) return(0);
95: PetscLogEventBegin(SVD_Solve,svd,0,0,0);
97: /* call setup */
98: SVDSetUp(svd);
99: svd->its = 0;
100: svd->nconv = 0;
101: for (i=0;i<svd->ncv;i++) {
102: svd->sigma[i] = 0.0;
103: svd->errest[i] = 0.0;
104: svd->perm[i] = i;
105: }
106: SVDViewFromOptions(svd,NULL,"-svd_view_pre");
108: (*svd->ops->solve)(svd);
109: svd->state = (svd->leftbasis)? SVD_STATE_VECTORS: SVD_STATE_SOLVED;
111: /* sort singular triplets */
112: if (svd->which == SVD_SMALLEST) {
113: PetscSortRealWithPermutation(svd->nconv,svd->sigma,svd->perm);
114: } else {
115: PetscMalloc(sizeof(PetscInt)*svd->nconv,&workperm);
116: for (i=0;i<svd->nconv;i++) workperm[i] = i;
117: PetscSortRealWithPermutation(svd->nconv,svd->sigma,workperm);
118: for (i=0;i<svd->nconv;i++) svd->perm[i] = workperm[svd->nconv-i-1];
119: PetscFree(workperm);
120: }
121: PetscLogEventEnd(SVD_Solve,svd,0,0,0);
123: /* various viewers */
124: SVDViewFromOptions(svd,NULL,"-svd_view");
125: SVDReasonViewFromOptions(svd);
126: SVDErrorViewFromOptions(svd);
127: SVDValuesViewFromOptions(svd);
128: SVDVectorsViewFromOptions(svd);
129: MatViewFromOptions(svd->OP,(PetscObject)svd,"-svd_view_mat");
131: /* Remove the initial subspaces */
132: svd->nini = 0;
133: svd->ninil = 0;
134: return(0);
135: }
139: /*@
140: SVDGetIterationNumber - Gets the current iteration number. If the
141: call to SVDSolve() is complete, then it returns the number of iterations
142: carried out by the solution method.
144: Not Collective
146: Input Parameter:
147: . svd - the singular value solver context
149: Output Parameter:
150: . its - number of iterations
152: Level: intermediate
154: Note:
155: During the i-th iteration this call returns i-1. If SVDSolve() is
156: complete, then parameter "its" contains either the iteration number at
157: which convergence was successfully reached, or failure was detected.
158: Call SVDGetConvergedReason() to determine if the solver converged or
159: failed and why.
161: .seealso: SVDGetConvergedReason(), SVDSetTolerances()
162: @*/
163: PetscErrorCode SVDGetIterationNumber(SVD svd,PetscInt *its)164: {
168: *its = svd->its;
169: return(0);
170: }
174: /*@
175: SVDGetConvergedReason - Gets the reason why the SVDSolve() iteration was
176: stopped.
178: Not Collective
180: Input Parameter:
181: . svd - the singular value solver context
183: Output Parameter:
184: . reason - negative value indicates diverged, positive value converged
185: (see SVDConvergedReason)
187: Notes:
189: Possible values for reason are
190: + SVD_CONVERGED_TOL - converged up to tolerance
191: . SVD_CONVERGED_USER - converged due to a user-defined condition
192: . SVD_DIVERGED_ITS - required more than max_it iterations to reach convergence
193: - SVD_DIVERGED_BREAKDOWN - generic breakdown in method
195: Can only be called after the call to SVDSolve() is complete.
197: Level: intermediate
199: .seealso: SVDSetTolerances(), SVDSolve(), SVDConvergedReason200: @*/
201: PetscErrorCode SVDGetConvergedReason(SVD svd,SVDConvergedReason *reason)202: {
206: SVDCheckSolved(svd,1);
207: *reason = svd->reason;
208: return(0);
209: }
213: /*@
214: SVDGetConverged - Gets the number of converged singular values.
216: Not Collective
218: Input Parameter:
219: . svd - the singular value solver context
221: Output Parameter:
222: . nconv - number of converged singular values
224: Note:
225: This function should be called after SVDSolve() has finished.
227: Level: beginner
229: @*/
230: PetscErrorCode SVDGetConverged(SVD svd,PetscInt *nconv)231: {
235: SVDCheckSolved(svd,1);
236: *nconv = svd->nconv;
237: return(0);
238: }
242: /*@
243: SVDGetSingularTriplet - Gets the i-th triplet of the singular value decomposition
244: as computed by SVDSolve(). The solution consists in the singular value and its left
245: and right singular vectors.
247: Not Collective, but vectors are shared by all processors that share the SVD249: Input Parameters:
250: + svd - singular value solver context
251: - i - index of the solution
253: Output Parameters:
254: + sigma - singular value
255: . u - left singular vector
256: - v - right singular vector
258: Note:
259: Both U or V can be NULL if singular vectors are not required.
260: Otherwise, the caller must provide valid Vec objects, i.e.,
261: they must be created by the calling program with e.g. MatCreateVecs().
263: The index i should be a value between 0 and nconv-1 (see SVDGetConverged()).
264: Singular triplets are indexed according to the ordering criterion established
265: with SVDSetWhichSingularTriplets().
267: Level: beginner
269: .seealso: SVDSolve(), SVDGetConverged(), SVDSetWhichSingularTriplets()
270: @*/
271: PetscErrorCode SVDGetSingularTriplet(SVD svd,PetscInt i,PetscReal *sigma,Vec u,Vec v)272: {
274: PetscInt M,N;
275: Vec w;
280: SVDCheckSolved(svd,1);
283: if (i<0 || i>=svd->nconv) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
284: *sigma = svd->sigma[svd->perm[i]];
285: MatGetSize(svd->OP,&M,&N);
286: if (M<N) { w = u; u = v; v = w; }
287: if (u) {
288: SVDComputeVectors(svd);
289: BVCopyVec(svd->U,svd->perm[i],u);
290: }
291: if (v) {
292: BVCopyVec(svd->V,svd->perm[i],v);
293: }
294: return(0);
295: }
299: /*
300: SVDComputeResidualNorms_Private - Computes the norms of the left and
301: right residuals associated with the i-th computed singular triplet.
302: @*/
303: static PetscErrorCode SVDComputeResidualNorms_Private(SVD svd,PetscInt i,PetscReal *norm1,PetscReal *norm2)304: {
306: Vec u,v,x = NULL,y = NULL;
307: PetscReal sigma;
308: PetscInt M,N;
311: MatCreateVecs(svd->OP,&v,&u);
312: SVDGetSingularTriplet(svd,i,&sigma,u,v);
313: /* norm1 = ||A*v-sigma*u||_2 */
314: if (norm1) {
315: VecDuplicate(u,&x);
316: MatMult(svd->OP,v,x);
317: VecAXPY(x,-sigma,u);
318: VecNorm(x,NORM_2,norm1);
319: }
320: /* norm2 = ||A^T*u-sigma*v||_2 */
321: if (norm2) {
322: VecDuplicate(v,&y);
323: if (svd->A && svd->AT) {
324: MatGetSize(svd->OP,&M,&N);
325: if (M<N) {
326: MatMult(svd->A,u,y);
327: } else {
328: MatMult(svd->AT,u,y);
329: }
330: } else {
331: #if defined(PETSC_USE_COMPLEX)
332: MatMultHermitianTranspose(svd->OP,u,y);
333: #else
334: MatMultTranspose(svd->OP,u,y);
335: #endif
336: }
337: VecAXPY(y,-sigma,v);
338: VecNorm(y,NORM_2,norm2);
339: }
341: VecDestroy(&v);
342: VecDestroy(&u);
343: VecDestroy(&x);
344: VecDestroy(&y);
345: return(0);
346: }
350: /*@
351: SVDComputeError - Computes the error (based on the residual norm) associated
352: with the i-th singular triplet.
354: Collective on SVD356: Input Parameter:
357: + svd - the singular value solver context
358: . i - the solution index
359: - type - the type of error to compute
361: Output Parameter:
362: . error - the error
364: Notes:
365: The error can be computed in various ways, all of them based on the residual
366: norm obtained as sqrt(n1^2+n2^2) with n1 = ||A*v-sigma*u||_2 and
367: n2 = ||A^T*u-sigma*v||_2, where sigma is the singular value, u is the left
368: singular vector and v is the right singular vector.
370: Level: beginner
372: .seealso: SVDErrorType, SVDSolve()
373: @*/
374: PetscErrorCode SVDComputeError(SVD svd,PetscInt i,SVDErrorType type,PetscReal *error)375: {
377: PetscReal sigma,norm1,norm2;
384: SVDCheckSolved(svd,1);
385: SVDGetSingularTriplet(svd,i,&sigma,NULL,NULL);
386: SVDComputeResidualNorms_Private(svd,i,&norm1,&norm2);
387: *error = PetscSqrtReal(norm1*norm1+norm2*norm2);
388: switch (type) {
389: case SVD_ERROR_ABSOLUTE:
390: break;
391: case SVD_ERROR_RELATIVE:
392: *error /= sigma;
393: break;
394: default:395: SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Invalid error type");
396: }
397: return(0);
398: }