Actual source code: svdbasic.c

  1: /*
  2:      The basic SVD routines, Create, View, etc. are here.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2013, 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;
 30: static PetscBool  SVDPackageInitialized = PETSC_FALSE;

 34: /*@C
 35:    SVDFinalizePackage - This function destroys everything in the Slepc interface
 36:    to the SVD package. It is called from SlepcFinalize().

 38:    Level: developer

 40: .seealso: SlepcFinalize()
 41: @*/
 42: PetscErrorCode SVDFinalizePackage(void)
 43: {

 47:   PetscFunctionListDestroy(&SVDList);
 48:   SVDPackageInitialized = PETSC_FALSE;
 49:   SVDRegisterAllCalled  = PETSC_FALSE;
 50:   return(0);
 51: }

 55: /*@C
 56:    SVDInitializePackage - This function initializes everything in the SVD package. It is called
 57:    from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to SVDCreate()
 58:    when using static libraries.

 60:    Level: developer

 62: .seealso: SlepcInitialize()
 63: @*/
 64: PetscErrorCode SVDInitializePackage(void)
 65: {
 66:   char           logList[256];
 67:   char           *className;
 68:   PetscBool      opt;

 72:   if (SVDPackageInitialized) return(0);
 73:   SVDPackageInitialized = PETSC_TRUE;
 74:   /* Register Classes */
 75:   PetscClassIdRegister("Singular Value Solver",&SVD_CLASSID);
 76:   /* Register Constructors */
 77:   SVDRegisterAll();
 78:   /* Register Events */
 79:   PetscLogEventRegister("SVDSetUp",SVD_CLASSID,&SVD_SetUp);
 80:   PetscLogEventRegister("SVDSolve",SVD_CLASSID,&SVD_Solve);
 81:   /* Process info exclusions */
 82:   PetscOptionsGetString(NULL,"-info_exclude",logList,256,&opt);
 83:   if (opt) {
 84:     PetscStrstr(logList,"svd",&className);
 85:     if (className) {
 86:       PetscInfoDeactivateClass(SVD_CLASSID);
 87:     }
 88:   }
 89:   /* Process summary exclusions */
 90:   PetscOptionsGetString(NULL,"-log_summary_exclude",logList,256,&opt);
 91:   if (opt) {
 92:     PetscStrstr(logList,"svd",&className);
 93:     if (className) {
 94:       PetscLogEventDeactivateClass(SVD_CLASSID);
 95:     }
 96:   }
 97:   PetscRegisterFinalize(SVDFinalizePackage);
 98:   return(0);
 99: }

103: /*@C
104:    SVDView - Prints the SVD data structure.

106:    Collective on SVD

108:    Input Parameters:
109: +  svd - the singular value solver context
110: -  viewer - optional visualization context

112:    Options Database Key:
113: .  -svd_view -  Calls SVDView() at end of SVDSolve()

115:    Note:
116:    The available visualization contexts include
117: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
118: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
119:          output where only the first processor opens
120:          the file.  All other processors send their
121:          data to the first processor to print.

123:    The user can open an alternative visualization context with
124:    PetscViewerASCIIOpen() - output to a specified file.

126:    Level: beginner

128: .seealso: STView(), PetscViewerASCIIOpen()
129: @*/
130: PetscErrorCode SVDView(SVD svd,PetscViewer viewer)
131: {
133:   PetscBool      isascii,isshell;

137:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)svd));

141:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
142:   if (isascii) {
143:     PetscObjectPrintClassNamePrefixType((PetscObject)svd,viewer,"SVD Object");
144:     if (svd->ops->view) {
145:       PetscViewerASCIIPushTab(viewer);
146:       (*svd->ops->view)(svd,viewer);
147:       PetscViewerASCIIPopTab(viewer);
148:     }
149:     switch (svd->transmode) {
150:       case SVD_TRANSPOSE_EXPLICIT:
151:         PetscViewerASCIIPrintf(viewer,"  transpose mode: explicit\n");
152:         break;
153:       case SVD_TRANSPOSE_IMPLICIT:
154:         PetscViewerASCIIPrintf(viewer,"  transpose mode: implicit\n");
155:         break;
156:       default:
157:         PetscViewerASCIIPrintf(viewer,"  transpose mode: not yet set\n");
158:     }
159:     if (svd->which == SVD_LARGEST) {
160:       PetscViewerASCIIPrintf(viewer,"  selected portion of the spectrum: largest\n");
161:     } else {
162:       PetscViewerASCIIPrintf(viewer,"  selected portion of the spectrum: smallest\n");
163:     }
164:     PetscViewerASCIIPrintf(viewer,"  number of singular values (nsv): %D\n",svd->nsv);
165:     PetscViewerASCIIPrintf(viewer,"  number of column vectors (ncv): %D\n",svd->ncv);
166:     PetscViewerASCIIPrintf(viewer,"  maximum dimension of projected problem (mpd): %D\n",svd->mpd);
167:     PetscViewerASCIIPrintf(viewer,"  maximum number of iterations: %D\n",svd->max_it);
168:     PetscViewerASCIIPrintf(viewer,"  tolerance: %G\n",svd->tol);
169:     if (svd->nini) {
170:       PetscViewerASCIIPrintf(viewer,"  dimension of user-provided initial space: %D\n",PetscAbs(svd->nini));
171:     }
172:     if (svd->ninil) {
173:       PetscViewerASCIIPrintf(viewer,"  dimension of user-provided initial left space: %D\n",PetscAbs(svd->ninil));
174:     }
175:   } else {
176:     if (svd->ops->view) {
177:       (*svd->ops->view)(svd,viewer);
178:     }
179:   }
180:   PetscObjectTypeCompareAny((PetscObject)svd,&isshell,SVDCROSS,SVDCYCLIC,"");
181:   if (!isshell) {
182:     if (!svd->ip) { SVDGetIP(svd,&svd->ip); }
183:     IPView(svd->ip,viewer);
184:     if (!svd->ds) { SVDGetDS(svd,&svd->ds); }
185:     PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
186:     DSView(svd->ds,viewer);
187:     PetscViewerPopFormat(viewer);
188:   }
189:   return(0);
190: }

194: /*@
195:    SVDPrintSolution - Prints the computed singular values.

197:    Collective on SVD

199:    Input Parameters:
200: +  svd - the singular value solver context
201: -  viewer - optional visualization context

203:    Options Database Key:
204: .  -svd_terse - print only minimal information

206:    Note:
207:    By default, this function prints a table with singular values and associated
208:    relative errors. With -svd_terse only the singular values are printed.

210:    Level: intermediate

212: .seealso: PetscViewerASCIIOpen()
213: @*/
214: PetscErrorCode SVDPrintSolution(SVD svd,PetscViewer viewer)
215: {
216:   PetscBool      terse,errok,isascii;
217:   PetscReal      error,sigma;
218:   PetscInt       i,j;

223:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)svd));
226:   if (!svd->sigma) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_WRONGSTATE,"SVDSolve must be called first");
227:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
228:   if (!isascii) return(0);

230:   PetscOptionsHasName(NULL,"-svd_terse",&terse);
231:   if (terse) {
232:     if (svd->nconv<svd->nsv) {
233:       PetscViewerASCIIPrintf(viewer," Problem: less than %D singular values converged\n\n",svd->nsv);
234:     } else {
235:       errok = PETSC_TRUE;
236:       for (i=0;i<svd->nsv;i++) {
237:         SVDComputeRelativeError(svd,i,&error);
238:         errok = (errok && error<svd->tol)? PETSC_TRUE: PETSC_FALSE;
239:       }
240:       if (errok) {
241:         PetscViewerASCIIPrintf(viewer," All requested singular values computed up to the required tolerance:");
242:         for (i=0;i<=(svd->nsv-1)/8;i++) {
243:           PetscViewerASCIIPrintf(viewer,"\n     ");
244:           for (j=0;j<PetscMin(8,svd->nsv-8*i);j++) {
245:             SVDGetSingularTriplet(svd,8*i+j,&sigma,NULL,NULL);
246:             PetscViewerASCIIPrintf(viewer,"%.5F",sigma);
247:             if (8*i+j+1<svd->nsv) { PetscViewerASCIIPrintf(viewer,", "); }
248:           }
249:         }
250:         PetscViewerASCIIPrintf(viewer,"\n\n");
251:       } else {
252:         PetscViewerASCIIPrintf(viewer," Problem: some of the first %D relative errors are higher than the tolerance\n\n",svd->nsv);
253:       }
254:     }
255:   } else {
256:     PetscViewerASCIIPrintf(viewer," Number of converged approximate singular triplets: %D\n\n",svd->nconv);
257:     if (svd->nconv>0) {
258:       PetscViewerASCIIPrintf(viewer,
259:            "          sigma            relative error\n"
260:            "   --------------------- ------------------\n");
261:       for (i=0;i<svd->nconv;i++) {
262:         SVDGetSingularTriplet(svd,i,&sigma,NULL,NULL);
263:         SVDComputeRelativeError(svd,i,&error);
264:         PetscViewerASCIIPrintf(viewer,"       % 6F          %12G\n",sigma,error);
265:       }
266:       PetscViewerASCIIPrintf(viewer,"\n");
267:     }
268:   }
269:   return(0);
270: }

274: /*@C
275:    SVDCreate - Creates the default SVD context.

277:    Collective on MPI_Comm

279:    Input Parameter:
280: .  comm - MPI communicator

282:    Output Parameter:
283: .  svd - location to put the SVD context

285:    Note:
286:    The default SVD type is SVDCROSS

288:    Level: beginner

290: .seealso: SVDSetUp(), SVDSolve(), SVDDestroy(), SVD
291: @*/
292: PetscErrorCode SVDCreate(MPI_Comm comm,SVD *outsvd)
293: {
295:   SVD            svd;

299:   *outsvd = 0;
300: #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
301:   SVDInitializePackage();
302: #endif

304:   SlepcHeaderCreate(svd,_p_SVD,struct _SVDOps,SVD_CLASSID,"SVD","Singular Value Decomposition","SVD",comm,SVDDestroy,SVDView);

306:   svd->OP             = NULL;
307:   svd->ip             = NULL;
308:   svd->ds             = NULL;
309:   svd->A              = NULL;
310:   svd->AT             = NULL;
311:   svd->transmode      = (SVDTransposeMode)PETSC_DECIDE;
312:   svd->sigma          = NULL;
313:   svd->perm           = NULL;
314:   svd->U              = NULL;
315:   svd->V              = NULL;
316:   svd->IS             = NULL;
317:   svd->ISL            = NULL;
318:   svd->tl             = NULL;
319:   svd->tr             = NULL;
320:   svd->rand           = NULL;
321:   svd->which          = SVD_LARGEST;
322:   svd->n              = 0;
323:   svd->nconv          = 0;
324:   svd->nsv            = 1;
325:   svd->ncv            = 0;
326:   svd->mpd            = 0;
327:   svd->nini           = 0;
328:   svd->ninil          = 0;
329:   svd->its            = 0;
330:   svd->max_it         = 0;
331:   svd->tol            = PETSC_DEFAULT;
332:   svd->errest         = NULL;
333:   svd->data           = NULL;
334:   svd->setupcalled    = 0;
335:   svd->reason         = SVD_CONVERGED_ITERATING;
336:   svd->numbermonitors = 0;
337:   svd->matvecs        = 0;
338:   svd->trackall       = PETSC_FALSE;

340:   PetscRandomCreate(comm,&svd->rand);
341:   PetscRandomSetSeed(svd->rand,0x12345678);
342:   PetscLogObjectParent(svd,svd->rand);
343:   *outsvd = svd;
344:   return(0);
345: }

349: /*@
350:    SVDReset - Resets the SVD context to the setupcalled=0 state and removes any
351:    allocated objects.

353:    Collective on SVD

355:    Input Parameter:
356: .  svd - singular value solver context obtained from SVDCreate()

358:    Level: advanced

360: .seealso: SVDDestroy()
361: @*/
362: PetscErrorCode SVDReset(SVD svd)
363: {

368:   if (svd->ops->reset) { (svd->ops->reset)(svd); }
369:   if (svd->ip) { IPReset(svd->ip); }
370:   if (svd->ds) { DSReset(svd->ds); }
371:   MatDestroy(&svd->OP);
372:   MatDestroy(&svd->A);
373:   MatDestroy(&svd->AT);
374:   VecDestroy(&svd->tl);
375:   VecDestroy(&svd->tr);
376:   if (svd->n) {
377:     PetscFree(svd->sigma);
378:     PetscFree(svd->perm);
379:     PetscFree(svd->errest);
380:     VecDestroyVecs(svd->n,&svd->U);
381:     VecDestroyVecs(svd->n,&svd->V);
382:   }
383:   svd->transmode   = (SVDTransposeMode)PETSC_DECIDE;
384:   svd->matvecs     = 0;
385:   svd->setupcalled = 0;
386:   return(0);
387: }

391: /*@C
392:    SVDDestroy - Destroys the SVD context.

394:    Collective on SVD

396:    Input Parameter:
397: .  svd - singular value solver context obtained from SVDCreate()

399:    Level: beginner

401: .seealso: SVDCreate(), SVDSetUp(), SVDSolve()
402: @*/
403: PetscErrorCode SVDDestroy(SVD *svd)
404: {

408:   if (!*svd) return(0);
410:   if (--((PetscObject)(*svd))->refct > 0) { *svd = 0; return(0); }
411:   SVDReset(*svd);
412:   if ((*svd)->ops->destroy) { (*(*svd)->ops->destroy)(*svd); }
413:   IPDestroy(&(*svd)->ip);
414:   DSDestroy(&(*svd)->ds);
415:   /* just in case the initial vectors have not been used */
416:   SlepcBasisDestroy_Private(&(*svd)->nini,&(*svd)->IS);
417:   SlepcBasisDestroy_Private(&(*svd)->ninil,&(*svd)->ISL);
418:   PetscRandomDestroy(&(*svd)->rand);
419:   SVDMonitorCancel(*svd);
420:   PetscHeaderDestroy(svd);
421:   return(0);
422: }

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

429:    Logically Collective on SVD

431:    Input Parameters:
432: +  svd      - the singular value solver context
433: -  type     - a known method

435:    Options Database Key:
436: .  -svd_type <method> - Sets the method; use -help for a list
437:     of available methods

439:    Notes:
440:    See "slepc/include/slepcsvd.h" for available methods. The default
441:    is SVDCROSS.

443:    Normally, it is best to use the SVDSetFromOptions() command and
444:    then set the SVD type from the options database rather than by using
445:    this routine.  Using the options database provides the user with
446:    maximum flexibility in evaluating the different available methods.
447:    The SVDSetType() routine is provided for those situations where it
448:    is necessary to set the iterative solver independently of the command
449:    line or options database.

451:    Level: intermediate

453: .seealso: SVDType
454: @*/
455: PetscErrorCode SVDSetType(SVD svd,SVDType type)
456: {
457:   PetscErrorCode ierr,(*r)(SVD);
458:   PetscBool      match;


464:   PetscObjectTypeCompare((PetscObject)svd,type,&match);
465:   if (match) return(0);

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

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

473:   svd->setupcalled = 0;
474:   PetscObjectChangeTypeName((PetscObject)svd,type);
475:   (*r)(svd);
476:   return(0);
477: }

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

484:    Not Collective

486:    Input Parameter:
487: .  svd - the singular value solver context

489:    Output Parameter:
490: .  name - name of SVD method

492:    Level: intermediate

494: .seealso: SVDSetType()
495: @*/
496: PetscErrorCode SVDGetType(SVD svd,SVDType *type)
497: {
501:   *type = ((PetscObject)svd)->type_name;
502:   return(0);
503: }

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

510:    Not Collective

512:    Input Parameters:
513: +  name - name of a new user-defined solver
514: -  function - routine to create the solver context

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

519:    Sample usage:
520: .vb
521:    SVDRegister("my_solver",MySolverCreate);
522: .ve

524:    Then, your solver can be chosen with the procedural interface via
525: $     SVDSetType(svd,"my_solver")
526:    or at runtime via the option
527: $     -svd_type my_solver

529:    Level: advanced

531: .seealso: SVDRegisterAll()
532: @*/
533: PetscErrorCode SVDRegister(const char *name,PetscErrorCode (*function)(SVD))
534: {

538:   PetscFunctionListAdd(&SVDList,name,function);
539:   return(0);
540: }

544: /*@
545:    SVDSetIP - Associates an inner product object to the
546:    singular value solver.

548:    Collective on SVD

550:    Input Parameters:
551: +  svd - singular value solver context obtained from SVDCreate()
552: -  ip  - the inner product object

554:    Note:
555:    Use SVDGetIP() to retrieve the inner product context (for example,
556:    to free it at the end of the computations).

558:    Level: advanced

560: .seealso: SVDGetIP()
561: @*/
562: PetscErrorCode SVDSetIP(SVD svd,IP ip)
563: {

570:   PetscObjectReference((PetscObject)ip);
571:   IPDestroy(&svd->ip);
572:   svd->ip = ip;
573:   PetscLogObjectParent(svd,svd->ip);
574:   return(0);
575: }

579: /*@C
580:    SVDGetIP - Obtain the inner product object associated
581:    to the singular value solver object.

583:    Not Collective

585:    Input Parameters:
586: .  svd - singular value solver context obtained from SVDCreate()

588:    Output Parameter:
589: .  ip - inner product context

591:    Level: advanced

593: .seealso: SVDSetIP()
594: @*/
595: PetscErrorCode SVDGetIP(SVD svd,IP *ip)
596: {

602:   if (!svd->ip) {
603:     IPCreate(PetscObjectComm((PetscObject)svd),&svd->ip);
604:     PetscLogObjectParent(svd,svd->ip);
605:   }
606:   *ip = svd->ip;
607:   return(0);
608: }

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

615:    Collective on SVD

617:    Input Parameters:
618: +  svd - singular value solver context obtained from SVDCreate()
619: -  ds  - the direct solver object

621:    Note:
622:    Use SVDGetDS() to retrieve the direct solver context (for example,
623:    to free it at the end of the computations).

625:    Level: advanced

627: .seealso: SVDGetDS()
628: @*/
629: PetscErrorCode SVDSetDS(SVD svd,DS ds)
630: {

637:   PetscObjectReference((PetscObject)ds);
638:   DSDestroy(&svd->ds);
639:   svd->ds = ds;
640:   PetscLogObjectParent(svd,svd->ds);
641:   return(0);
642: }

646: /*@C
647:    SVDGetDS - Obtain the direct solver object associated to the singular value
648:    solver object.

650:    Not Collective

652:    Input Parameters:
653: .  svd - singular value solver context obtained from SVDCreate()

655:    Output Parameter:
656: .  ds - direct solver context

658:    Level: advanced

660: .seealso: SVDSetDS()
661: @*/
662: PetscErrorCode SVDGetDS(SVD svd,DS *ds)
663: {

669:   if (!svd->ds) {
670:     DSCreate(PetscObjectComm((PetscObject)svd),&svd->ds);
671:     PetscLogObjectParent(svd,svd->ds);
672:   }
673:   *ds = svd->ds;
674:   return(0);
675: }