Actual source code: svdbasic.c

slepc-3.15.2 2021-09-20
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    Basic SVD routines
 12: */

 14: #include <slepc/private/svdimpl.h>

 16: /* Logging support */
 17: PetscClassId      SVD_CLASSID = 0;
 18: PetscLogEvent     SVD_SetUp = 0,SVD_Solve = 0;

 20: /* List of registered SVD routines */
 21: PetscFunctionList SVDList = 0;
 22: PetscBool         SVDRegisterAllCalled = PETSC_FALSE;

 24: /* List of registered SVD monitors */
 25: PetscFunctionList SVDMonitorList              = NULL;
 26: PetscFunctionList SVDMonitorCreateList        = NULL;
 27: PetscFunctionList SVDMonitorDestroyList       = NULL;
 28: PetscBool         SVDMonitorRegisterAllCalled = PETSC_FALSE;

 30: /*@
 31:    SVDCreate - Creates the default SVD context.

 33:    Collective

 35:    Input Parameter:
 36: .  comm - MPI communicator

 38:    Output Parameter:
 39: .  svd - location to put the SVD context

 41:    Note:
 42:    The default SVD type is SVDCROSS

 44:    Level: beginner

 46: .seealso: SVDSetUp(), SVDSolve(), SVDDestroy(), SVD
 47: @*/
 48: PetscErrorCode SVDCreate(MPI_Comm comm,SVD *outsvd)
 49: {
 51:   SVD            svd;

 55:   *outsvd = 0;
 56:   SVDInitializePackage();
 57:   SlepcHeaderCreate(svd,SVD_CLASSID,"SVD","Singular Value Decomposition","SVD",comm,SVDDestroy,SVDView);

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

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

 85:   svd->ds               = NULL;
 86:   svd->U                = NULL;
 87:   svd->V                = NULL;
 88:   svd->A                = NULL;
 89:   svd->B                = NULL;
 90:   svd->AT               = NULL;
 91:   svd->BT               = NULL;
 92:   svd->IS               = NULL;
 93:   svd->ISL              = NULL;
 94:   svd->sigma            = NULL;
 95:   svd->errest           = NULL;
 96:   svd->perm             = NULL;
 97:   svd->nworkl           = 0;
 98:   svd->nworkr           = 0;
 99:   svd->workl            = NULL;
100:   svd->workr            = NULL;
101:   svd->data             = NULL;

103:   svd->state            = SVD_STATE_INITIAL;
104:   svd->nconv            = 0;
105:   svd->its              = 0;
106:   svd->leftbasis        = PETSC_FALSE;
107:   svd->swapped          = PETSC_FALSE;
108:   svd->expltrans        = PETSC_FALSE;
109:   svd->isgeneralized    = PETSC_FALSE;
110:   svd->reason           = SVD_CONVERGED_ITERATING;

112:   PetscNewLog(svd,&svd->sc);
113:   *outsvd = svd;
114:   return(0);
115: }

117: /*@
118:    SVDReset - Resets the SVD context to the initial state (prior to setup)
119:    and destroys any allocated Vecs and Mats.

121:    Collective on svd

123:    Input Parameter:
124: .  svd - singular value solver context obtained from SVDCreate()

126:    Level: advanced

128: .seealso: SVDDestroy()
129: @*/
130: PetscErrorCode SVDReset(SVD svd)
131: {

136:   if (!svd) return(0);
137:   if (svd->ops->reset) { (svd->ops->reset)(svd); }
138:   MatDestroy(&svd->OP);
139:   MatDestroy(&svd->OPb);
140:   MatDestroy(&svd->A);
141:   MatDestroy(&svd->B);
142:   MatDestroy(&svd->AT);
143:   MatDestroy(&svd->BT);
144:   BVDestroy(&svd->U);
145:   BVDestroy(&svd->V);
146:   VecDestroyVecs(svd->nworkl,&svd->workl);
147:   svd->nworkl = 0;
148:   VecDestroyVecs(svd->nworkr,&svd->workr);
149:   svd->nworkr = 0;
150:   svd->state = SVD_STATE_INITIAL;
151:   return(0);
152: }

154: /*@C
155:    SVDDestroy - Destroys the SVD context.

157:    Collective on svd

159:    Input Parameter:
160: .  svd - singular value solver context obtained from SVDCreate()

162:    Level: beginner

164: .seealso: SVDCreate(), SVDSetUp(), SVDSolve()
165: @*/
166: PetscErrorCode SVDDestroy(SVD *svd)
167: {

171:   if (!*svd) return(0);
173:   if (--((PetscObject)(*svd))->refct > 0) { *svd = 0; return(0); }
174:   SVDReset(*svd);
175:   if ((*svd)->ops->destroy) { (*(*svd)->ops->destroy)(*svd); }
176:   if ((*svd)->sigma) {
177:     PetscFree3((*svd)->sigma,(*svd)->perm,(*svd)->errest);
178:   }
179:   DSDestroy(&(*svd)->ds);
180:   PetscFree((*svd)->sc);
181:   /* just in case the initial vectors have not been used */
182:   SlepcBasisDestroy_Private(&(*svd)->nini,&(*svd)->IS);
183:   SlepcBasisDestroy_Private(&(*svd)->ninil,&(*svd)->ISL);
184:   SVDMonitorCancel(*svd);
185:   PetscHeaderDestroy(svd);
186:   return(0);
187: }

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

192:    Logically Collective on svd

194:    Input Parameters:
195: +  svd      - the singular value solver context
196: -  type     - a known method

198:    Options Database Key:
199: .  -svd_type <method> - Sets the method; use -help for a list
200:     of available methods

202:    Notes:
203:    See "slepc/include/slepcsvd.h" for available methods. The default
204:    is SVDCROSS.

206:    Normally, it is best to use the SVDSetFromOptions() command and
207:    then set the SVD type from the options database rather than by using
208:    this routine.  Using the options database provides the user with
209:    maximum flexibility in evaluating the different available methods.
210:    The SVDSetType() routine is provided for those situations where it
211:    is necessary to set the iterative solver independently of the command
212:    line or options database.

214:    Level: intermediate

216: .seealso: SVDType
217: @*/
218: PetscErrorCode SVDSetType(SVD svd,SVDType type)
219: {
220:   PetscErrorCode ierr,(*r)(SVD);
221:   PetscBool      match;


227:   PetscObjectTypeCompare((PetscObject)svd,type,&match);
228:   if (match) return(0);

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

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

236:   svd->state = SVD_STATE_INITIAL;
237:   PetscObjectChangeTypeName((PetscObject)svd,type);
238:   (*r)(svd);
239:   return(0);
240: }

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

245:    Not Collective

247:    Input Parameter:
248: .  svd - the singular value solver context

250:    Output Parameter:
251: .  name - name of SVD method

253:    Level: intermediate

255: .seealso: SVDSetType()
256: @*/
257: PetscErrorCode SVDGetType(SVD svd,SVDType *type)
258: {
262:   *type = ((PetscObject)svd)->type_name;
263:   return(0);
264: }

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

269:    Not Collective

271:    Input Parameters:
272: +  name - name of a new user-defined solver
273: -  function - routine to create the solver context

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

278:    Sample usage:
279: .vb
280:     SVDRegister("my_solver",MySolverCreate);
281: .ve

283:    Then, your solver can be chosen with the procedural interface via
284: $     SVDSetType(svd,"my_solver")
285:    or at runtime via the option
286: $     -svd_type my_solver

288:    Level: advanced

290: .seealso: SVDRegisterAll()
291: @*/
292: PetscErrorCode SVDRegister(const char *name,PetscErrorCode (*function)(SVD))
293: {

297:   SVDInitializePackage();
298:   PetscFunctionListAdd(&SVDList,name,function);
299:   return(0);
300: }

302: /*@C
303:    SVDMonitorRegister - Adds SVD monitor routine.

305:    Not Collective

307:    Input Parameters:
308: +  name    - name of a new monitor routine
309: .  vtype   - a PetscViewerType for the output
310: .  format  - a PetscViewerFormat for the output
311: .  monitor - monitor routine
312: .  create  - creation routine, or NULL
313: -  destroy - destruction routine, or NULL

315:    Notes:
316:    SVDMonitorRegister() may be called multiple times to add several user-defined monitors.

318:    Sample usage:
319: .vb
320:    SVDMonitorRegister("my_monitor",PETSCVIEWERASCII,PETSC_VIEWER_ASCII_INFO_DETAIL,MyMonitor,NULL,NULL);
321: .ve

323:    Then, your monitor can be chosen with the procedural interface via
324: $      SVDMonitorSetFromOptions(svd,"-svd_monitor_my_monitor","my_monitor",NULL)
325:    or at runtime via the option
326: $      -svd_monitor_my_monitor

328:    Level: advanced

330: .seealso: SVDMonitorRegisterAll()
331: @*/
332: PetscErrorCode SVDMonitorRegister(const char name[],PetscViewerType vtype,PetscViewerFormat format,PetscErrorCode (*monitor)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,PetscViewerAndFormat*),PetscErrorCode (*create)(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**),PetscErrorCode (*destroy)(PetscViewerAndFormat**))
333: {
334:   char           key[PETSC_MAX_PATH_LEN];

338:   SVDInitializePackage();
339:   SlepcMonitorMakeKey_Internal(name,vtype,format,key);
340:   PetscFunctionListAdd(&SVDMonitorList,key,monitor);
341:   if (create)  { PetscFunctionListAdd(&SVDMonitorCreateList,key,create); }
342:   if (destroy) { PetscFunctionListAdd(&SVDMonitorDestroyList,key,destroy); }
343:   return(0);
344: }

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

349:    Collective on svd

351:    Input Parameters:
352: +  svd - singular value solver context obtained from SVDCreate()
353: .  V   - the basis vectors object for right singular vectors
354: -  U   - the basis vectors object for left singular vectors

356:    Note:
357:    Use SVDGetBV() to retrieve the basis vectors contexts (for example,
358:    to free them at the end of the computations).

360:    Level: advanced

362: .seealso: SVDGetBV()
363: @*/
364: PetscErrorCode SVDSetBV(SVD svd,BV V,BV U)
365: {

370:   if (V) {
373:     PetscObjectReference((PetscObject)V);
374:     BVDestroy(&svd->V);
375:     svd->V = V;
376:     PetscLogObjectParent((PetscObject)svd,(PetscObject)svd->V);
377:   }
378:   if (U) {
381:     PetscObjectReference((PetscObject)U);
382:     BVDestroy(&svd->U);
383:     svd->U = U;
384:     PetscLogObjectParent((PetscObject)svd,(PetscObject)svd->U);
385:   }
386:   return(0);
387: }

389: /*@
390:    SVDGetBV - Obtain the basis vectors objects associated to the singular
391:    value solver object.

393:    Not Collective

395:    Input Parameters:
396: .  svd - singular value solver context obtained from SVDCreate()

398:    Output Parameter:
399: +  V - basis vectors context for right singular vectors
400: -  U - basis vectors context for left singular vectors

402:    Level: advanced

404: .seealso: SVDSetBV()
405: @*/
406: PetscErrorCode SVDGetBV(SVD svd,BV *V,BV *U)
407: {

412:   if (V) {
413:     if (!svd->V) {
414:       BVCreate(PetscObjectComm((PetscObject)svd),&svd->V);
415:       PetscObjectIncrementTabLevel((PetscObject)svd->V,(PetscObject)svd,0);
416:       PetscLogObjectParent((PetscObject)svd,(PetscObject)svd->V);
417:       PetscObjectSetOptions((PetscObject)svd->V,((PetscObject)svd)->options);
418:     }
419:     *V = svd->V;
420:   }
421:   if (U) {
422:     if (!svd->U) {
423:       BVCreate(PetscObjectComm((PetscObject)svd),&svd->U);
424:       PetscObjectIncrementTabLevel((PetscObject)svd->U,(PetscObject)svd,0);
425:       PetscLogObjectParent((PetscObject)svd,(PetscObject)svd->U);
426:       PetscObjectSetOptions((PetscObject)svd->U,((PetscObject)svd)->options);
427:     }
428:     *U = svd->U;
429:   }
430:   return(0);
431: }

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

436:    Collective on svd

438:    Input Parameters:
439: +  svd - singular value solver context obtained from SVDCreate()
440: -  ds  - the direct solver object

442:    Note:
443:    Use SVDGetDS() to retrieve the direct solver context (for example,
444:    to free it at the end of the computations).

446:    Level: advanced

448: .seealso: SVDGetDS()
449: @*/
450: PetscErrorCode SVDSetDS(SVD svd,DS ds)
451: {

458:   PetscObjectReference((PetscObject)ds);
459:   DSDestroy(&svd->ds);
460:   svd->ds = ds;
461:   PetscLogObjectParent((PetscObject)svd,(PetscObject)svd->ds);
462:   return(0);
463: }

465: /*@
466:    SVDGetDS - Obtain the direct solver object associated to the singular value
467:    solver object.

469:    Not Collective

471:    Input Parameters:
472: .  svd - singular value solver context obtained from SVDCreate()

474:    Output Parameter:
475: .  ds - direct solver context

477:    Level: advanced

479: .seealso: SVDSetDS()
480: @*/
481: PetscErrorCode SVDGetDS(SVD svd,DS *ds)
482: {

488:   if (!svd->ds) {
489:     DSCreate(PetscObjectComm((PetscObject)svd),&svd->ds);
490:     PetscObjectIncrementTabLevel((PetscObject)svd->ds,(PetscObject)svd,0);
491:     PetscLogObjectParent((PetscObject)svd,(PetscObject)svd->ds);
492:     PetscObjectSetOptions((PetscObject)svd->ds,((PetscObject)svd)->options);
493:   }
494:   *ds = svd->ds;
495:   return(0);
496: }