Actual source code: svdsetup.c
1: /*
2: SVD routines for setting up the solver.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
9:
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 <private/svdimpl.h> /*I "slepcsvd.h" I*/
25: #include <private/ipimpl.h>
29: /*@
30: SVDSetOperator - Set the matrix associated with the singular value problem.
32: Collective on SVD and Mat
34: Input Parameters:
35: + svd - the singular value solver context
36: - A - the matrix associated with the singular value problem
38: Level: beginner
40: .seealso: SVDSolve(), SVDGetOperator()
41: @*/
42: PetscErrorCode SVDSetOperator(SVD svd,Mat mat)
43: {
50: if (svd->setupcalled) { SVDReset(svd); }
51: PetscObjectReference((PetscObject)mat);
52: MatDestroy(&svd->OP);
53: svd->OP = mat;
54: return(0);
55: }
59: /*@
60: SVDGetOperator - Get the matrix associated with the singular value problem.
62: Not collective, though parallel Mats are returned if the SVD is parallel
64: Input Parameter:
65: . svd - the singular value solver context
67: Output Parameters:
68: . A - the matrix associated with the singular value problem
70: Level: advanced
72: .seealso: SVDSolve(), SVDSetOperator()
73: @*/
74: PetscErrorCode SVDGetOperator(SVD svd,Mat *A)
75: {
79: *A = svd->OP;
80: return(0);
81: }
85: /*@
86: SVDSetUp - Sets up all the internal data structures necessary for the
87: execution of the singular value solver.
89: Collective on SVD
91: Input Parameter:
92: . svd - singular value solver context
94: Level: advanced
96: Notes:
97: This function need not be called explicitly in most cases, since SVDSolve()
98: calls it. It can be useful when one wants to measure the set-up time
99: separately from the solve time.
101: .seealso: SVDCreate(), SVDSolve(), SVDDestroy()
102: @*/
103: PetscErrorCode SVDSetUp(SVD svd)
104: {
106: PetscBool flg,lindep;
107: PetscInt i,k,M,N;
108: PetscReal norm;
109:
112: if (svd->setupcalled) return(0);
113: PetscLogEventBegin(SVD_SetUp,svd,0,0,0);
115: /* Set default solver type (SVDSetFromOptions was not called) */
116: if (!((PetscObject)svd)->type_name) {
117: SVDSetType(svd,SVDCROSS);
118: }
119: if (!svd->ip) { SVDGetIP(svd,&svd->ip); }
120: if (!((PetscObject)svd->ip)->type_name) {
121: IPSetDefaultType_Private(svd->ip);
122: }
123: if (!((PetscObject)svd->rand)->type_name) {
124: PetscRandomSetFromOptions(svd->rand);
125: }
127: /* check matrix */
128: if (!svd->OP) SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_WRONGSTATE,"SVDSetOperator must be called first");
129:
130: /* determine how to build the transpose */
131: if (svd->transmode == PETSC_DECIDE) {
132: MatHasOperation(svd->OP,MATOP_TRANSPOSE,&flg);
133: if (flg) svd->transmode = SVD_TRANSPOSE_EXPLICIT;
134: else svd->transmode = SVD_TRANSPOSE_IMPLICIT;
135: }
136:
137: /* build transpose matrix */
138: MatDestroy(&svd->A);
139: MatDestroy(&svd->AT);
140: MatGetSize(svd->OP,&M,&N);
141: PetscObjectReference((PetscObject)svd->OP);
142: switch (svd->transmode) {
143: case SVD_TRANSPOSE_EXPLICIT:
144: MatHasOperation(svd->OP,MATOP_TRANSPOSE,&flg);
145: if (!flg) SETERRQ(((PetscObject)svd)->comm,1,"Matrix has not defined the MatTranpose operation");
146: if (M>=N) {
147: svd->A = svd->OP;
148: MatTranspose(svd->OP,MAT_INITIAL_MATRIX,&svd->AT);
149: } else {
150: MatTranspose(svd->OP,MAT_INITIAL_MATRIX,&svd->A);
151: svd->AT = svd->OP;
152: }
153: break;
154: case SVD_TRANSPOSE_IMPLICIT:
155: MatHasOperation(svd->OP,MATOP_MULT_TRANSPOSE,&flg);
156: if (!flg) SETERRQ(((PetscObject)svd)->comm,1,"Matrix has not defined the MatMultTranpose operation");
157: if (M>=N) {
158: svd->A = svd->OP;
159: svd->AT = PETSC_NULL;
160: } else {
161: svd->A = PETSC_NULL;
162: svd->AT = svd->OP;
163: }
164: break;
165: default:
166: SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid transpose mode");
167: }
169: VecDestroy(&svd->tr);
170: VecDestroy(&svd->tl);
171: if (svd->A) {
172: SlepcMatGetVecsTemplate(svd->A,&svd->tr,&svd->tl);
173: } else {
174: SlepcMatGetVecsTemplate(svd->AT,&svd->tl,&svd->tr);
175: }
177: /* call specific solver setup */
178: (*svd->ops->setup)(svd);
180: /* set tolerance if not yet set */
181: if (svd->tol==PETSC_DEFAULT) svd->tol = SLEPC_DEFAULT_TOL;
183: if (svd->ncv > M || svd->ncv > N)
184: SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_OUTOFRANGE,"ncv bigger than matrix dimensions");
185: if (svd->nsv > svd->ncv)
186: SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nsv bigger than ncv");
188: if (svd->ncv != svd->n) {
189: /* free memory for previous solution */
190: if (svd->n) {
191: PetscFree(svd->sigma);
192: PetscFree(svd->perm);
193: PetscFree(svd->errest);
194: VecDestroyVecs(svd->n,&svd->V);
195: }
196: /* allocate memory for next solution */
197: PetscMalloc(svd->ncv*sizeof(PetscReal),&svd->sigma);
198: PetscMalloc(svd->ncv*sizeof(PetscInt),&svd->perm);
199: PetscMalloc(svd->ncv*sizeof(PetscReal),&svd->errest);
200: VecDuplicateVecs(svd->tr,svd->ncv,&svd->V);
201: svd->n = svd->ncv;
202: }
204: /* process initial vectors */
205: if (svd->nini<0) {
206: svd->nini = -svd->nini;
207: if (svd->nini>svd->ncv) SETERRQ(((PetscObject)svd)->comm,1,"The number of initial vectors is larger than ncv");
208: k = 0;
209: for (i=0;i<svd->nini;i++) {
210: VecCopy(svd->IS[i],svd->V[k]);
211: VecDestroy(&svd->IS[i]);
212: IPOrthogonalize(svd->ip,0,PETSC_NULL,k,PETSC_NULL,svd->V,svd->V[k],PETSC_NULL,&norm,&lindep);
213: if (norm==0.0 || lindep) PetscInfo(svd,"Linearly dependent initial vector found, removing...\n");
214: else {
215: VecScale(svd->V[k],1.0/norm);
216: k++;
217: }
218: }
219: svd->nini = k;
220: PetscFree(svd->IS);
221: }
223: PetscLogEventEnd(SVD_SetUp,svd,0,0,0);
224: svd->setupcalled = 1;
225: return(0);
226: }
230: /*@
231: SVDSetInitialSpace - Specify a basis of vectors that constitute the initial
232: space, that is, the subspace from which the solver starts to iterate.
234: Collective on SVD and Vec
236: Input Parameter:
237: + svd - the singular value solver context
238: . n - number of vectors
239: - is - set of basis vectors of the initial space
241: Notes:
242: Some solvers start to iterate on a single vector (initial vector). In that case,
243: the other vectors are ignored.
245: These vectors do not persist from one SVDSolve() call to the other, so the
246: initial space should be set every time.
248: The vectors do not need to be mutually orthonormal, since they are explicitly
249: orthonormalized internally.
251: Common usage of this function is when the user can provide a rough approximation
252: of the wanted singular space. Then, convergence may be faster.
254: Level: intermediate
255: @*/
256: PetscErrorCode SVDSetInitialSpace(SVD svd,PetscInt n,Vec *is)
257: {
259: PetscInt i;
260:
264: if (n<0) SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
266: /* free previous non-processed vectors */
267: if (svd->nini<0) {
268: for (i=0;i<-svd->nini;i++) {
269: VecDestroy(&svd->IS[i]);
270: }
271: PetscFree(svd->IS);
272: }
274: /* get references of passed vectors */
275: PetscMalloc(n*sizeof(Vec),&svd->IS);
276: for (i=0;i<n;i++) {
277: PetscObjectReference((PetscObject)is[i]);
278: svd->IS[i] = is[i];
279: }
281: svd->nini = -n;
282: svd->setupcalled = 0;
283: return(0);
284: }