Actual source code: svdsetup.c

slepc-3.7.3 2016-09-29
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-2016, 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);

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

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

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

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

171:   if (svd->ncv > PetscMin(M,N)) svd->ncv = PetscMin(M,N);
172:   if (svd->nsv > PetscMin(M,N)) svd->nsv = PetscMin(M,N);
173:   if (svd->ncv && svd->nsv > svd->ncv) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"nsv bigger than ncv");

175:   /* call specific solver setup */
176:   (*svd->ops->setup)(svd);

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

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

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

209:   PetscLogEventEnd(SVD_SetUp,svd,0,0,0);
210:   svd->state = SVD_STATE_SETUP;
211:   return(0);
212: }

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

221:    Collective on SVD and Vec

223:    Input Parameter:
224: +  svd   - the singular value solver context
225: .  n     - number of vectors
226: -  is    - set of basis vectors of the initial space

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

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

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

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

241:    Level: intermediate

243: .seealso: SVDSetInitialSpaceLeft()
244: @*/
245: PetscErrorCode SVDSetInitialSpace(SVD svd,PetscInt n,Vec *is)
246: {

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

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

265:    Collective on SVD and Vec

267:    Input Parameter:
268: +  svd   - the singular value solver context
269: .  n     - number of vectors
270: -  is    - set of basis vectors of the initial space

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

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

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

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

285:    Level: intermediate

287: .seealso: SVDSetInitialSpace()
288: @*/
289: PetscErrorCode SVDSetInitialSpaceLeft(SVD svd,PetscInt n,Vec *is)
290: {

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

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

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

332: /*@
333:    SVDAllocateSolution - Allocate memory storage for common variables such
334:    as the singular values and the basis vectors.

336:    Collective on SVD

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

343:    Developers Notes:
344:    This is PETSC_EXTERN because it may be required by user plugin SVD
345:    implementations.

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

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

358:   requested = svd->ncv + extra;

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

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