Actual source code: svdsolve.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: SVD routines related to the solution process
12: */
14: #include <slepc/private/svdimpl.h>
15: #include <slepc/private/bvimpl.h>
17: /*
18: SVDComputeVectors_Left - Compute left singular vectors as U=A*V.
19: Only done if the leftbasis flag is false. Assumes V is available.
20: */
21: PetscErrorCode SVDComputeVectors_Left(SVD svd)
22: {
23: Vec tl,omega2,u,v,w;
24: PetscInt i,oldsize;
25: VecType vtype;
26: const PetscScalar* varray;
28: PetscFunctionBegin;
29: if (!svd->leftbasis) {
30: /* generate left singular vectors on U */
31: if (!svd->U) PetscCall(SVDGetBV(svd,NULL,&svd->U));
32: PetscCall(BVGetSizes(svd->U,NULL,NULL,&oldsize));
33: if (!oldsize) {
34: if (!((PetscObject)svd->U)->type_name) PetscCall(BVSetType(svd->U,((PetscObject)svd->V)->type_name));
35: PetscCall(MatCreateVecsEmpty(svd->A,NULL,&tl));
36: PetscCall(BVSetSizesFromVec(svd->U,tl,svd->ncv));
37: PetscCall(VecDestroy(&tl));
38: }
39: PetscCall(BVSetActiveColumns(svd->V,0,svd->nconv));
40: PetscCall(BVSetActiveColumns(svd->U,0,svd->nconv));
41: if (!svd->ishyperbolic) PetscCall(BVMatMult(svd->V,svd->A,svd->U));
42: else if (svd->swapped) { /* compute right singular vectors as V=A'*Omega*U */
43: PetscCall(MatCreateVecs(svd->A,&w,NULL));
44: for (i=0;i<svd->nconv;i++) {
45: PetscCall(BVGetColumn(svd->V,i,&v));
46: PetscCall(BVGetColumn(svd->U,i,&u));
47: PetscCall(VecPointwiseMult(w,v,svd->omega));
48: PetscCall(MatMult(svd->A,w,u));
49: PetscCall(BVRestoreColumn(svd->V,i,&v));
50: PetscCall(BVRestoreColumn(svd->U,i,&u));
51: }
52: PetscCall(VecDestroy(&w));
53: } else { /* compute left singular vectors as usual U=A*V, and set-up Omega-orthogonalization of U */
54: PetscCall(BV_SetMatrixDiagonal(svd->U,svd->omega,svd->A));
55: PetscCall(BVMatMult(svd->V,svd->A,svd->U));
56: }
57: PetscCall(BVOrthogonalize(svd->U,NULL));
58: if (svd->ishyperbolic && !svd->swapped) { /* store signature after Omega-orthogonalization */
59: PetscCall(MatGetVecType(svd->A,&vtype));
60: PetscCall(VecCreate(PETSC_COMM_SELF,&omega2));
61: PetscCall(VecSetSizes(omega2,svd->nconv,svd->nconv));
62: PetscCall(VecSetType(omega2,vtype));
63: PetscCall(BVGetSignature(svd->U,omega2));
64: PetscCall(VecGetArrayRead(omega2,&varray));
65: for (i=0;i<svd->nconv;i++) {
66: svd->sign[i] = PetscRealPart(varray[i]);
67: if (PetscRealPart(varray[i])<0.0) PetscCall(BVScaleColumn(svd->U,i,-1.0));
68: }
69: PetscCall(VecRestoreArrayRead(omega2,&varray));
70: PetscCall(VecDestroy(&omega2));
71: }
72: }
73: PetscFunctionReturn(PETSC_SUCCESS);
74: }
76: PetscErrorCode SVDComputeVectors(SVD svd)
77: {
78: PetscFunctionBegin;
79: SVDCheckSolved(svd,1);
80: if (svd->state==SVD_STATE_SOLVED) PetscTryTypeMethod(svd,computevectors);
81: svd->state = SVD_STATE_VECTORS;
82: PetscFunctionReturn(PETSC_SUCCESS);
83: }
85: /*@
86: SVDSolve - Solves the singular value problem.
88: Collective
90: Input Parameter:
91: . svd - the singular value solver context
93: Options Database Keys:
94: + -svd_view - print information about the solver used
95: . -svd_view_mat0 - view the first matrix (A)
96: . -svd_view_mat1 - view the second matrix (B)
97: . -svd_view_signature - view the signature matrix (omega)
98: . -svd_view_vectors - view the computed singular vectors
99: . -svd_view_values - view the computed singular values
100: . -svd_converged_reason - print reason for convergence, and number of iterations
101: . -svd_error_absolute - print absolute errors of each singular triplet
102: . -svd_error_relative - print relative errors of each singular triplet
103: - -svd_error_norm - print errors relative to the matrix norms of each singular triplet
105: Notes:
106: All the command-line options listed above admit an optional argument specifying
107: the viewer type and options. For instance, use '-svd_view_mat0 binary:amatrix.bin'
108: to save the A matrix to a binary file, '-svd_view_values draw' to draw the computed
109: singular values graphically, or '-svd_error_relative :myerr.m:ascii_matlab' to save
110: the errors in a file that can be executed in Matlab.
112: Level: beginner
114: .seealso: [](ch:svd), `SVDCreate()`, `SVDSetUp()`, `SVDDestroy()`
115: @*/
116: PetscErrorCode SVDSolve(SVD svd)
117: {
118: PetscInt i,m,n,*workperm;
120: PetscFunctionBegin;
122: if (svd->state>=SVD_STATE_SOLVED) PetscFunctionReturn(PETSC_SUCCESS);
123: PetscCall(PetscLogEventBegin(SVD_Solve,svd,0,0,0));
125: /* call setup */
126: PetscCall(SVDSetUp(svd));
128: /* safeguard for matrices with zero rows or columns */
129: PetscCall(MatGetSize(svd->OP,&m,&n));
130: if (m == 0 || n == 0) {
131: svd->nconv = 0;
132: svd->reason = SVD_CONVERGED_TOL;
133: svd->state = SVD_STATE_SOLVED;
134: PetscFunctionReturn(PETSC_SUCCESS);
135: }
137: svd->its = 0;
138: svd->nconv = 0;
139: for (i=0;i<svd->ncv;i++) {
140: svd->sigma[i] = 0.0;
141: svd->errest[i] = 0.0;
142: svd->perm[i] = i;
143: }
144: PetscCall(SVDViewFromOptions(svd,NULL,"-svd_view_pre"));
146: switch (svd->problem_type) {
147: case SVD_STANDARD:
148: PetscUseTypeMethod(svd,solve);
149: break;
150: case SVD_GENERALIZED:
151: PetscUseTypeMethod(svd,solveg);
152: break;
153: case SVD_HYPERBOLIC:
154: PetscUseTypeMethod(svd,solveh);
155: break;
156: }
157: svd->state = SVD_STATE_SOLVED;
159: /* sort singular triplets */
160: if (svd->which == SVD_SMALLEST) PetscCall(PetscSortRealWithPermutation(svd->nconv,svd->sigma,svd->perm));
161: else {
162: PetscCall(PetscMalloc1(svd->nconv,&workperm));
163: for (i=0;i<svd->nconv;i++) workperm[i] = i;
164: PetscCall(PetscSortRealWithPermutation(svd->nconv,svd->sigma,workperm));
165: for (i=0;i<svd->nconv;i++) svd->perm[i] = workperm[svd->nconv-i-1];
166: PetscCall(PetscFree(workperm));
167: }
168: PetscCall(PetscLogEventEnd(SVD_Solve,svd,0,0,0));
170: /* various viewers */
171: PetscCall(SVDViewFromOptions(svd,NULL,"-svd_view"));
172: PetscCall(SVDConvergedReasonViewFromOptions(svd));
173: PetscCall(SVDErrorViewFromOptions(svd));
174: PetscCall(SVDValuesViewFromOptions(svd));
175: PetscCall(SVDVectorsViewFromOptions(svd));
176: PetscCall(MatViewFromOptions(svd->OP,(PetscObject)svd,"-svd_view_mat0"));
177: if (svd->isgeneralized) PetscCall(MatViewFromOptions(svd->OPb,(PetscObject)svd,"-svd_view_mat1"));
178: if (svd->ishyperbolic) PetscCall(VecViewFromOptions(svd->omega,(PetscObject)svd,"-svd_view_signature"));
180: /* Remove the initial subspaces */
181: svd->nini = 0;
182: svd->ninil = 0;
183: PetscFunctionReturn(PETSC_SUCCESS);
184: }
186: /*@
187: SVDGetIterationNumber - Gets the current iteration number. If the
188: call to SVDSolve() is complete, then it returns the number of iterations
189: carried out by the solution method.
191: Not Collective
193: Input Parameter:
194: . svd - the singular value solver context
196: Output Parameter:
197: . its - number of iterations
199: Note:
200: During the i-th iteration this call returns i-1. If SVDSolve() is
201: complete, then parameter "its" contains either the iteration number at
202: which convergence was successfully reached, or failure was detected.
203: Call SVDGetConvergedReason() to determine if the solver converged or
204: failed and why.
206: Level: intermediate
208: .seealso: [](ch:svd), `SVDGetConvergedReason()`, `SVDSetTolerances()`
209: @*/
210: PetscErrorCode SVDGetIterationNumber(SVD svd,PetscInt *its)
211: {
212: PetscFunctionBegin;
214: PetscAssertPointer(its,2);
215: *its = svd->its;
216: PetscFunctionReturn(PETSC_SUCCESS);
217: }
219: /*@
220: SVDGetConvergedReason - Gets the reason why the `SVDSolve()` iteration was
221: stopped.
223: Not Collective
225: Input Parameter:
226: . svd - the singular value solver context
228: Output Parameter:
229: . reason - negative value indicates diverged, positive value converged, see
230: `SVDConvergedReason` for the possible values
232: Options Database Key:
233: . -svd_converged_reason - print reason for convergence/divergence, and number of iterations
235: Note:
236: If this routine is called before or doing the `SVDSolve()` the value of
237: `SVD_CONVERGED_ITERATING` is returned.
239: Level: intermediate
241: .seealso: [](ch:svd), `SVDSetTolerances()`, `SVDSolve()`, `SVDConvergedReason`
242: @*/
243: PetscErrorCode SVDGetConvergedReason(SVD svd,SVDConvergedReason *reason)
244: {
245: PetscFunctionBegin;
247: PetscAssertPointer(reason,2);
248: SVDCheckSolved(svd,1);
249: *reason = svd->reason;
250: PetscFunctionReturn(PETSC_SUCCESS);
251: }
253: /*@
254: SVDGetConverged - Gets the number of converged singular values.
256: Not Collective
258: Input Parameter:
259: . svd - the singular value solver context
261: Output Parameter:
262: . nconv - number of converged singular values
264: Note:
265: This function should be called after SVDSolve() has finished.
267: Level: beginner
269: .seealso: [](ch:svd), `SVDSetDimensions()`, `SVDSolve()`, `SVDGetSingularTriplet()`
270: @*/
271: PetscErrorCode SVDGetConverged(SVD svd,PetscInt *nconv)
272: {
273: PetscFunctionBegin;
275: PetscAssertPointer(nconv,2);
276: SVDCheckSolved(svd,1);
277: *nconv = svd->nconv;
278: PetscFunctionReturn(PETSC_SUCCESS);
279: }
281: /*@
282: SVDGetSingularTriplet - Gets the i-th triplet of the singular value decomposition
283: as computed by SVDSolve(). The solution consists in the singular value and its left
284: and right singular vectors.
286: Collective
288: Input Parameters:
289: + svd - the singular value solver context
290: - i - index of the solution
292: Output Parameters:
293: + sigma - singular value
294: . u - left singular vector
295: - v - right singular vector
297: Note:
298: Both u or v can be NULL if singular vectors are not required.
299: Otherwise, the caller must provide valid Vec objects, i.e.,
300: they must be created by the calling program with e.g. MatCreateVecs().
302: The index i should be a value between 0 and nconv-1 (see SVDGetConverged()).
303: Singular triplets are indexed according to the ordering criterion established
304: with SVDSetWhichSingularTriplets().
306: In the case of GSVD, the solution consists in three vectors u,v,x that are
307: returned as follows. Vector x is returned in the right singular vector
308: (argument v) and has length equal to the number of columns of A and B.
309: The other two vectors are returned stacked on top of each other [u;v] in
310: the left singular vector argument, with length equal to m+n (number of rows
311: of A plus number of rows of B).
313: Level: beginner
315: .seealso: [](ch:svd), `SVDSolve()`, `SVDGetConverged()`, `SVDSetWhichSingularTriplets()`
316: @*/
317: PetscErrorCode SVDGetSingularTriplet(SVD svd,PetscInt i,PetscReal *sigma,Vec u,Vec v)
318: {
319: PetscInt M,N;
320: Vec w;
322: PetscFunctionBegin;
325: SVDCheckSolved(svd,1);
328: PetscCheck(i>=0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"The index cannot be negative");
329: PetscCheck(i<svd->nconv,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"The index can be nconv-1 at most, see SVDGetConverged()");
330: if (sigma) *sigma = svd->sigma[svd->perm[i]];
331: if (u || v) {
332: if (!svd->isgeneralized) {
333: PetscCall(MatGetSize(svd->OP,&M,&N));
334: if (M<N) { w = u; u = v; v = w; }
335: }
336: PetscCall(SVDComputeVectors(svd));
337: if (u) PetscCall(BVCopyVec(svd->U,svd->perm[i],u));
338: if (v) PetscCall(BVCopyVec(svd->V,svd->perm[i],v));
339: }
340: PetscFunctionReturn(PETSC_SUCCESS);
341: }
343: /*
344: SVDComputeResidualNorms_Standard - Computes the norms of the left and
345: right residuals associated with the i-th computed singular triplet.
347: Input Parameters:
348: sigma - singular value
349: u,v - singular vectors
350: x,y - two work vectors with the same dimensions as u,v
351: */
352: static PetscErrorCode SVDComputeResidualNorms_Standard(SVD svd,PetscReal sigma,Vec u,Vec v,Vec x,Vec y,PetscReal *norm1,PetscReal *norm2)
353: {
354: PetscInt M,N;
356: PetscFunctionBegin;
357: /* norm1 = ||A*v-sigma*u||_2 */
358: if (norm1) {
359: PetscCall(MatMult(svd->OP,v,x));
360: PetscCall(VecAXPY(x,-sigma,u));
361: PetscCall(VecNorm(x,NORM_2,norm1));
362: }
363: /* norm2 = ||A^T*u-sigma*v||_2 */
364: if (norm2) {
365: PetscCall(MatGetSize(svd->OP,&M,&N));
366: if (M<N) PetscCall(MatMult(svd->A,u,y));
367: else PetscCall(MatMult(svd->AT,u,y));
368: PetscCall(VecAXPY(y,-sigma,v));
369: PetscCall(VecNorm(y,NORM_2,norm2));
370: }
371: PetscFunctionReturn(PETSC_SUCCESS);
372: }
374: /*
375: SVDComputeResidualNorms_Generalized - In GSVD, compute the residual norms
376: norm1 = ||s^2*A'*u-c*B'*B*x||_2 and norm2 = ||c^2*B'*v-s*A'*A*x||_2.
378: Input Parameters:
379: sigma - singular value
380: uv - left singular vectors [u;v]
381: x - right singular vector
382: y,z - two work vectors with the same dimension as x
383: */
384: static PetscErrorCode SVDComputeResidualNorms_Generalized(SVD svd,PetscReal sigma,Vec uv,Vec x,Vec y,Vec z,PetscReal *norm1,PetscReal *norm2)
385: {
386: Vec u,v,ax,bx,nest,aux[2];
387: PetscReal c,s;
389: PetscFunctionBegin;
390: PetscCall(MatCreateVecs(svd->OP,NULL,&u));
391: PetscCall(MatCreateVecs(svd->OPb,NULL,&v));
392: aux[0] = u;
393: aux[1] = v;
394: PetscCall(VecCreateNest(PetscObjectComm((PetscObject)svd),2,NULL,aux,&nest));
395: PetscCall(VecCopy(uv,nest));
397: s = 1.0/PetscSqrtReal(1.0+sigma*sigma);
398: c = sigma*s;
400: /* norm1 = ||s^2*A'*u-c*B'*B*x||_2 */
401: if (norm1) {
402: PetscCall(VecDuplicate(v,&bx));
403: PetscCall(MatMultHermitianTranspose(svd->OP,u,z));
404: PetscCall(MatMult(svd->OPb,x,bx));
405: PetscCall(MatMultHermitianTranspose(svd->OPb,bx,y));
406: PetscCall(VecAXPBY(y,s*s,-c,z));
407: PetscCall(VecNorm(y,NORM_2,norm1));
408: PetscCall(VecDestroy(&bx));
409: }
410: /* norm2 = ||c^2*B'*v-s*A'*A*x||_2 */
411: if (norm2) {
412: PetscCall(VecDuplicate(u,&ax));
413: PetscCall(MatMultHermitianTranspose(svd->OPb,v,z));
414: PetscCall(MatMult(svd->OP,x,ax));
415: PetscCall(MatMultHermitianTranspose(svd->OP,ax,y));
416: PetscCall(VecAXPBY(y,c*c,-s,z));
417: PetscCall(VecNorm(y,NORM_2,norm2));
418: PetscCall(VecDestroy(&ax));
419: }
421: PetscCall(VecDestroy(&v));
422: PetscCall(VecDestroy(&u));
423: PetscCall(VecDestroy(&nest));
424: PetscFunctionReturn(PETSC_SUCCESS);
425: }
427: /*
428: SVDComputeResidualNorms_Hyperbolic - Computes the norms of the left and
429: right residuals associated with the i-th computed singular triplet.
431: Input Parameters:
432: sigma - singular value
433: sign - corresponding element of the signature Omega2
434: u,v - singular vectors
435: x,y,z - three work vectors with the same dimensions as u,v,u
436: */
437: static PetscErrorCode SVDComputeResidualNorms_Hyperbolic(SVD svd,PetscReal sigma,PetscReal sign,Vec u,Vec v,Vec x,Vec y,Vec z,PetscReal *norm1,PetscReal *norm2)
438: {
439: PetscInt M,N;
441: PetscFunctionBegin;
442: /* norm1 = ||A*v-sigma*u||_2 */
443: if (norm1) {
444: PetscCall(MatMult(svd->OP,v,x));
445: PetscCall(VecAXPY(x,-sigma,u));
446: PetscCall(VecNorm(x,NORM_2,norm1));
447: }
448: /* norm2 = ||A^T*Omega*u-sigma*sign*v||_2 */
449: if (norm2) {
450: PetscCall(MatGetSize(svd->OP,&M,&N));
451: PetscCall(VecPointwiseMult(z,u,svd->omega));
452: if (M<N) PetscCall(MatMult(svd->A,z,y));
453: else PetscCall(MatMult(svd->AT,z,y));
454: PetscCall(VecAXPY(y,-sigma*sign,v));
455: PetscCall(VecNorm(y,NORM_2,norm2));
456: }
457: PetscFunctionReturn(PETSC_SUCCESS);
458: }
460: /*@
461: SVDComputeError - Computes the error (based on the residual norm) associated
462: with the i-th singular triplet.
464: Collective
466: Input Parameters:
467: + svd - the singular value solver context
468: . i - the solution index
469: - type - the type of error to compute
471: Output Parameter:
472: . error - the error
474: Notes:
475: The error can be computed in various ways, all of them based on the residual
476: norm obtained as sqrt(n1^2+n2^2) with n1 = ||A*v-sigma*u||_2 and
477: n2 = ||A^T*u-sigma*v||_2, where sigma is the singular value, u is the left
478: singular vector and v is the right singular vector.
480: In the case of the GSVD, the two components of the residual norm are
481: n1 = ||s^2*A'*u-c*B'*B*x||_2 and n2 = ||c^2*B'*v-s*A'*A*x||_2, where [u;v]
482: are the left singular vectors and x is the right singular vector, with
483: sigma=c/s.
485: Level: beginner
487: .seealso: [](ch:svd), `SVDErrorType`, `SVDSolve()`
488: @*/
489: PetscErrorCode SVDComputeError(SVD svd,PetscInt i,SVDErrorType type,PetscReal *error)
490: {
491: PetscReal sigma,norm1,norm2,c,s;
492: Vec u=NULL,v=NULL,x=NULL,y=NULL,z=NULL;
493: PetscReal vecnorm=1.0;
495: PetscFunctionBegin;
499: PetscAssertPointer(error,4);
500: SVDCheckSolved(svd,1);
502: /* allocate work vectors */
503: switch (svd->problem_type) {
504: case SVD_STANDARD:
505: PetscCall(SVDSetWorkVecs(svd,2,2));
506: u = svd->workl[0];
507: v = svd->workr[0];
508: x = svd->workl[1];
509: y = svd->workr[1];
510: break;
511: case SVD_GENERALIZED:
512: PetscCall(SVDSetWorkVecs(svd,1,3));
513: u = svd->workl[0];
514: v = svd->workr[0];
515: x = svd->workr[1];
516: y = svd->workr[2];
517: break;
518: case SVD_HYPERBOLIC:
519: PetscCall(SVDSetWorkVecs(svd,3,2));
520: u = svd->workl[0];
521: v = svd->workr[0];
522: x = svd->workl[1];
523: y = svd->workr[1];
524: z = svd->workl[2];
525: break;
526: }
528: /* compute residual norm */
529: PetscCall(SVDGetSingularTriplet(svd,i,&sigma,u,v));
530: switch (svd->problem_type) {
531: case SVD_STANDARD:
532: PetscCall(SVDComputeResidualNorms_Standard(svd,sigma,u,v,x,y,&norm1,&norm2));
533: break;
534: case SVD_GENERALIZED:
535: PetscCall(SVDComputeResidualNorms_Generalized(svd,sigma,u,v,x,y,&norm1,&norm2));
536: break;
537: case SVD_HYPERBOLIC:
538: PetscCall(SVDComputeResidualNorms_Hyperbolic(svd,sigma,svd->sign[svd->perm[i]],u,v,x,y,z,&norm1,&norm2));
539: break;
540: }
541: *error = SlepcAbs(norm1,norm2);
543: /* compute 2-norm of eigenvector of the cyclic form */
544: if (type!=SVD_ERROR_ABSOLUTE) {
545: switch (svd->problem_type) {
546: case SVD_STANDARD:
547: vecnorm = PETSC_SQRT2;
548: break;
549: case SVD_GENERALIZED:
550: PetscCall(VecNorm(v,NORM_2,&vecnorm));
551: vecnorm = PetscSqrtReal(1.0+vecnorm*vecnorm);
552: break;
553: case SVD_HYPERBOLIC:
554: PetscCall(VecNorm(u,NORM_2,&vecnorm));
555: vecnorm = PetscSqrtReal(1.0+vecnorm*vecnorm);
556: break;
557: }
558: }
560: /* compute error */
561: switch (type) {
562: case SVD_ERROR_ABSOLUTE:
563: break;
564: case SVD_ERROR_RELATIVE:
565: if (svd->isgeneralized) {
566: s = 1.0/PetscSqrtReal(1.0+sigma*sigma);
567: c = sigma*s;
568: norm1 /= c*vecnorm;
569: norm2 /= s*vecnorm;
570: *error = PetscMax(norm1,norm2);
571: } else *error /= sigma*vecnorm;
572: break;
573: case SVD_ERROR_NORM:
574: if (!svd->nrma) PetscCall(MatNorm(svd->OP,NORM_INFINITY,&svd->nrma));
575: if (svd->isgeneralized && !svd->nrmb) PetscCall(MatNorm(svd->OPb,NORM_INFINITY,&svd->nrmb));
576: *error /= PetscMax(svd->nrma,svd->nrmb)*vecnorm;
577: break;
578: default:
579: SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Invalid error type");
580: }
581: PetscFunctionReturn(PETSC_SUCCESS);
582: }