Actual source code: svdsetup.c
slepc-3.18.3 2023-03-24
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 for setting up the solver
12: */
14: #include <slepc/private/svdimpl.h>
16: /*@
17: SVDSetOperators - Set the matrices associated with the singular value problem.
19: Collective on svd
21: Input Parameters:
22: + svd - the singular value solver context
23: . A - the matrix associated with the singular value problem
24: - B - the second matrix in the case of GSVD
26: Level: beginner
28: .seealso: SVDSolve(), SVDGetOperators()
29: @*/
30: PetscErrorCode SVDSetOperators(SVD svd,Mat A,Mat B)
31: {
32: PetscInt Ma,Na,Mb,Nb,ma,na,mb,nb,M0,N0,m0,n0;
33: PetscBool samesize=PETSC_TRUE;
41: /* Check matrix sizes */
42: MatGetSize(A,&Ma,&Na);
43: MatGetLocalSize(A,&ma,&na);
44: if (svd->OP) {
45: MatGetSize(svd->OP,&M0,&N0);
46: MatGetLocalSize(svd->OP,&m0,&n0);
47: if (M0!=Ma || N0!=Na || m0!=ma || n0!=na) samesize = PETSC_FALSE;
48: }
49: if (B) {
50: MatGetSize(B,&Mb,&Nb);
51: MatGetLocalSize(B,&mb,&nb);
54: if (svd->OPb) {
55: MatGetSize(svd->OPb,&M0,&N0);
56: MatGetLocalSize(svd->OPb,&m0,&n0);
57: if (M0!=Mb || N0!=Nb || m0!=mb || n0!=nb) samesize = PETSC_FALSE;
58: }
59: }
61: PetscObjectReference((PetscObject)A);
62: if (B) PetscObjectReference((PetscObject)B);
63: if (svd->state && !samesize) SVDReset(svd);
64: else {
65: MatDestroy(&svd->OP);
66: MatDestroy(&svd->OPb);
67: MatDestroy(&svd->A);
68: MatDestroy(&svd->B);
69: MatDestroy(&svd->AT);
70: MatDestroy(&svd->BT);
71: }
72: svd->nrma = 0.0;
73: svd->nrmb = 0.0;
74: svd->OP = A;
75: svd->OPb = B;
76: svd->state = SVD_STATE_INITIAL;
77: return 0;
78: }
80: /*@
81: SVDGetOperators - Get the matrices associated with the singular value problem.
83: Collective on svd
85: Input Parameter:
86: . svd - the singular value solver context
88: Output Parameters:
89: + A - the matrix associated with the singular value problem
90: - B - the second matrix in the case of GSVD
92: Level: intermediate
94: .seealso: SVDSolve(), SVDSetOperators()
95: @*/
96: PetscErrorCode SVDGetOperators(SVD svd,Mat *A,Mat *B)
97: {
99: if (A) *A = svd->OP;
100: if (B) *B = svd->OPb;
101: return 0;
102: }
104: /*@
105: SVDSetSignature - Set the signature matrix defining a hyperbolic singular value problem.
107: Collective on svd
109: Input Parameters:
110: + svd - the singular value solver context
111: - omega - a vector containing the diagonal elements of the signature matrix (or NULL)
113: Notes:
114: The signature matrix is relevant only for hyperbolic problems (HSVD).
115: Use NULL to reset a previously set signature.
117: Level: intermediate
119: .seealso: SVDSetProblemType(), SVDSetOperators(), SVDGetSignature()
120: @*/
121: PetscErrorCode SVDSetSignature(SVD svd,Vec omega)
122: {
123: PetscInt N,Ma,n,ma;
126: if (omega) {
129: }
131: if (omega && svd->OP) { /* Check sizes */
132: VecGetSize(omega,&N);
133: VecGetLocalSize(omega,&n);
134: MatGetSize(svd->OP,&Ma,NULL);
135: MatGetLocalSize(svd->OP,&ma,NULL);
138: }
140: if (omega) PetscObjectReference((PetscObject)omega);
141: VecDestroy(&svd->omega);
142: svd->omega = omega;
143: svd->state = SVD_STATE_INITIAL;
144: return 0;
145: }
147: /*@
148: SVDGetSignature - Get the signature matrix defining a hyperbolic singular value problem.
150: Collective on svd
152: Input Parameter:
153: . svd - the singular value solver context
155: Output Parameter:
156: . omega - a vector containing the diagonal elements of the signature matrix
158: Level: intermediate
160: .seealso: SVDSetSignature()
161: @*/
162: PetscErrorCode SVDGetSignature(SVD svd,Vec *omega)
163: {
166: *omega = svd->omega;
167: return 0;
168: }
170: /*@
171: SVDSetUp - Sets up all the internal data structures necessary for the
172: execution of the singular value solver.
174: Collective on svd
176: Input Parameter:
177: . svd - singular value solver context
179: Notes:
180: This function need not be called explicitly in most cases, since SVDSolve()
181: calls it. It can be useful when one wants to measure the set-up time
182: separately from the solve time.
184: Level: developer
186: .seealso: SVDCreate(), SVDSolve(), SVDDestroy()
187: @*/
188: PetscErrorCode SVDSetUp(SVD svd)
189: {
190: PetscBool flg;
191: PetscInt M,N,P=0,k,maxnsol;
192: SlepcSC sc;
193: Vec *T;
194: BV bv;
197: if (svd->state) return 0;
198: PetscLogEventBegin(SVD_SetUp,svd,0,0,0);
200: /* reset the convergence flag from the previous solves */
201: svd->reason = SVD_CONVERGED_ITERATING;
203: /* set default solver type (SVDSetFromOptions was not called) */
204: if (!((PetscObject)svd)->type_name) SVDSetType(svd,SVDCROSS);
205: if (!svd->ds) SVDGetDS(svd,&svd->ds);
207: /* check matrices */
210: if (!svd->problem_type) { /* set default problem type */
211: if (svd->OPb) {
213: SVDSetProblemType(svd,SVD_GENERALIZED);
214: } else {
215: if (svd->omega) SVDSetProblemType(svd,SVD_HYPERBOLIC);
216: else SVDSetProblemType(svd,SVD_STANDARD);
217: }
218: } else { /* check consistency of problem type set by user */
219: if (svd->OPb) {
222: } else {
226: }
227: }
229: /* determine how to handle the transpose */
230: svd->expltrans = PETSC_TRUE;
231: if (svd->impltrans) svd->expltrans = PETSC_FALSE;
232: else {
233: MatHasOperation(svd->OP,MATOP_TRANSPOSE,&flg);
234: if (!flg) svd->expltrans = PETSC_FALSE;
235: else {
236: PetscObjectTypeCompareAny((PetscObject)svd,&flg,SVDLAPACK,SVDSCALAPACK,SVDELEMENTAL,"");
237: if (flg) svd->expltrans = PETSC_FALSE;
238: }
239: }
241: /* get matrix dimensions */
242: MatGetSize(svd->OP,&M,&N);
243: if (svd->isgeneralized) {
244: MatGetSize(svd->OPb,&P,NULL);
246: }
248: /* build transpose matrix */
249: MatDestroy(&svd->A);
250: MatDestroy(&svd->AT);
251: PetscObjectReference((PetscObject)svd->OP);
252: if (svd->expltrans) {
253: if (svd->isgeneralized || M>=N) {
254: svd->A = svd->OP;
255: MatHermitianTranspose(svd->OP,MAT_INITIAL_MATRIX,&svd->AT);
256: } else {
257: MatHermitianTranspose(svd->OP,MAT_INITIAL_MATRIX,&svd->A);
258: svd->AT = svd->OP;
259: }
260: } else {
261: if (svd->isgeneralized || M>=N) {
262: svd->A = svd->OP;
263: MatCreateHermitianTranspose(svd->OP,&svd->AT);
264: } else {
265: MatCreateHermitianTranspose(svd->OP,&svd->A);
266: svd->AT = svd->OP;
267: }
268: }
270: /* build transpose matrix B for GSVD */
271: if (svd->isgeneralized) {
272: MatDestroy(&svd->B);
273: MatDestroy(&svd->BT);
274: PetscObjectReference((PetscObject)svd->OPb);
275: if (svd->expltrans) {
276: svd->B = svd->OPb;
277: MatHermitianTranspose(svd->OPb,MAT_INITIAL_MATRIX,&svd->BT);
278: } else {
279: svd->B = svd->OPb;
280: MatCreateHermitianTranspose(svd->OPb,&svd->BT);
281: }
282: }
284: if (!svd->isgeneralized && M<N) {
285: /* swap initial vectors */
286: if (svd->nini || svd->ninil) {
287: T=svd->ISL; svd->ISL=svd->IS; svd->IS=T;
288: k=svd->ninil; svd->ninil=svd->nini; svd->nini=k;
289: }
290: /* swap basis vectors */
291: if (!svd->swapped) { /* only the first time in case of multiple calls */
292: bv=svd->V; svd->V=svd->U; svd->U=bv;
293: svd->swapped = PETSC_TRUE;
294: }
295: }
297: maxnsol = svd->isgeneralized? PetscMin(PetscMin(M,N),P): PetscMin(M,N);
298: svd->ncv = PetscMin(svd->ncv,maxnsol);
299: svd->nsv = PetscMin(svd->nsv,maxnsol);
302: /* relative convergence criterion is not allowed in GSVD */
303: if (svd->conv==(SVDConv)-1) SVDSetConvergenceTest(svd,svd->isgeneralized?SVD_CONV_NORM:SVD_CONV_REL);
306: /* initialization of matrix norm (stardard case only, for GSVD it is done inside setup()) */
307: if (!svd->isgeneralized && svd->conv==SVD_CONV_NORM && !svd->nrma) MatNorm(svd->OP,NORM_INFINITY,&svd->nrma);
309: /* call specific solver setup */
310: PetscUseTypeMethod(svd,setup);
312: /* set tolerance if not yet set */
313: if (svd->tol==PETSC_DEFAULT) svd->tol = SLEPC_DEFAULT_TOL;
315: /* fill sorting criterion context */
316: DSGetSlepcSC(svd->ds,&sc);
317: sc->comparison = (svd->which==SVD_LARGEST)? SlepcCompareLargestReal: SlepcCompareSmallestReal;
318: sc->comparisonctx = NULL;
319: sc->map = NULL;
320: sc->mapobj = NULL;
322: /* process initial vectors */
323: if (svd->nini<0) {
324: k = -svd->nini;
326: BVInsertVecs(svd->V,0,&k,svd->IS,PETSC_TRUE);
327: SlepcBasisDestroy_Private(&svd->nini,&svd->IS);
328: svd->nini = k;
329: }
330: if (svd->ninil<0) {
331: k = 0;
332: if (svd->leftbasis) {
333: k = -svd->ninil;
335: BVInsertVecs(svd->U,0,&k,svd->ISL,PETSC_TRUE);
336: } else PetscInfo(svd,"Ignoring initial left vectors\n");
337: SlepcBasisDestroy_Private(&svd->ninil,&svd->ISL);
338: svd->ninil = k;
339: }
341: PetscLogEventEnd(SVD_SetUp,svd,0,0,0);
342: svd->state = SVD_STATE_SETUP;
343: return 0;
344: }
346: /*@C
347: SVDSetInitialSpaces - Specify two basis of vectors that constitute the initial
348: right and/or left spaces.
350: Collective on svd
352: Input Parameters:
353: + svd - the singular value solver context
354: . nr - number of right vectors
355: . isr - set of basis vectors of the right initial space
356: . nl - number of left vectors
357: - isl - set of basis vectors of the left initial space
359: Notes:
360: The initial right and left spaces are rough approximations to the right and/or
361: left singular subspaces from which the solver starts to iterate.
362: It is not necessary to provide both sets of vectors.
364: Some solvers start to iterate on a single vector (initial vector). In that case,
365: the other vectors are ignored.
367: These vectors do not persist from one SVDSolve() call to the other, so the
368: initial space should be set every time.
370: The vectors do not need to be mutually orthonormal, since they are explicitly
371: orthonormalized internally.
373: Common usage of this function is when the user can provide a rough approximation
374: of the wanted singular space. Then, convergence may be faster.
376: Level: intermediate
378: .seealso: SVDSetUp()
379: @*/
380: PetscErrorCode SVDSetInitialSpaces(SVD svd,PetscInt nr,Vec isr[],PetscInt nl,Vec isl[])
381: {
386: if (nr>0) {
389: }
391: if (nl>0) {
394: }
395: SlepcBasisReference_Private(nr,isr,&svd->nini,&svd->IS);
396: SlepcBasisReference_Private(nl,isl,&svd->ninil,&svd->ISL);
397: if (nr>0 || nl>0) svd->state = SVD_STATE_INITIAL;
398: return 0;
399: }
401: /*
402: SVDSetDimensions_Default - Set reasonable values for ncv, mpd if not set
403: by the user. This is called at setup.
404: */
405: PetscErrorCode SVDSetDimensions_Default(SVD svd)
406: {
407: PetscInt N,M,P,maxnsol;
409: MatGetSize(svd->OP,&M,&N);
410: maxnsol = PetscMin(M,N);
411: if (svd->isgeneralized) {
412: MatGetSize(svd->OPb,&P,NULL);
413: maxnsol = PetscMin(maxnsol,P);
414: }
415: if (svd->ncv!=PETSC_DEFAULT) { /* ncv set */
417: } else if (svd->mpd!=PETSC_DEFAULT) { /* mpd set */
418: svd->ncv = PetscMin(maxnsol,svd->nsv+svd->mpd);
419: } else { /* neither set: defaults depend on nsv being small or large */
420: if (svd->nsv<500) svd->ncv = PetscMin(maxnsol,PetscMax(2*svd->nsv,10));
421: else {
422: svd->mpd = 500;
423: svd->ncv = PetscMin(maxnsol,svd->nsv+svd->mpd);
424: }
425: }
426: if (svd->mpd==PETSC_DEFAULT) svd->mpd = svd->ncv;
427: return 0;
428: }
430: /*@
431: SVDAllocateSolution - Allocate memory storage for common variables such
432: as the singular values and the basis vectors.
434: Collective on svd
436: Input Parameters:
437: + svd - eigensolver context
438: - extra - number of additional positions, used for methods that require a
439: working basis slightly larger than ncv
441: Developer Notes:
442: This is SLEPC_EXTERN because it may be required by user plugin SVD
443: implementations.
445: This is called at setup after setting the value of ncv and the flag leftbasis.
447: Level: developer
449: .seealso: SVDSetUp()
450: @*/
451: PetscErrorCode SVDAllocateSolution(SVD svd,PetscInt extra)
452: {
453: PetscInt oldsize,requested;
454: Vec tr,tl;
456: requested = svd->ncv + extra;
458: /* oldsize is zero if this is the first time setup is called */
459: BVGetSizes(svd->V,NULL,NULL,&oldsize);
461: /* allocate sigma */
462: if (requested != oldsize || !svd->sigma) {
463: PetscFree3(svd->sigma,svd->perm,svd->errest);
464: if (svd->sign) PetscFree(svd->sign);
465: PetscMalloc3(requested,&svd->sigma,requested,&svd->perm,requested,&svd->errest);
466: if (svd->ishyperbolic) PetscMalloc1(requested,&svd->sign);
467: }
468: /* allocate V */
469: if (!svd->V) SVDGetBV(svd,&svd->V,NULL);
470: if (!oldsize) {
471: if (!((PetscObject)(svd->V))->type_name) BVSetType(svd->V,BVSVEC);
472: MatCreateVecsEmpty(svd->A,&tr,NULL);
473: BVSetSizesFromVec(svd->V,tr,requested);
474: VecDestroy(&tr);
475: } else BVResize(svd->V,requested,PETSC_FALSE);
476: /* allocate U */
477: if (svd->leftbasis && !svd->isgeneralized) {
478: if (!svd->U) SVDGetBV(svd,NULL,&svd->U);
479: if (!oldsize) {
480: if (!((PetscObject)(svd->U))->type_name) BVSetType(svd->U,((PetscObject)(svd->V))->type_name);
481: MatCreateVecsEmpty(svd->A,NULL,&tl);
482: BVSetSizesFromVec(svd->U,tl,requested);
483: VecDestroy(&tl);
484: } else BVResize(svd->U,requested,PETSC_FALSE);
485: } else if (svd->isgeneralized) { /* left basis for the GSVD */
486: if (!svd->U) SVDGetBV(svd,NULL,&svd->U);
487: if (!oldsize) {
488: if (!((PetscObject)(svd->U))->type_name) BVSetType(svd->U,((PetscObject)(svd->V))->type_name);
489: SVDCreateLeftTemplate(svd,&tl);
490: BVSetSizesFromVec(svd->U,tl,requested);
491: VecDestroy(&tl);
492: } else BVResize(svd->U,requested,PETSC_FALSE);
493: }
494: return 0;
495: }