Actual source code: svdbasic.c

slepc-3.7.2 2016-07-19
Report Typos and Errors
  1: /*
  2:    The basic SVD routines, Create, Destroy, etc. are here.

  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*/

 26: PetscFunctionList SVDList = 0;
 27: PetscBool         SVDRegisterAllCalled = PETSC_FALSE;
 28: PetscClassId      SVD_CLASSID = 0;
 29: PetscLogEvent     SVD_SetUp = 0,SVD_Solve = 0;

 33: /*@
 34:    SVDCreate - Creates the default SVD context.

 36:    Collective on MPI_Comm

 38:    Input Parameter:
 39: .  comm - MPI communicator

 41:    Output Parameter:
 42: .  svd - location to put the SVD context

 44:    Note:
 45:    The default SVD type is SVDCROSS

 47:    Level: beginner

 49: .seealso: SVDSetUp(), SVDSolve(), SVDDestroy(), SVD
 50: @*/
 51: PetscErrorCode SVDCreate(MPI_Comm comm,SVD *outsvd)
 52: {
 54:   SVD            svd;

 58:   *outsvd = 0;
 59:   SVDInitializePackage();
 60:   SlepcHeaderCreate(svd,SVD_CLASSID,"SVD","Singular Value Decomposition","SVD",comm,SVDDestroy,SVDView);

 62:   svd->OP               = NULL;
 63:   svd->max_it           = 0;
 64:   svd->nsv              = 1;
 65:   svd->ncv              = 0;
 66:   svd->mpd              = 0;
 67:   svd->nini             = 0;
 68:   svd->ninil            = 0;
 69:   svd->tol              = PETSC_DEFAULT;
 70:   svd->conv             = SVD_CONV_REL;
 71:   svd->stop             = SVD_STOP_BASIC;
 72:   svd->which            = SVD_LARGEST;
 73:   svd->impltrans        = PETSC_FALSE;
 74:   svd->trackall         = PETSC_FALSE;

 76:   svd->converged        = SVDConvergedRelative;
 77:   svd->convergeddestroy = NULL;
 78:   svd->stopping         = SVDStoppingBasic;
 79:   svd->stoppingdestroy  = NULL;
 80:   svd->convergedctx     = NULL;
 81:   svd->stoppingctx      = NULL;
 82:   svd->numbermonitors   = 0;

 84:   svd->ds               = NULL;
 85:   svd->U                = NULL;
 86:   svd->V                = NULL;
 87:   svd->A                = NULL;
 88:   svd->AT               = NULL;
 89:   svd->IS               = NULL;
 90:   svd->ISL              = NULL;
 91:   svd->sigma            = NULL;
 92:   svd->perm             = NULL;
 93:   svd->errest           = NULL;
 94:   svd->data             = NULL;

 96:   svd->state            = SVD_STATE_INITIAL;
 97:   svd->nconv            = 0;
 98:   svd->its              = 0;
 99:   svd->leftbasis        = PETSC_FALSE;
100:   svd->reason           = SVD_CONVERGED_ITERATING;

102:   PetscNewLog(svd,&svd->sc);
103:   *outsvd = svd;
104:   return(0);
105: }

109: /*@
110:    SVDReset - Resets the SVD context to the initial state and removes any
111:    allocated objects.

113:    Collective on SVD

115:    Input Parameter:
116: .  svd - singular value solver context obtained from SVDCreate()

118:    Level: advanced

120: .seealso: SVDDestroy()
121: @*/
122: PetscErrorCode SVDReset(SVD svd)
123: {
125:   PetscInt       ncols;

129:   if (svd->ops->reset) { (svd->ops->reset)(svd); }
130:   if (svd->ds) { DSReset(svd->ds); }
131:   MatDestroy(&svd->OP);
132:   MatDestroy(&svd->A);
133:   MatDestroy(&svd->AT);
134:   BVGetSizes(svd->V,NULL,NULL,&ncols);
135:   if (ncols) {
136:     PetscFree3(svd->sigma,svd->perm,svd->errest);
137:   }
138:   BVDestroy(&svd->U);
139:   BVDestroy(&svd->V);
140:   svd->state = SVD_STATE_INITIAL;
141:   return(0);
142: }

146: /*@
147:    SVDDestroy - Destroys the SVD context.

149:    Collective on SVD

151:    Input Parameter:
152: .  svd - singular value solver context obtained from SVDCreate()

154:    Level: beginner

156: .seealso: SVDCreate(), SVDSetUp(), SVDSolve()
157: @*/
158: PetscErrorCode SVDDestroy(SVD *svd)
159: {

163:   if (!*svd) return(0);
165:   if (--((PetscObject)(*svd))->refct > 0) { *svd = 0; return(0); }
166:   SVDReset(*svd);
167:   if ((*svd)->ops->destroy) { (*(*svd)->ops->destroy)(*svd); }
168:   DSDestroy(&(*svd)->ds);
169:   PetscFree((*svd)->sc);
170:   /* just in case the initial vectors have not been used */
171:   SlepcBasisDestroy_Private(&(*svd)->nini,&(*svd)->IS);
172:   SlepcBasisDestroy_Private(&(*svd)->ninil,&(*svd)->ISL);
173:   SVDMonitorCancel(*svd);
174:   PetscHeaderDestroy(svd);
175:   return(0);
176: }

180: /*@C
181:    SVDSetType - Selects the particular solver to be used in the SVD object.

183:    Logically Collective on SVD

185:    Input Parameters:
186: +  svd      - the singular value solver context
187: -  type     - a known method

189:    Options Database Key:
190: .  -svd_type <method> - Sets the method; use -help for a list
191:     of available methods

193:    Notes:
194:    See "slepc/include/slepcsvd.h" for available methods. The default
195:    is SVDCROSS.

197:    Normally, it is best to use the SVDSetFromOptions() command and
198:    then set the SVD type from the options database rather than by using
199:    this routine.  Using the options database provides the user with
200:    maximum flexibility in evaluating the different available methods.
201:    The SVDSetType() routine is provided for those situations where it
202:    is necessary to set the iterative solver independently of the command
203:    line or options database.

205:    Level: intermediate

207: .seealso: SVDType
208: @*/
209: PetscErrorCode SVDSetType(SVD svd,SVDType type)
210: {
211:   PetscErrorCode ierr,(*r)(SVD);
212:   PetscBool      match;


218:   PetscObjectTypeCompare((PetscObject)svd,type,&match);
219:   if (match) return(0);

221:   PetscFunctionListFind(SVDList,type,&r);
222:   if (!r) SETERRQ1(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown SVD type given: %s",type);

224:   if (svd->ops->destroy) { (*svd->ops->destroy)(svd); }
225:   PetscMemzero(svd->ops,sizeof(struct _SVDOps));

227:   svd->state = SVD_STATE_INITIAL;
228:   PetscObjectChangeTypeName((PetscObject)svd,type);
229:   (*r)(svd);
230:   return(0);
231: }

235: /*@C
236:    SVDGetType - Gets the SVD type as a string from the SVD object.

238:    Not Collective

240:    Input Parameter:
241: .  svd - the singular value solver context

243:    Output Parameter:
244: .  name - name of SVD method

246:    Level: intermediate

248: .seealso: SVDSetType()
249: @*/
250: PetscErrorCode SVDGetType(SVD svd,SVDType *type)
251: {
255:   *type = ((PetscObject)svd)->type_name;
256:   return(0);
257: }

261: /*@C
262:    SVDRegister - Adds a method to the singular value solver package.

264:    Not Collective

266:    Input Parameters:
267: +  name - name of a new user-defined solver
268: -  function - routine to create the solver context

270:    Notes:
271:    SVDRegister() may be called multiple times to add several user-defined solvers.

273:    Sample usage:
274: .vb
275:     SVDRegister("my_solver",MySolverCreate);
276: .ve

278:    Then, your solver can be chosen with the procedural interface via
279: $     SVDSetType(svd,"my_solver")
280:    or at runtime via the option
281: $     -svd_type my_solver

283:    Level: advanced

285: .seealso: SVDRegisterAll()
286: @*/
287: PetscErrorCode SVDRegister(const char *name,PetscErrorCode (*function)(SVD))
288: {

292:   PetscFunctionListAdd(&SVDList,name,function);
293:   return(0);
294: }

298: /*@
299:    SVDSetBV - Associates basis vectors objects to the singular value solver.

301:    Collective on SVD

303:    Input Parameters:
304: +  svd - singular value solver context obtained from SVDCreate()
305: .  V   - the basis vectors object for right singular vectors
306: -  U   - the basis vectors object for left singular vectors

308:    Note:
309:    Use SVDGetBV() to retrieve the basis vectors contexts (for example,
310:    to free them at the end of the computations).

312:    Level: advanced

314: .seealso: SVDGetBV()
315: @*/
316: PetscErrorCode SVDSetBV(SVD svd,BV V,BV U)
317: {

322:   if (V) {
325:     PetscObjectReference((PetscObject)V);
326:     BVDestroy(&svd->V);
327:     svd->V = V;
328:     PetscLogObjectParent((PetscObject)svd,(PetscObject)svd->V);
329:   }
330:   if (U) {
333:     PetscObjectReference((PetscObject)U);
334:     BVDestroy(&svd->U);
335:     svd->U = U;
336:     PetscLogObjectParent((PetscObject)svd,(PetscObject)svd->U);
337:   }
338:   return(0);
339: }

343: /*@
344:    SVDGetBV - Obtain the basis vectors objects associated to the singular
345:    value solver object.

347:    Not Collective

349:    Input Parameters:
350: .  svd - singular value solver context obtained from SVDCreate()

352:    Output Parameter:
353: +  V - basis vectors context for right singular vectors
354: -  U - basis vectors context for left singular vectors

356:    Level: advanced

358: .seealso: SVDSetBV()
359: @*/
360: PetscErrorCode SVDGetBV(SVD svd,BV *V,BV *U)
361: {

366:   if (V) {
367:     if (!svd->V) {
368:       BVCreate(PetscObjectComm((PetscObject)svd),&svd->V);
369:       PetscLogObjectParent((PetscObject)svd,(PetscObject)svd->V);
370:     }
371:     *V = svd->V;
372:   }
373:   if (U) {
374:     if (!svd->U) {
375:       BVCreate(PetscObjectComm((PetscObject)svd),&svd->U);
376:       PetscLogObjectParent((PetscObject)svd,(PetscObject)svd->U);
377:     }
378:     *U = svd->U;
379:   }
380:   return(0);
381: }

385: /*@
386:    SVDSetDS - Associates a direct solver object to the singular value solver.

388:    Collective on SVD

390:    Input Parameters:
391: +  svd - singular value solver context obtained from SVDCreate()
392: -  ds  - the direct solver object

394:    Note:
395:    Use SVDGetDS() to retrieve the direct solver context (for example,
396:    to free it at the end of the computations).

398:    Level: advanced

400: .seealso: SVDGetDS()
401: @*/
402: PetscErrorCode SVDSetDS(SVD svd,DS ds)
403: {

410:   PetscObjectReference((PetscObject)ds);
411:   DSDestroy(&svd->ds);
412:   svd->ds = ds;
413:   PetscLogObjectParent((PetscObject)svd,(PetscObject)svd->ds);
414:   return(0);
415: }

419: /*@
420:    SVDGetDS - Obtain the direct solver object associated to the singular value
421:    solver object.

423:    Not Collective

425:    Input Parameters:
426: .  svd - singular value solver context obtained from SVDCreate()

428:    Output Parameter:
429: .  ds - direct solver context

431:    Level: advanced

433: .seealso: SVDSetDS()
434: @*/
435: PetscErrorCode SVDGetDS(SVD svd,DS *ds)
436: {

442:   if (!svd->ds) {
443:     DSCreate(PetscObjectComm((PetscObject)svd),&svd->ds);
444:     PetscLogObjectParent((PetscObject)svd,(PetscObject)svd->ds);
445:   }
446:   *ds = svd->ds;
447:   return(0);
448: }