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: Possible values for reason:
188: + SVD_CONVERGED_TOL - converged up to tolerance
189: . SVD_CONVERGED_USER - converged due to a user-defined condition
190: . SVD_DIVERGED_ITS - required more than max_it iterations to reach convergence
191: - SVD_DIVERGED_BREAKDOWN - generic breakdown in method
193: Note:
194: Can only be called after the call to SVDSolve() is complete.
196: Level: intermediate
198: .seealso: SVDSetTolerances(), SVDSolve(), SVDConvergedReason199: @*/
200: PetscErrorCode SVDGetConvergedReason(SVD svd,SVDConvergedReason *reason)201: {
205: SVDCheckSolved(svd,1);
206: *reason = svd->reason;
207: return(0);
208: }
212: /*@
213: SVDGetConverged - Gets the number of converged singular values.
215: Not Collective
217: Input Parameter:
218: . svd - the singular value solver context
220: Output Parameter:
221: . nconv - number of converged singular values
223: Note:
224: This function should be called after SVDSolve() has finished.
226: Level: beginner
228: @*/
229: PetscErrorCode SVDGetConverged(SVD svd,PetscInt *nconv)230: {
234: SVDCheckSolved(svd,1);
235: *nconv = svd->nconv;
236: return(0);
237: }
241: /*@
242: SVDGetSingularTriplet - Gets the i-th triplet of the singular value decomposition
243: as computed by SVDSolve(). The solution consists in the singular value and its left
244: and right singular vectors.
246: Not Collective, but vectors are shared by all processors that share the SVD248: Input Parameters:
249: + svd - singular value solver context
250: - i - index of the solution
252: Output Parameters:
253: + sigma - singular value
254: . u - left singular vector
255: - v - right singular vector
257: Note:
258: Both U or V can be NULL if singular vectors are not required.
259: Otherwise, the caller must provide valid Vec objects, i.e.,
260: they must be created by the calling program with e.g. MatCreateVecs().
262: The index i should be a value between 0 and nconv-1 (see SVDGetConverged()).
263: Singular triplets are indexed according to the ordering criterion established
264: with SVDSetWhichSingularTriplets().
266: Level: beginner
268: .seealso: SVDSolve(), SVDGetConverged(), SVDSetWhichSingularTriplets()
269: @*/
270: PetscErrorCode SVDGetSingularTriplet(SVD svd,PetscInt i,PetscReal *sigma,Vec u,Vec v)271: {
273: PetscInt M,N;
274: Vec w;
279: SVDCheckSolved(svd,1);
282: if (i<0 || i>=svd->nconv) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
283: *sigma = svd->sigma[svd->perm[i]];
284: MatGetSize(svd->OP,&M,&N);
285: if (M<N) { w = u; u = v; v = w; }
286: if (u) {
287: SVDComputeVectors(svd);
288: BVCopyVec(svd->U,svd->perm[i],u);
289: }
290: if (v) {
291: BVCopyVec(svd->V,svd->perm[i],v);
292: }
293: return(0);
294: }
298: /*
299: SVDComputeResidualNorms_Private - Computes the norms of the left and
300: right residuals associated with the i-th computed singular triplet.
301: @*/
302: static PetscErrorCode SVDComputeResidualNorms_Private(SVD svd,PetscInt i,PetscReal *norm1,PetscReal *norm2)303: {
305: Vec u,v,x = NULL,y = NULL;
306: PetscReal sigma;
307: PetscInt M,N;
310: MatCreateVecs(svd->OP,&v,&u);
311: SVDGetSingularTriplet(svd,i,&sigma,u,v);
312: /* norm1 = ||A*v-sigma*u||_2 */
313: if (norm1) {
314: VecDuplicate(u,&x);
315: MatMult(svd->OP,v,x);
316: VecAXPY(x,-sigma,u);
317: VecNorm(x,NORM_2,norm1);
318: }
319: /* norm2 = ||A^T*u-sigma*v||_2 */
320: if (norm2) {
321: VecDuplicate(v,&y);
322: if (svd->A && svd->AT) {
323: MatGetSize(svd->OP,&M,&N);
324: if (M<N) {
325: MatMult(svd->A,u,y);
326: } else {
327: MatMult(svd->AT,u,y);
328: }
329: } else {
330: #if defined(PETSC_USE_COMPLEX)
331: MatMultHermitianTranspose(svd->OP,u,y);
332: #else
333: MatMultTranspose(svd->OP,u,y);
334: #endif
335: }
336: VecAXPY(y,-sigma,v);
337: VecNorm(y,NORM_2,norm2);
338: }
340: VecDestroy(&v);
341: VecDestroy(&u);
342: VecDestroy(&x);
343: VecDestroy(&y);
344: return(0);
345: }
349: /*@
350: SVDComputeError - Computes the error (based on the residual norm) associated
351: with the i-th singular triplet.
353: Collective on SVD355: Input Parameter:
356: + svd - the singular value solver context
357: . i - the solution index
358: - type - the type of error to compute
360: Output Parameter:
361: . error - the error
363: Notes:
364: The error can be computed in various ways, all of them based on the residual
365: norm obtained as sqrt(n1^2+n2^2) with n1 = ||A*v-sigma*u||_2 and
366: n2 = ||A^T*u-sigma*v||_2, where sigma is the singular value, u is the left
367: singular vector and v is the right singular vector.
369: Level: beginner
371: .seealso: SVDErrorType, SVDSolve()
372: @*/
373: PetscErrorCode SVDComputeError(SVD svd,PetscInt i,SVDErrorType type,PetscReal *error)374: {
376: PetscReal sigma,norm1,norm2;
383: SVDCheckSolved(svd,1);
384: SVDGetSingularTriplet(svd,i,&sigma,NULL,NULL);
385: SVDComputeResidualNorms_Private(svd,i,&norm1,&norm2);
386: *error = PetscSqrtReal(norm1*norm1+norm2*norm2);
387: switch (type) {
388: case SVD_ERROR_ABSOLUTE:
389: break;
390: case SVD_ERROR_RELATIVE:
391: *error /= sigma;
392: break;
393: default:394: SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Invalid error type");
395: }
396: return(0);
397: }