Actual source code: svdsetup.c

slepc-3.6.1 2015-09-03
Report Typos and Errors
  1: /*
  2:      SVD routines for setting up the solver.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

 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 <slepc/private/svdimpl.h>      /*I "slepcsvd.h" I*/

 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:   if (svd->state) { SVDReset(svd); }
 50:   PetscObjectReference((PetscObject)mat);
 51:   MatDestroy(&svd->OP);
 52:   svd->OP = mat;
 53:   return(0);
 54: }

 58: /*@
 59:    SVDGetOperator - Get the matrix associated with the singular value problem.

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

 63:    Input Parameter:
 64: .  svd - the singular value solver context

 66:    Output Parameters:
 67: .  A    - the matrix associated with the singular value problem

 69:    Level: advanced

 71: .seealso: SVDSolve(), SVDSetOperator()
 72: @*/
 73: PetscErrorCode SVDGetOperator(SVD svd,Mat *A)
 74: {
 78:   *A = svd->OP;
 79:   return(0);
 80: }

 84: /*@
 85:    SVDSetUp - Sets up all the internal data structures necessary for the
 86:    execution of the singular value solver.

 88:    Collective on SVD

 90:    Input Parameter:
 91: .  svd   - singular value solver context

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

 98:    Level: developer

100: .seealso: SVDCreate(), SVDSolve(), SVDDestroy()
101: @*/
102: PetscErrorCode SVDSetUp(SVD svd)
103: {
105:   PetscBool      expltrans,flg;
106:   PetscInt       M,N,k;
107:   SlepcSC        sc;
108:   Vec            *T;

112:   if (svd->state) return(0);
113:   PetscLogEventBegin(SVD_SetUp,svd,0,0,0);

115:   /* reset the convergence flag from the previous solves */
116:   svd->reason = SVD_CONVERGED_ITERATING;

118:   /* Set default solver type (SVDSetFromOptions was not called) */
119:   if (!((PetscObject)svd)->type_name) {
120:     SVDSetType(svd,SVDCROSS);
121:   }
122:   if (!svd->ds) { SVDGetDS(svd,&svd->ds); }
123:   DSReset(svd->ds);
124:   if (!((PetscObject)svd->rand)->type_name) {
125:     PetscRandomSetFromOptions(svd->rand);
126:   }

128:   /* check matrix */
129:   if (!svd->OP) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_WRONGSTATE,"SVDSetOperator must be called first");

131:   /* determine how to handle the transpose */
132:   expltrans = PETSC_TRUE;
133:   if (svd->impltrans) expltrans = PETSC_FALSE;
134:   else {
135:     MatHasOperation(svd->OP,MATOP_TRANSPOSE,&flg);
136:     if (!flg) expltrans = PETSC_FALSE;
137:     else {
138:       PetscObjectTypeCompare((PetscObject)svd,SVDLAPACK,&flg);
139:       if (flg) expltrans = PETSC_FALSE;
140:     }
141:   }

143:   /* build transpose matrix */
144:   MatDestroy(&svd->A);
145:   MatDestroy(&svd->AT);
146:   MatGetSize(svd->OP,&M,&N);
147:   PetscObjectReference((PetscObject)svd->OP);
148:   if (expltrans) {
149:     if (M>=N) {
150:       svd->A = svd->OP;
151:       MatTranspose(svd->OP,MAT_INITIAL_MATRIX,&svd->AT);
152:       MatConjugate(svd->AT);
153:     } else {
154:       MatTranspose(svd->OP,MAT_INITIAL_MATRIX,&svd->A);
155:       MatConjugate(svd->A);
156:       svd->AT = svd->OP;
157:     }
158:   } else {
159:     if (M>=N) {
160:       svd->A = svd->OP;
161:       svd->AT = NULL;
162:     } else {
163:       svd->A = NULL;
164:       svd->AT = svd->OP;
165:     }
166:   }

168:   /* swap initial vectors if necessary */
169:   if (M<N) {
170:     T=svd->ISL; svd->ISL=svd->IS; svd->IS=T;
171:     k=svd->ninil; svd->ninil=svd->nini; svd->nini=k;
172:   }

174:   if (svd->ncv > PetscMin(M,N)) svd->ncv = PetscMin(M,N);
175:   if (svd->nsv > PetscMin(M,N)) svd->nsv = PetscMin(M,N);
176:   if (svd->ncv && svd->nsv > svd->ncv) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"nsv bigger than ncv");

178:   /* call specific solver setup */
179:   (*svd->ops->setup)(svd);

181:   /* set tolerance if not yet set */
182:   if (svd->tol==PETSC_DEFAULT) svd->tol = SLEPC_DEFAULT_TOL;

184:   /* fill sorting criterion context */
185:   DSGetSlepcSC(svd->ds,&sc);
186:   sc->comparison    = (svd->which==SVD_LARGEST)? SlepcCompareLargestReal: SlepcCompareSmallestReal;
187:   sc->comparisonctx = NULL;
188:   sc->map           = NULL;
189:   sc->mapobj        = NULL;

191:   /* process initial vectors */
192:   if (svd->nini<0) {
193:     k = -svd->nini;
194:     if (k>svd->ncv) SETERRQ(PetscObjectComm((PetscObject)svd),1,"The number of initial vectors is larger than ncv");
195:     BVInsertVecs(svd->V,0,&k,svd->IS,PETSC_TRUE);
196:     SlepcBasisDestroy_Private(&svd->nini,&svd->IS);
197:     svd->nini = k;
198:   }
199:   if (svd->ninil<0) {
200:     k = 0;
201:     if (svd->leftbasis) {
202:       k = -svd->ninil;
203:       if (k>svd->ncv) SETERRQ(PetscObjectComm((PetscObject)svd),1,"The number of left initial vectors is larger than ncv");
204:       BVInsertVecs(svd->U,0,&k,svd->ISL,PETSC_TRUE);
205:     } else {
206:       PetscInfo(svd,"Ignoring initial left vectors\n");
207:     }
208:     SlepcBasisDestroy_Private(&svd->ninil,&svd->ISL);
209:     svd->ninil = k;
210:   }

212:   PetscLogEventEnd(SVD_SetUp,svd,0,0,0);
213:   svd->state = SVD_STATE_SETUP;
214:   return(0);
215: }

219: /*@
220:    SVDSetInitialSpace - Specify a basis of vectors that constitute the initial
221:    (right) space, that is, a rough approximation to the right singular subspace
222:    from which the solver starts to iterate.

224:    Collective on SVD and Vec

226:    Input Parameter:
227: +  svd   - the singular value solver context
228: .  n     - number of vectors
229: -  is    - set of basis vectors of the initial space

231:    Notes:
232:    Some solvers start to iterate on a single vector (initial vector). In that case,
233:    the other vectors are ignored.

235:    These vectors do not persist from one SVDSolve() call to the other, so the
236:    initial space should be set every time.

238:    The vectors do not need to be mutually orthonormal, since they are explicitly
239:    orthonormalized internally.

241:    Common usage of this function is when the user can provide a rough approximation
242:    of the wanted singular space. Then, convergence may be faster.

244:    Level: intermediate

246: .seealso: SVDSetInitialSpaceLeft()
247: @*/
248: PetscErrorCode SVDSetInitialSpace(SVD svd,PetscInt n,Vec *is)
249: {

255:   if (n<0) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
256:   SlepcBasisReference_Private(n,is,&svd->nini,&svd->IS);
257:   if (n>0) svd->state = SVD_STATE_INITIAL;
258:   return(0);
259: }

263: /*@
264:    SVDSetInitialSpaceLeft - Specify a basis of vectors that constitute the initial
265:    left space, that is, a rough approximation to the left singular subspace
266:    from which the solver starts to iterate.

268:    Collective on SVD and Vec

270:    Input Parameter:
271: +  svd   - the singular value solver context
272: .  n     - number of vectors
273: -  is    - set of basis vectors of the initial space

275:    Notes:
276:    Some solvers start to iterate on a single vector (initial vector). In that case,
277:    the other vectors are ignored.

279:    These vectors do not persist from one SVDSolve() call to the other, so the
280:    initial space should be set every time.

282:    The vectors do not need to be mutually orthonormal, since they are explicitly
283:    orthonormalized internally.

285:    Common usage of this function is when the user can provide a rough approximation
286:    of the wanted singular space. Then, convergence may be faster.

288:    Level: intermediate

290: .seealso: SVDSetInitialSpace()
291: @*/
292: PetscErrorCode SVDSetInitialSpaceLeft(SVD svd,PetscInt n,Vec *is)
293: {

299:   if (n<0) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
300:   SlepcBasisReference_Private(n,is,&svd->ninil,&svd->ISL);
301:   if (n>0) svd->state = SVD_STATE_INITIAL;
302:   return(0);
303: }

307: /*
308:   SVDSetDimensions_Default - Set reasonable values for ncv, mpd if not set
309:   by the user. This is called at setup.
310:  */
311: PetscErrorCode SVDSetDimensions_Default(SVD svd)
312: {
314:   PetscInt       N;

317:   SVDMatGetSize(svd,NULL,&N);
318:   if (svd->ncv) { /* ncv set */
319:     if (svd->ncv<svd->nsv) SETERRQ(PetscObjectComm((PetscObject)svd),1,"The value of ncv must be at least nsv");
320:   } else if (svd->mpd) { /* mpd set */
321:     svd->ncv = PetscMin(N,svd->nsv+svd->mpd);
322:   } else { /* neither set: defaults depend on nsv being small or large */
323:     if (svd->nsv<500) svd->ncv = PetscMin(N,PetscMax(2*svd->nsv,10));
324:     else {
325:       svd->mpd = 500;
326:       svd->ncv = PetscMin(N,svd->nsv+svd->mpd);
327:     }
328:   }
329:   if (!svd->mpd) svd->mpd = svd->ncv;
330:   return(0);
331: }

335: /*@
336:    SVDAllocateSolution - Allocate memory storage for common variables such
337:    as the singular values and the basis vectors.

339:    Collective on SVD

341:    Input Parameters:
342: +  svd   - eigensolver context
343: -  extra - number of additional positions, used for methods that require a
344:            working basis slightly larger than ncv

346:    Developers Notes:
347:    This is PETSC_EXTERN because it may be required by user plugin SVD
348:    implementations.

350:    This is called at setup after setting the value of ncv and the flag leftbasis.

352:    Level: developer
353: @*/
354: PetscErrorCode SVDAllocateSolution(SVD svd,PetscInt extra)
355: {
357:   PetscInt       oldsize,requested;
358:   Vec            tr,tl;

361:   requested = svd->ncv + extra;

363:   /* oldsize is zero if this is the first time setup is called */
364:   BVGetSizes(svd->V,NULL,NULL,&oldsize);

366:   /* allocate sigma */
367:   if (requested != oldsize || !svd->sigma) {
368:     if (oldsize) {
369:       PetscFree3(svd->sigma,svd->perm,svd->errest);
370:     }
371:     PetscMalloc3(requested,&svd->sigma,requested,&svd->perm,requested,&svd->errest);
372:     PetscLogObjectMemory((PetscObject)svd,PetscMax(0,requested-oldsize)*(2*sizeof(PetscReal)+sizeof(PetscInt)));
373:   }
374:   /* allocate V */
375:   if (!svd->V) { SVDGetBV(svd,&svd->V,NULL); }
376:   if (!oldsize) {
377:     if (!((PetscObject)(svd->V))->type_name) {
378:       BVSetType(svd->V,BVSVEC);
379:     }
380:     SVDMatCreateVecs(svd,&tr,NULL);
381:     BVSetSizesFromVec(svd->V,tr,requested);
382:     VecDestroy(&tr);
383:   } else {
384:     BVResize(svd->V,requested,PETSC_FALSE);
385:   }
386:   /* allocate U */
387:   if (svd->leftbasis) {
388:     if (!svd->U) { SVDGetBV(svd,NULL,&svd->U); }
389:     if (!oldsize) {
390:       if (!((PetscObject)(svd->U))->type_name) {
391:         BVSetType(svd->U,BVSVEC);
392:       }
393:       SVDMatCreateVecs(svd,NULL,&tl);
394:       BVSetSizesFromVec(svd->U,tl,requested);
395:       VecDestroy(&tl);
396:     } else {
397:       BVResize(svd->U,requested,PETSC_FALSE);
398:     }
399:   }
400:   return(0);
401: }