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-2007, Universidad Politecnica de Valencia, Spain
8: This file is part of SLEPc. See the README file for conditions of use
9: and additional information.
10: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
11: */
13: #include src/svd/svdimpl.h
17: /*@
18: SVDSetOperator - Set the matrix associated with the singular value problem.
20: Collective on SVD and Mat
22: Input Parameters:
23: + svd - the singular value solver context
24: - A - the matrix associated with the singular value problem
26: Level: beginner
28: .seealso: SVDSolve(), SVDGetOperator()
29: @*/
30: PetscErrorCode SVDSetOperator(SVD svd,Mat mat)
31: {
38: PetscObjectReference((PetscObject)mat);
39: if (svd->OP) {
40: MatDestroy(svd->OP);
41: }
42: svd->OP = mat;
43: if (svd->vec_initial) {
44: VecDestroy(svd->vec_initial);
45: svd->vec_initial = PETSC_NULL;
46: }
47: svd->setupcalled = 0;
48: return(0);
49: }
53: /*@
54: SVDGetOperator - Get the matrix associated with the singular value problem.
56: Not collective, though parallel Mats are returned if the SVD is parallel
58: Input Parameter:
59: . svd - the singular value solver context
61: Output Parameters:
62: . A - the matrix associated with the singular value problem
64: Level: advanced
66: .seealso: SVDSolve(), SVDSetOperator()
67: @*/
68: PetscErrorCode SVDGetOperator(SVD svd,Mat *A)
69: {
73: *A = svd->OP;
74: return(0);
75: }
79: /*@
80: SVDSetInitialVector - Sets the initial vector from which the
81: singular value solver starts to iterate.
83: Collective on SVD and Vec
85: Input Parameters:
86: + svd - the singular value solver context
87: - vec - the vector
89: Level: intermediate
91: .seealso: SVDGetInitialVector()
93: @*/
94: PetscErrorCode SVDSetInitialVector(SVD svd,Vec vec)
95: {
97:
102: PetscObjectReference((PetscObject)vec);
103: if (svd->vec_initial) {
104: VecDestroy(svd->vec_initial);
105: }
106: svd->vec_initial = vec;
107: return(0);
108: }
112: /*@
113: SVDGetInitialVector - Gets the initial vector associated with the
114: singular value solver; if the vector was not set it will return a 0
115: pointer or a vector randomly generated by SVDSetUp().
117: Not collective, but vector is shared by all processors that share the SVD
119: Input Parameter:
120: . svd - the singular value solver context
122: Output Parameter:
123: . vec - the vector
125: Level: intermediate
127: .seealso: SVDSetInitialVector()
129: @*/
130: PetscErrorCode SVDGetInitialVector(SVD svd,Vec *vec)
131: {
135: *vec = svd->vec_initial;
136: return(0);
137: }
141: /*@
142: SVDSetUp - Sets up all the internal data structures necessary for the
143: execution of the singular value solver.
145: Collective on SVD
147: Input Parameter:
148: . SVD - singular value solver context
150: Level: advanced
152: Notes:
153: This function need not be called explicitly in most cases, since SVDSolve()
154: calls it. It can be useful when one wants to measure the set-up time
155: separately from the solve time.
157: .seealso: SVDCreate(), SVDSolve(), SVDDestroy()
158: @*/
159: PetscErrorCode SVDSetUp(SVD svd)
160: {
162: int i;
163: PetscTruth flg;
164: PetscInt M,N;
165:
168: if (svd->setupcalled) return(0);
169: PetscLogEventBegin(SVD_SetUp,svd,0,0,0);
171: /* Set default solver type */
172: if (!svd->type_name) {
173: SVDSetType(svd,SVDCROSS);
174: }
176: /* check matrix */
177: if (!svd->OP)
178: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "SVDSetOperator must be called first");
179:
180: /* determine how to build the transpose */
181: if (svd->transmode == PETSC_DECIDE) {
182: MatHasOperation(svd->OP,MATOP_TRANSPOSE,&flg);
183: if (flg) svd->transmode = SVD_TRANSPOSE_EXPLICIT;
184: else svd->transmode = SVD_TRANSPOSE_IMPLICIT;
185: }
186:
187: /* build transpose matrix */
188: if (svd->A) { MatDestroy(svd->A); }
189: if (svd->AT) { MatDestroy(svd->AT); }
190: MatGetSize(svd->OP,&M,&N);
191: PetscObjectReference((PetscObject)svd->OP);
192: switch (svd->transmode) {
193: case SVD_TRANSPOSE_EXPLICIT:
194: MatHasOperation(svd->OP,MATOP_TRANSPOSE,&flg);
195: if (!flg) SETERRQ(1,"Matrix has not defined the MatTranpose operation");
196: if (M>=N) {
197: svd->A = svd->OP;
198: MatTranspose(svd->OP,&svd->AT);
199: } else {
200: MatTranspose(svd->OP,&svd->A);
201: svd->AT = svd->OP;
202: }
203: break;
204: case SVD_TRANSPOSE_IMPLICIT:
205: MatHasOperation(svd->OP,MATOP_MULT_TRANSPOSE,&flg);
206: if (!flg) SETERRQ(1,"Matrix has not defined the MatMultTranpose operation");
207: if (M>=N) {
208: svd->A = svd->OP;
209: svd->AT = PETSC_NULL;
210: } else {
211: svd->A = PETSC_NULL;
212: svd->AT = svd->OP;
213: }
214: break;
215: default:
216: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid transpose mode");
217: }
219: /* set initial vector */
220: if (!svd->vec_initial) {
221: SVDMatGetVecs(svd,&svd->vec_initial,PETSC_NULL);
222: SlepcVecSetRandom(svd->vec_initial);
223: }
225: /* call specific solver setup */
226: (*svd->ops->setup)(svd);
228: if (svd->ncv > M || svd->ncv > N)
229: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"ncv bigger than matrix dimensions");
230: if (svd->nsv > svd->ncv)
231: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"nsv bigger than ncv");
233: if (svd->ncv != svd->n) {
234: /* free memory for previous solution */
235: if (svd->n) {
236: PetscFree(svd->sigma);
237: PetscFree(svd->errest);
238: for (i=0;i<svd->n;i++) {
239: VecDestroy(svd->V[i]);
240: }
241: PetscFree(svd->V);
242: }
243: /* allocate memory for next solution */
244: PetscMalloc(svd->ncv*sizeof(PetscReal),&svd->sigma);
245: PetscMalloc(svd->ncv*sizeof(PetscReal),&svd->errest);
246: PetscMalloc(svd->ncv*sizeof(Vec),&svd->V);
247: for (i=0;i<svd->ncv;i++) {
248: SVDMatGetVecs(svd,svd->V+i,PETSC_NULL);
249: }
250: svd->n = svd->ncv;
251: }
253: PetscLogEventEnd(SVD_SetUp,svd,0,0,0);
254: svd->setupcalled = 1;
255: return(0);
256: }