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-2009, Universidad 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

 28: /*@
 29:    SVDSetOperator - Set the matrix associated with the singular value problem.

 31:    Collective on SVD and Mat

 33:    Input Parameters:
 34: +  svd - the singular value solver context
 35: -  A  - the matrix associated with the singular value problem

 37:    Level: beginner

 39: .seealso: SVDSolve(), SVDGetOperator()
 40: @*/
 41: PetscErrorCode SVDSetOperator(SVD svd,Mat mat)
 42: {

 49:   PetscObjectReference((PetscObject)mat);
 50:   if (svd->OP) {
 51:     MatDestroy(svd->OP);
 52:   }
 53:   svd->OP = mat;
 54:   if (svd->vec_initial) {
 55:     VecDestroy(svd->vec_initial);
 56:     svd->vec_initial = PETSC_NULL;
 57:   }
 58:   svd->setupcalled = 0;
 59:   return(0);
 60: }

 64: /*@
 65:    SVDGetOperator - Get the matrix associated with the singular value problem.

 67:    Not collective, though parallel Mats are returned if the SVD is parallel

 69:    Input Parameter:
 70: .  svd - the singular value solver context

 72:    Output Parameters:
 73: .  A    - the matrix associated with the singular value problem

 75:    Level: advanced

 77: .seealso: SVDSolve(), SVDSetOperator()
 78: @*/
 79: PetscErrorCode SVDGetOperator(SVD svd,Mat *A)
 80: {
 84:   *A = svd->OP;
 85:   return(0);
 86: }

 90: /*@
 91:    SVDSetInitialVector - Sets the initial vector from which the 
 92:    singular value solver starts to iterate.

 94:    Collective on SVD and Vec

 96:    Input Parameters:
 97: +  svd - the singular value solver context
 98: -  vec - the vector

100:    Level: intermediate

102: .seealso: SVDGetInitialVector()

104: @*/
105: PetscErrorCode SVDSetInitialVector(SVD svd,Vec vec)
106: {
108: 
113:   PetscObjectReference((PetscObject)vec);
114:   if (svd->vec_initial) {
115:     VecDestroy(svd->vec_initial);
116:   }
117:   svd->vec_initial = vec;
118:   return(0);
119: }

123: /*@
124:    SVDGetInitialVector - Gets the initial vector associated with the 
125:    singular value solver; if the vector was not set it will return a 0 
126:    pointer or a vector randomly generated by SVDSetUp().

128:    Not collective, but vector is shared by all processors that share the SVD

130:    Input Parameter:
131: .  svd - the singular value solver context

133:    Output Parameter:
134: .  vec - the vector

136:    Level: intermediate

138: .seealso: SVDSetInitialVector()

140: @*/
141: PetscErrorCode SVDGetInitialVector(SVD svd,Vec *vec)
142: {
146:   *vec = svd->vec_initial;
147:   return(0);
148: }

152: /*@
153:    SVDSetUp - Sets up all the internal data structures necessary for the
154:    execution of the singular value solver.

156:    Collective on SVD

158:    Input Parameter:
159: .  svd   - singular value solver context

161:    Level: advanced

163:    Notes:
164:    This function need not be called explicitly in most cases, since SVDSolve()
165:    calls it. It can be useful when one wants to measure the set-up time 
166:    separately from the solve time.

168: .seealso: SVDCreate(), SVDSolve(), SVDDestroy()
169: @*/
170: PetscErrorCode SVDSetUp(SVD svd)
171: {
173:   PetscTruth     flg;
174:   PetscInt       i,M,N,nloc;
175:   PetscScalar    *pV;
176: 
179:   if (svd->setupcalled) return(0);
180:   PetscLogEventBegin(SVD_SetUp,svd,0,0,0);

182:   /* Set default solver type */
183:   if (!((PetscObject)svd)->type_name) {
184:     SVDSetType(svd,SVDCROSS);
185:   }

187:   /* check matrix */
188:   if (!svd->OP)
189:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "SVDSetOperator must be called first");
190: 
191:   /* determine how to build the transpose */
192:   if (svd->transmode == PETSC_DECIDE) {
193:     MatHasOperation(svd->OP,MATOP_TRANSPOSE,&flg);
194:     if (flg) svd->transmode = SVD_TRANSPOSE_EXPLICIT;
195:     else svd->transmode = SVD_TRANSPOSE_IMPLICIT;
196:   }
197: 
198:   /* build transpose matrix */
199:   if (svd->A) { MatDestroy(svd->A); }
200:   if (svd->AT) { MatDestroy(svd->AT); }
201:   MatGetSize(svd->OP,&M,&N);
202:   PetscObjectReference((PetscObject)svd->OP);
203:   switch (svd->transmode) {
204:     case SVD_TRANSPOSE_EXPLICIT:
205:       MatHasOperation(svd->OP,MATOP_TRANSPOSE,&flg);
206:       if (!flg) SETERRQ(1,"Matrix has not defined the MatTranpose operation");
207:       if (M>=N) {
208:         svd->A = svd->OP;
209:         MatTranspose(svd->OP, MAT_INITIAL_MATRIX,&svd->AT);
210:       } else {
211:         MatTranspose(svd->OP, MAT_INITIAL_MATRIX,&svd->A);
212:         svd->AT = svd->OP;
213:       }
214:       break;
215:     case SVD_TRANSPOSE_IMPLICIT:
216:       MatHasOperation(svd->OP,MATOP_MULT_TRANSPOSE,&flg);
217:       if (!flg) SETERRQ(1,"Matrix has not defined the MatMultTranpose operation");
218:       if (M>=N) {
219:         svd->A = svd->OP;
220:         svd->AT = PETSC_NULL;
221:       } else {
222:         svd->A = PETSC_NULL;
223:         svd->AT = svd->OP;
224:       }
225:       break;
226:     default:
227:       SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid transpose mode");
228:   }

230:   /* set initial vector */
231:   if (!svd->vec_initial) {
232:     SVDMatGetVecs(svd,&svd->vec_initial,PETSC_NULL);
233:     SlepcVecSetRandom(svd->vec_initial);
234:   }

236:   /* call specific solver setup */
237:   (*svd->ops->setup)(svd);

239:   if (svd->ncv > M || svd->ncv > N)
240:     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"ncv bigger than matrix dimensions");
241:   if (svd->nsv > svd->ncv)
242:     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"nsv bigger than ncv");

244:   if (svd->ncv != svd->n) {
245:     /* free memory for previous solution  */
246:     if (svd->n) {
247:       PetscFree(svd->sigma);
248:       PetscFree(svd->perm);
249:       PetscFree(svd->errest);
250:       VecGetArray(svd->V[0],&pV);
251:       for (i=0;i<svd->n;i++) {
252:         VecDestroy(svd->V[i]);
253:       }
254:       PetscFree(pV);
255:       PetscFree(svd->V);
256:     }
257:     /* allocate memory for next solution */
258:     PetscMalloc(svd->ncv*sizeof(PetscReal),&svd->sigma);
259:     PetscMalloc(svd->ncv*sizeof(PetscReal),&svd->perm);
260:     PetscMalloc(svd->ncv*sizeof(PetscReal),&svd->errest);
261:     PetscMalloc(svd->ncv*sizeof(Vec),&svd->V);
262:     VecGetLocalSize(svd->vec_initial,&nloc);
263:     PetscMalloc(svd->ncv*nloc*sizeof(PetscScalar),&pV);
264:     for (i=0;i<svd->ncv;i++) {
265:       VecCreateMPIWithArray(((PetscObject)svd)->comm,nloc,PETSC_DECIDE,pV+i*nloc,&svd->V[i]);
266:     }
267:     svd->n = svd->ncv;
268:   }

270:   PetscLogEventEnd(SVD_SetUp,svd,0,0,0);
271:   svd->setupcalled = 1;
272:   return(0);
273: }