Actual source code: svdsolve.c
slepc-3.19.0 2023-03-31
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 - singular value solver context obtained from SVDCreate()
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: SVDCreate(), SVDSetUp(), SVDDestroy()
115: @*/
116: PetscErrorCode SVDSolve(SVD svd)
117: {
118: PetscInt i,*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));
127: svd->its = 0;
128: svd->nconv = 0;
129: for (i=0;i<svd->ncv;i++) {
130: svd->sigma[i] = 0.0;
131: svd->errest[i] = 0.0;
132: svd->perm[i] = i;
133: }
134: PetscCall(SVDViewFromOptions(svd,NULL,"-svd_view_pre"));
136: switch (svd->problem_type) {
137: case SVD_STANDARD:
138: PetscUseTypeMethod(svd,solve);
139: break;
140: case SVD_GENERALIZED:
141: PetscUseTypeMethod(svd,solveg);
142: break;
143: case SVD_HYPERBOLIC:
144: PetscUseTypeMethod(svd,solveh);
145: break;
146: }
147: svd->state = SVD_STATE_SOLVED;
149: /* sort singular triplets */
150: if (svd->which == SVD_SMALLEST) PetscCall(PetscSortRealWithPermutation(svd->nconv,svd->sigma,svd->perm));
151: else {
152: PetscCall(PetscMalloc1(svd->nconv,&workperm));
153: for (i=0;i<svd->nconv;i++) workperm[i] = i;
154: PetscCall(PetscSortRealWithPermutation(svd->nconv,svd->sigma,workperm));
155: for (i=0;i<svd->nconv;i++) svd->perm[i] = workperm[svd->nconv-i-1];
156: PetscCall(PetscFree(workperm));
157: }
158: PetscCall(PetscLogEventEnd(SVD_Solve,svd,0,0,0));
160: /* various viewers */
161: PetscCall(SVDViewFromOptions(svd,NULL,"-svd_view"));
162: PetscCall(SVDConvergedReasonViewFromOptions(svd));
163: PetscCall(SVDErrorViewFromOptions(svd));
164: PetscCall(SVDValuesViewFromOptions(svd));
165: PetscCall(SVDVectorsViewFromOptions(svd));
166: PetscCall(MatViewFromOptions(svd->OP,(PetscObject)svd,"-svd_view_mat0"));
167: if (svd->isgeneralized) PetscCall(MatViewFromOptions(svd->OPb,(PetscObject)svd,"-svd_view_mat1"));
168: if (svd->ishyperbolic) PetscCall(VecViewFromOptions(svd->omega,(PetscObject)svd,"-svd_view_signature"));
170: /* Remove the initial subspaces */
171: svd->nini = 0;
172: svd->ninil = 0;
173: PetscFunctionReturn(PETSC_SUCCESS);
174: }
176: /*@
177: SVDGetIterationNumber - Gets the current iteration number. If the
178: call to SVDSolve() is complete, then it returns the number of iterations
179: carried out by the solution method.
181: Not Collective
183: Input Parameter:
184: . svd - the singular value solver context
186: Output Parameter:
187: . its - number of iterations
189: Note:
190: During the i-th iteration this call returns i-1. If SVDSolve() is
191: complete, then parameter "its" contains either the iteration number at
192: which convergence was successfully reached, or failure was detected.
193: Call SVDGetConvergedReason() to determine if the solver converged or
194: failed and why.
196: Level: intermediate
198: .seealso: SVDGetConvergedReason(), SVDSetTolerances()
199: @*/
200: PetscErrorCode SVDGetIterationNumber(SVD svd,PetscInt *its)
201: {
202: PetscFunctionBegin;
205: *its = svd->its;
206: PetscFunctionReturn(PETSC_SUCCESS);
207: }
209: /*@
210: SVDGetConvergedReason - Gets the reason why the SVDSolve() iteration was
211: stopped.
213: Not Collective
215: Input Parameter:
216: . svd - the singular value solver context
218: Output Parameter:
219: . reason - negative value indicates diverged, positive value converged
220: (see SVDConvergedReason)
222: Options Database Key:
223: . -svd_converged_reason - print the reason to a viewer
225: Notes:
226: Possible values for reason are
227: + SVD_CONVERGED_TOL - converged up to tolerance
228: . SVD_CONVERGED_USER - converged due to a user-defined condition
229: . SVD_CONVERGED_MAXIT - reached the maximum number of iterations with SVD_CONV_MAXIT criterion
230: . SVD_DIVERGED_ITS - required more than max_it iterations to reach convergence
231: - SVD_DIVERGED_BREAKDOWN - generic breakdown in method
233: Can only be called after the call to SVDSolve() is complete.
235: Level: intermediate
237: .seealso: SVDSetTolerances(), SVDSolve(), SVDConvergedReason
238: @*/
239: PetscErrorCode SVDGetConvergedReason(SVD svd,SVDConvergedReason *reason)
240: {
241: PetscFunctionBegin;
244: SVDCheckSolved(svd,1);
245: *reason = svd->reason;
246: PetscFunctionReturn(PETSC_SUCCESS);
247: }
249: /*@
250: SVDGetConverged - Gets the number of converged singular values.
252: Not Collective
254: Input Parameter:
255: . svd - the singular value solver context
257: Output Parameter:
258: . nconv - number of converged singular values
260: Note:
261: This function should be called after SVDSolve() has finished.
263: Level: beginner
265: .seealso: SVDSetDimensions(), SVDSolve(), SVDGetSingularTriplet()
266: @*/
267: PetscErrorCode SVDGetConverged(SVD svd,PetscInt *nconv)
268: {
269: PetscFunctionBegin;
272: SVDCheckSolved(svd,1);
273: *nconv = svd->nconv;
274: PetscFunctionReturn(PETSC_SUCCESS);
275: }
277: /*@C
278: SVDGetSingularTriplet - Gets the i-th triplet of the singular value decomposition
279: as computed by SVDSolve(). The solution consists in the singular value and its left
280: and right singular vectors.
282: Collective
284: Input Parameters:
285: + svd - singular value solver context
286: - i - index of the solution
288: Output Parameters:
289: + sigma - singular value
290: . u - left singular vector
291: - v - right singular vector
293: Note:
294: Both u or v can be NULL if singular vectors are not required.
295: Otherwise, the caller must provide valid Vec objects, i.e.,
296: they must be created by the calling program with e.g. MatCreateVecs().
298: The index i should be a value between 0 and nconv-1 (see SVDGetConverged()).
299: Singular triplets are indexed according to the ordering criterion established
300: with SVDSetWhichSingularTriplets().
302: In the case of GSVD, the solution consists in three vectors u,v,x that are
303: returned as follows. Vector x is returned in the right singular vector
304: (argument v) and has length equal to the number of columns of A and B.
305: The other two vectors are returned stacked on top of each other [u;v] in
306: the left singular vector argument, with length equal to m+n (number of rows
307: of A plus number of rows of B).
309: Level: beginner
311: .seealso: SVDSolve(), SVDGetConverged(), SVDSetWhichSingularTriplets()
312: @*/
313: PetscErrorCode SVDGetSingularTriplet(SVD svd,PetscInt i,PetscReal *sigma,Vec u,Vec v)
314: {
315: PetscInt M,N;
316: Vec w;
318: PetscFunctionBegin;
321: SVDCheckSolved(svd,1);
324: PetscCheck(i>=0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"The index cannot be negative");
325: PetscCheck(i<svd->nconv,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"The index can be nconv-1 at most, see SVDGetConverged()");
326: if (sigma) *sigma = svd->sigma[svd->perm[i]];
327: if (u || v) {
328: if (!svd->isgeneralized) {
329: PetscCall(MatGetSize(svd->OP,&M,&N));
330: if (M<N) { w = u; u = v; v = w; }
331: }
332: PetscCall(SVDComputeVectors(svd));
333: if (u) PetscCall(BVCopyVec(svd->U,svd->perm[i],u));
334: if (v) PetscCall(BVCopyVec(svd->V,svd->perm[i],v));
335: }
336: PetscFunctionReturn(PETSC_SUCCESS);
337: }
339: /*
340: SVDComputeResidualNorms_Standard - Computes the norms of the left and
341: right residuals associated with the i-th computed singular triplet.
343: Input Parameters:
344: sigma - singular value
345: u,v - singular vectors
346: x,y - two work vectors with the same dimensions as u,v
347: */
348: static PetscErrorCode SVDComputeResidualNorms_Standard(SVD svd,PetscReal sigma,Vec u,Vec v,Vec x,Vec y,PetscReal *norm1,PetscReal *norm2)
349: {
350: PetscInt M,N;
352: PetscFunctionBegin;
353: /* norm1 = ||A*v-sigma*u||_2 */
354: if (norm1) {
355: PetscCall(MatMult(svd->OP,v,x));
356: PetscCall(VecAXPY(x,-sigma,u));
357: PetscCall(VecNorm(x,NORM_2,norm1));
358: }
359: /* norm2 = ||A^T*u-sigma*v||_2 */
360: if (norm2) {
361: PetscCall(MatGetSize(svd->OP,&M,&N));
362: if (M<N) PetscCall(MatMult(svd->A,u,y));
363: else PetscCall(MatMult(svd->AT,u,y));
364: PetscCall(VecAXPY(y,-sigma,v));
365: PetscCall(VecNorm(y,NORM_2,norm2));
366: }
367: PetscFunctionReturn(PETSC_SUCCESS);
368: }
370: /*
371: SVDComputeResidualNorms_Generalized - In GSVD, compute the residual norms
372: norm1 = ||s^2*A'*u-c*B'*B*x||_2 and norm2 = ||c^2*B'*v-s*A'*A*x||_2.
374: Input Parameters:
375: sigma - singular value
376: uv - left singular vectors [u;v]
377: x - right singular vector
378: y,z - two work vectors with the same dimension as x
379: */
380: static PetscErrorCode SVDComputeResidualNorms_Generalized(SVD svd,PetscReal sigma,Vec uv,Vec x,Vec y,Vec z,PetscReal *norm1,PetscReal *norm2)
381: {
382: Vec u,v,ax,bx,nest,aux[2];
383: PetscReal c,s;
385: PetscFunctionBegin;
386: PetscCall(MatCreateVecs(svd->OP,NULL,&u));
387: PetscCall(MatCreateVecs(svd->OPb,NULL,&v));
388: aux[0] = u;
389: aux[1] = v;
390: PetscCall(VecCreateNest(PetscObjectComm((PetscObject)svd),2,NULL,aux,&nest));
391: PetscCall(VecCopy(uv,nest));
393: s = 1.0/PetscSqrtReal(1.0+sigma*sigma);
394: c = sigma*s;
396: /* norm1 = ||s^2*A'*u-c*B'*B*x||_2 */
397: if (norm1) {
398: PetscCall(VecDuplicate(v,&bx));
399: PetscCall(MatMultHermitianTranspose(svd->OP,u,z));
400: PetscCall(MatMult(svd->OPb,x,bx));
401: PetscCall(MatMultHermitianTranspose(svd->OPb,bx,y));
402: PetscCall(VecAXPBY(y,s*s,-c,z));
403: PetscCall(VecNorm(y,NORM_2,norm1));
404: PetscCall(VecDestroy(&bx));
405: }
406: /* norm2 = ||c^2*B'*v-s*A'*A*x||_2 */
407: if (norm2) {
408: PetscCall(VecDuplicate(u,&ax));
409: PetscCall(MatMultHermitianTranspose(svd->OPb,v,z));
410: PetscCall(MatMult(svd->OP,x,ax));
411: PetscCall(MatMultHermitianTranspose(svd->OP,ax,y));
412: PetscCall(VecAXPBY(y,c*c,-s,z));
413: PetscCall(VecNorm(y,NORM_2,norm2));
414: PetscCall(VecDestroy(&ax));
415: }
417: PetscCall(VecDestroy(&v));
418: PetscCall(VecDestroy(&u));
419: PetscCall(VecDestroy(&nest));
420: PetscFunctionReturn(PETSC_SUCCESS);
421: }
423: /*
424: SVDComputeResidualNorms_Hyperbolic - Computes the norms of the left and
425: right residuals associated with the i-th computed singular triplet.
427: Input Parameters:
428: sigma - singular value
429: sign - corresponding element of the signature Omega2
430: u,v - singular vectors
431: x,y,z - three work vectors with the same dimensions as u,v,u
432: */
433: static PetscErrorCode SVDComputeResidualNorms_Hyperbolic(SVD svd,PetscReal sigma,PetscReal sign,Vec u,Vec v,Vec x,Vec y,Vec z,PetscReal *norm1,PetscReal *norm2)
434: {
435: PetscInt M,N;
437: PetscFunctionBegin;
438: /* norm1 = ||A*v-sigma*u||_2 */
439: if (norm1) {
440: PetscCall(MatMult(svd->OP,v,x));
441: PetscCall(VecAXPY(x,-sigma,u));
442: PetscCall(VecNorm(x,NORM_2,norm1));
443: }
444: /* norm2 = ||A^T*Omega*u-sigma*sign*v||_2 */
445: if (norm2) {
446: PetscCall(MatGetSize(svd->OP,&M,&N));
447: PetscCall(VecPointwiseMult(z,u,svd->omega));
448: if (M<N) PetscCall(MatMult(svd->A,z,y));
449: else PetscCall(MatMult(svd->AT,z,y));
450: PetscCall(VecAXPY(y,-sigma*sign,v));
451: PetscCall(VecNorm(y,NORM_2,norm2));
452: }
453: PetscFunctionReturn(PETSC_SUCCESS);
454: }
456: /*@
457: SVDComputeError - Computes the error (based on the residual norm) associated
458: with the i-th singular triplet.
460: Collective
462: Input Parameters:
463: + svd - the singular value solver context
464: . i - the solution index
465: - type - the type of error to compute
467: Output Parameter:
468: . error - the error
470: Notes:
471: The error can be computed in various ways, all of them based on the residual
472: norm obtained as sqrt(n1^2+n2^2) with n1 = ||A*v-sigma*u||_2 and
473: n2 = ||A^T*u-sigma*v||_2, where sigma is the singular value, u is the left
474: singular vector and v is the right singular vector.
476: In the case of the GSVD, the two components of the residual norm are
477: n1 = ||s^2*A'*u-c*B'*B*x||_2 and n2 = ||c^2*B'*v-s*A'*A*x||_2, where [u;v]
478: are the left singular vectors and x is the right singular vector, with
479: sigma=c/s.
481: Level: beginner
483: .seealso: SVDErrorType, SVDSolve()
484: @*/
485: PetscErrorCode SVDComputeError(SVD svd,PetscInt i,SVDErrorType type,PetscReal *error)
486: {
487: PetscReal sigma,norm1,norm2;
488: Vec u=NULL,v=NULL,x=NULL,y=NULL,z=NULL;
490: PetscFunctionBegin;
495: SVDCheckSolved(svd,1);
497: /* allocate work vectors */
498: switch (svd->problem_type) {
499: case SVD_STANDARD:
500: PetscCall(SVDSetWorkVecs(svd,2,2));
501: u = svd->workl[0];
502: v = svd->workr[0];
503: x = svd->workl[1];
504: y = svd->workr[1];
505: break;
506: case SVD_GENERALIZED:
507: PetscCheck(type!=SVD_ERROR_RELATIVE,PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"In GSVD the error should be either absolute or relative to the norms");
508: PetscCall(SVDSetWorkVecs(svd,1,3));
509: u = svd->workl[0];
510: v = svd->workr[0];
511: x = svd->workr[1];
512: y = svd->workr[2];
513: break;
514: case SVD_HYPERBOLIC:
515: PetscCall(SVDSetWorkVecs(svd,3,2));
516: u = svd->workl[0];
517: v = svd->workr[0];
518: x = svd->workl[1];
519: y = svd->workr[1];
520: z = svd->workl[2];
521: break;
522: }
524: /* compute residual norm and error */
525: PetscCall(SVDGetSingularTriplet(svd,i,&sigma,u,v));
526: switch (svd->problem_type) {
527: case SVD_STANDARD:
528: PetscCall(SVDComputeResidualNorms_Standard(svd,sigma,u,v,x,y,&norm1,&norm2));
529: break;
530: case SVD_GENERALIZED:
531: PetscCall(SVDComputeResidualNorms_Generalized(svd,sigma,u,v,x,y,&norm1,&norm2));
532: break;
533: case SVD_HYPERBOLIC:
534: PetscCall(SVDComputeResidualNorms_Hyperbolic(svd,sigma,svd->sign[svd->perm[i]],u,v,x,y,z,&norm1,&norm2));
535: break;
536: }
537: *error = SlepcAbs(norm1,norm2);
538: switch (type) {
539: case SVD_ERROR_ABSOLUTE:
540: break;
541: case SVD_ERROR_RELATIVE:
542: *error /= sigma;
543: break;
544: case SVD_ERROR_NORM:
545: if (!svd->nrma) PetscCall(MatNorm(svd->OP,NORM_INFINITY,&svd->nrma));
546: if (svd->isgeneralized && !svd->nrmb) PetscCall(MatNorm(svd->OPb,NORM_INFINITY,&svd->nrmb));
547: *error /= PetscMax(svd->nrma,svd->nrmb);
548: break;
549: default:
550: SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Invalid error type");
551: }
552: PetscFunctionReturn(PETSC_SUCCESS);
553: }