1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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> /*I "slepcsvd.h" I*/
16: /*@
17: SVDSetOperator - Set the matrix 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
25: Level: beginner
27: .seealso: SVDSolve(), SVDGetOperator()
28: @*/
29: PetscErrorCode SVDSetOperator(SVD svd,Mat mat) 30: {
37: PetscObjectReference((PetscObject)mat);
38: if (svd->state) { SVDReset(svd); }
39: else { MatDestroy(&svd->OP); }
40: svd->OP = mat;
41: svd->state = SVD_STATE_INITIAL;
42: return(0);
43: }
45: /*@
46: SVDGetOperator - Get the matrix associated with the singular value problem.
48: Not collective, though parallel Mats are returned if the SVD is parallel
50: Input Parameter:
51: . svd - the singular value solver context
53: Output Parameters:
54: . A - the matrix associated with the singular value problem
56: Level: advanced
58: .seealso: SVDSolve(), SVDSetOperator()
59: @*/
60: PetscErrorCode SVDGetOperator(SVD svd,Mat *A) 61: {
65: *A = svd->OP;
66: return(0);
67: }
69: /*@
70: SVDSetUp - Sets up all the internal data structures necessary for the
71: execution of the singular value solver.
73: Collective on svd
75: Input Parameter:
76: . svd - singular value solver context
78: Notes:
79: This function need not be called explicitly in most cases, since SVDSolve()
80: calls it. It can be useful when one wants to measure the set-up time
81: separately from the solve time.
83: Level: developer
85: .seealso: SVDCreate(), SVDSolve(), SVDDestroy()
86: @*/
87: PetscErrorCode SVDSetUp(SVD svd) 88: {
90: PetscBool expltrans,flg;
91: PetscInt M,N,k;
92: SlepcSC sc;
93: Vec *T;
94: BV bv;
98: if (svd->state) return(0);
99: PetscLogEventBegin(SVD_SetUp,svd,0,0,0);
101: /* reset the convergence flag from the previous solves */
102: svd->reason = SVD_CONVERGED_ITERATING;
104: /* Set default solver type (SVDSetFromOptions was not called) */
105: if (!((PetscObject)svd)->type_name) {
106: SVDSetType(svd,SVDCROSS);
107: }
108: if (!svd->ds) { SVDGetDS(svd,&svd->ds); }
110: /* check matrix */
111: if (!svd->OP) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_WRONGSTATE,"SVDSetOperator must be called first");
113: /* determine how to handle the transpose */
114: expltrans = PETSC_TRUE;
115: if (svd->impltrans) expltrans = PETSC_FALSE;
116: else {
117: MatHasOperation(svd->OP,MATOP_TRANSPOSE,&flg);
118: if (!flg) expltrans = PETSC_FALSE;
119: else {
120: PetscObjectTypeCompare((PetscObject)svd,SVDLAPACK,&flg);
121: if (flg) expltrans = PETSC_FALSE;
122: }
123: }
125: /* build transpose matrix */
126: MatDestroy(&svd->A);
127: MatDestroy(&svd->AT);
128: MatGetSize(svd->OP,&M,&N);
129: PetscObjectReference((PetscObject)svd->OP);
130: if (expltrans) {
131: if (M>=N) {
132: svd->A = svd->OP;
133: MatTranspose(svd->OP,MAT_INITIAL_MATRIX,&svd->AT);
134: MatConjugate(svd->AT);
135: } else {
136: MatTranspose(svd->OP,MAT_INITIAL_MATRIX,&svd->A);
137: MatConjugate(svd->A);
138: svd->AT = svd->OP;
139: }
140: } else {
141: if (M>=N) {
142: svd->A = svd->OP;
143: svd->AT = NULL;
144: } else {
145: svd->A = NULL;
146: svd->AT = svd->OP;
147: }
148: }
150: if (M<N) {
151: /* swap initial vectors and basis vectors */
152: T=svd->ISL; svd->ISL=svd->IS; svd->IS=T;
153: k=svd->ninil; svd->ninil=svd->nini; svd->nini=k;
154: bv=svd->V; svd->V=svd->U; svd->U=bv;
155: }
157: if (svd->ncv > PetscMin(M,N)) svd->ncv = PetscMin(M,N);
158: if (svd->nsv > PetscMin(M,N)) svd->nsv = PetscMin(M,N);
159: if (svd->ncv && svd->nsv > svd->ncv) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"nsv bigger than ncv");
161: /* call specific solver setup */
162: (*svd->ops->setup)(svd);
164: /* set tolerance if not yet set */
165: if (svd->tol==PETSC_DEFAULT) svd->tol = SLEPC_DEFAULT_TOL;
167: /* fill sorting criterion context */
168: DSGetSlepcSC(svd->ds,&sc);
169: sc->comparison = (svd->which==SVD_LARGEST)? SlepcCompareLargestReal: SlepcCompareSmallestReal;
170: sc->comparisonctx = NULL;
171: sc->map = NULL;
172: sc->mapobj = NULL;
174: /* process initial vectors */
175: if (svd->nini<0) {
176: k = -svd->nini;
177: if (k>svd->ncv) SETERRQ(PetscObjectComm((PetscObject)svd),1,"The number of initial vectors is larger than ncv");
178: BVInsertVecs(svd->V,0,&k,svd->IS,PETSC_TRUE);
179: SlepcBasisDestroy_Private(&svd->nini,&svd->IS);
180: svd->nini = k;
181: }
182: if (svd->ninil<0) {
183: k = 0;
184: if (svd->leftbasis) {
185: k = -svd->ninil;
186: if (k>svd->ncv) SETERRQ(PetscObjectComm((PetscObject)svd),1,"The number of left initial vectors is larger than ncv");
187: BVInsertVecs(svd->U,0,&k,svd->ISL,PETSC_TRUE);
188: } else {
189: PetscInfo(svd,"Ignoring initial left vectors\n");
190: }
191: SlepcBasisDestroy_Private(&svd->ninil,&svd->ISL);
192: svd->ninil = k;
193: }
195: PetscLogEventEnd(SVD_SetUp,svd,0,0,0);
196: svd->state = SVD_STATE_SETUP;
197: return(0);
198: }
200: /*@C
201: SVDSetInitialSpaces - Specify two basis of vectors that constitute the initial
202: right and/or left spaces, that is, a rough approximation to the right and/or
203: left singular subspaces from which the solver starts to iterate.
205: Collective on svd
207: Input Parameter:
208: + svd - the singular value solver context
209: . nr - number of right vectors
210: . isr - set of basis vectors of the right initial space
211: . nl - number of left vectors
212: - isl - set of basis vectors of the left initial space
214: Notes:
215: It is not necessary to provide both sets of vectors.
217: Some solvers start to iterate on a single vector (initial vector). In that case,
218: the other vectors are ignored.
220: These vectors do not persist from one SVDSolve() call to the other, so the
221: initial space should be set every time.
223: The vectors do not need to be mutually orthonormal, since they are explicitly
224: orthonormalized internally.
226: Common usage of this function is when the user can provide a rough approximation
227: of the wanted singular space. Then, convergence may be faster.
229: Level: intermediate
230: @*/
231: PetscErrorCode SVDSetInitialSpaces(SVD svd,PetscInt nr,Vec *isr,PetscInt nl,Vec *isl)232: {
239: if (nr<0) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Argument nr cannot be negative");
240: if (nl<0) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Argument nl cannot be negative");
241: SlepcBasisReference_Private(nr,isr,&svd->nini,&svd->IS);
242: SlepcBasisReference_Private(nl,isl,&svd->ninil,&svd->ISL);
243: if (nr>0 || nl>0) svd->state = SVD_STATE_INITIAL;
244: return(0);
245: }
247: /*
248: SVDSetDimensions_Default - Set reasonable values for ncv, mpd if not set
249: by the user. This is called at setup.
250: */
251: PetscErrorCode SVDSetDimensions_Default(SVD svd)252: {
254: PetscInt N;
257: SVDMatGetSize(svd,NULL,&N);
258: if (svd->ncv) { /* ncv set */
259: if (svd->ncv<svd->nsv) SETERRQ(PetscObjectComm((PetscObject)svd),1,"The value of ncv must be at least nsv");
260: } else if (svd->mpd) { /* mpd set */
261: svd->ncv = PetscMin(N,svd->nsv+svd->mpd);
262: } else { /* neither set: defaults depend on nsv being small or large */
263: if (svd->nsv<500) svd->ncv = PetscMin(N,PetscMax(2*svd->nsv,10));
264: else {
265: svd->mpd = 500;
266: svd->ncv = PetscMin(N,svd->nsv+svd->mpd);
267: }
268: }
269: if (!svd->mpd) svd->mpd = svd->ncv;
270: return(0);
271: }
273: /*@
274: SVDAllocateSolution - Allocate memory storage for common variables such
275: as the singular values and the basis vectors.
277: Collective on svd
279: Input Parameters:
280: + svd - eigensolver context
281: - extra - number of additional positions, used for methods that require a
282: working basis slightly larger than ncv
284: Developers Notes:
285: This is SLEPC_EXTERN because it may be required by user plugin SVD286: implementations.
288: This is called at setup after setting the value of ncv and the flag leftbasis.
290: Level: developer
291: @*/
292: PetscErrorCode SVDAllocateSolution(SVD svd,PetscInt extra)293: {
295: PetscInt oldsize,requested;
296: Vec tr,tl;
299: requested = svd->ncv + extra;
301: /* oldsize is zero if this is the first time setup is called */
302: BVGetSizes(svd->V,NULL,NULL,&oldsize);
304: /* allocate sigma */
305: if (requested != oldsize || !svd->sigma) {
306: PetscFree3(svd->sigma,svd->perm,svd->errest);
307: PetscMalloc3(requested,&svd->sigma,requested,&svd->perm,requested,&svd->errest);
308: PetscLogObjectMemory((PetscObject)svd,PetscMax(0,requested-oldsize)*(2*sizeof(PetscReal)+sizeof(PetscInt)));
309: }
310: /* allocate V */
311: if (!svd->V) { SVDGetBV(svd,&svd->V,NULL); }
312: if (!oldsize) {
313: if (!((PetscObject)(svd->V))->type_name) {
314: BVSetType(svd->V,BVSVEC);
315: }
316: SVDMatCreateVecsEmpty(svd,&tr,NULL);
317: BVSetSizesFromVec(svd->V,tr,requested);
318: VecDestroy(&tr);
319: } else {
320: BVResize(svd->V,requested,PETSC_FALSE);
321: }
322: /* allocate U */
323: if (svd->leftbasis) {
324: if (!svd->U) { SVDGetBV(svd,NULL,&svd->U); }
325: if (!oldsize) {
326: if (!((PetscObject)(svd->U))->type_name) {
327: BVSetType(svd->U,BVSVEC);
328: }
329: SVDMatCreateVecsEmpty(svd,NULL,&tl);
330: BVSetSizesFromVec(svd->U,tl,requested);
331: VecDestroy(&tl);
332: } else {
333: BVResize(svd->U,requested,PETSC_FALSE);
334: }
335: }
336: return(0);
337: }