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-2007, Universidad Politecnica de Valencia, Spain

  8:       This file is part of SLEPc. See the README file for conditions of use
  9:       and additional information.
 10:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 11: */

 13:  #include src/svd/svdimpl.h

 15: PetscFList SVDList = 0;
 16: PetscCookie SVD_COOKIE = 0;
 17: PetscEvent SVD_SetUp = 0, SVD_Solve = 0, SVD_Dense = 0;

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

 26:   Input Parameter:
 27:   path - The dynamic library path, or PETSC_NULL

 29:   Level: developer

 31: .seealso: SlepcInitialize()
 32: @*/
 33: PetscErrorCode SVDInitializePackage(char *path)
 34: {
 35:   static PetscTruth initialized = PETSC_FALSE;
 36:   char              logList[256];
 37:   char              *className;
 38:   PetscTruth        opt;
 39:   PetscErrorCode    ierr;

 42:   if (initialized) return(0);
 43:   initialized = PETSC_TRUE;
 44:   /* Register Classes */
 45:   PetscLogClassRegister(&SVD_COOKIE,"Singular Value Solver");
 46:   /* Register Constructors */
 47:   SVDRegisterAll(path);
 48:   /* Register Events */
 49:   PetscLogEventRegister(&SVD_SetUp,"SVDSetUp",SVD_COOKIE);
 50:   PetscLogEventRegister(&SVD_Solve,"SVDSolve",SVD_COOKIE);
 51:   PetscLogEventRegister(&SVD_Dense,"SVDDense",SVD_COOKIE);
 52:   /* Process info exclusions */
 53:   PetscOptionsGetString(PETSC_NULL, "-log_info_exclude", logList, 256, &opt);
 54:   if (opt) {
 55:     PetscStrstr(logList, "svd", &className);
 56:     if (className) {
 57:       PetscInfoDeactivateClass(SVD_COOKIE);
 58:     }
 59:   }
 60:   /* Process summary exclusions */
 61:   PetscOptionsGetString(PETSC_NULL, "-log_summary_exclude", logList, 256, &opt);
 62:   if (opt) {
 63:     PetscStrstr(logList, "svd", &className);
 64:     if (className) {
 65:       PetscLogEventDeactivateClass(SVD_COOKIE);
 66:     }
 67:   }
 68:   return(0);
 69: }

 73: /*@C
 74:    SVDView - Prints the SVD data structure.

 76:    Collective on SVD

 78:    Input Parameters:
 79: +  svd - the singular value solver context
 80: -  viewer - optional visualization context

 82:    Options Database Key:
 83: .  -svd_view -  Calls SVDView() at end of SVDSolve()

 85:    Note:
 86:    The available visualization contexts include
 87: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
 88: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
 89:          output where only the first processor opens
 90:          the file.  All other processors send their 
 91:          data to the first processor to print. 

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

 96:    Level: beginner

 98: .seealso: STView(), PetscViewerASCIIOpen()
 99: @*/
100: PetscErrorCode SVDView(SVD svd,PetscViewer viewer)
101: {
103:   const char     *type;
104:   PetscTruth     isascii;

108:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(svd->comm);

112:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
113:   if (isascii) {
114:     PetscViewerASCIIPrintf(viewer,"SVD Object:\n");
115:     SVDGetType(svd,&type);
116:     if (type) {
117:       PetscViewerASCIIPrintf(viewer,"  method: %s\n",type);
118:     } else {
119:       PetscViewerASCIIPrintf(viewer,"  method: not yet set\n");
120:     }
121:     switch (svd->transmode) {
122:       case SVD_TRANSPOSE_EXPLICIT:
123:         PetscViewerASCIIPrintf(viewer,"  transpose mode: explicit\n");
124:         break;
125:       case SVD_TRANSPOSE_IMPLICIT:
126:         PetscViewerASCIIPrintf(viewer,"  transpose mode: implicit\n");
127:         break;
128:       default:
129:         PetscViewerASCIIPrintf(viewer,"  transpose mode: not yet set\n");
130:     }
131:     if (svd->which == SVD_LARGEST) {
132:       PetscViewerASCIIPrintf(viewer,"  selected portion of the spectrum: largest\n");
133:     } else {
134:       PetscViewerASCIIPrintf(viewer,"  selected portion of the spectrum: smallest\n");
135:     }
136:     PetscViewerASCIIPrintf(viewer,"  number of singular values (nsv): %d\n",svd->nsv);
137:     PetscViewerASCIIPrintf(viewer,"  number of column vectors (ncv): %d\n",svd->ncv);
138:     PetscViewerASCIIPrintf(viewer,"  maximum number of iterations: %d\n",svd->max_it);
139:     PetscViewerASCIIPrintf(viewer,"  tolerance: %g\n",svd->tol);
140:     if (svd->ops->view) {
141:       PetscViewerASCIIPushTab(viewer);
142:       (*svd->ops->view)(svd,viewer);
143:       PetscViewerASCIIPopTab(viewer);
144:     }
145:     PetscViewerASCIIPushTab(viewer);
146:     IPView(svd->ip,viewer);
147:     PetscViewerASCIIPopTab(viewer);
148:   } else {
149:     if (svd->ops->view) {
150:       (*svd->ops->view)(svd,viewer);
151:     }
152:   }
153:   return(0);
154: }

158: static PetscErrorCode SVDPublish_Petsc(PetscObject object)
159: {
161:   return(0);
162: }

166: /*@C
167:    SVDCreate - Creates the default SVD context.

169:    Collective on MPI_Comm

171:    Input Parameter:
172: .  comm - MPI communicator

174:    Output Parameter:
175: .  svd - location to put the SVD context

177:    Note:
178:    The default SVD type is SVDEIGENSOLVER

180:    Level: beginner

182: .seealso: SVDSetUp(), SVDSolve(), SVDDestroy(), SVD
183: @*/
184: PetscErrorCode SVDCreate(MPI_Comm comm,SVD *outsvd)
185: {
187:   SVD            svd;


192:   PetscHeaderCreate(svd,_p_SVD,struct _SVDOps,SVD_COOKIE,-1,"SVD",comm,SVDDestroy,SVDView);
193:   *outsvd = svd;

195:   svd->bops->publish   = SVDPublish_Petsc;
196:   PetscMemzero(svd->ops,sizeof(struct _SVDOps));

198:   svd->type_name   = PETSC_NULL;
199:   svd->OP          = PETSC_NULL;
200:   svd->A           = PETSC_NULL;
201:   svd->AT          = PETSC_NULL;
202:   svd->transmode   = (SVDTransposeMode)PETSC_DECIDE;
203:   svd->sigma       = PETSC_NULL;
204:   svd->U           = PETSC_NULL;
205:   svd->V           = PETSC_NULL;
206:   svd->vec_initial = PETSC_NULL;
207:   svd->which       = SVD_LARGEST;
208:   svd->n           = 0;
209:   svd->nconv       = 0;
210:   svd->nsv         = 1;
211:   svd->ncv         = PETSC_DECIDE;
212:   svd->its         = 0;
213:   svd->max_it      = PETSC_DECIDE;
214:   svd->tol         = 1e-7;
215:   svd->errest      = PETSC_NULL;
216:   svd->data        = PETSC_NULL;
217:   svd->setupcalled = 0;
218:   svd->reason      = SVD_CONVERGED_ITERATING;
219:   svd->numbermonitors = 0;
220:   svd->matvecs = 0;

222:   IPCreate(comm,&svd->ip);
223:   IPSetOptionsPrefix(svd->ip,svd->prefix);
224:   IPAppendOptionsPrefix(svd->ip,"svd_");
225:   PetscLogObjectParent(svd,svd->ip);

227:   PetscPublishAll(svd);
228:   return(0);
229: }
230: 
233: /*@
234:    SVDDestroy - Destroys the SVD context.

236:    Collective on SVD

238:    Input Parameter:
239: .  svd - singular value solver context obtained from SVDCreate()

241:    Level: beginner

243: .seealso: SVDCreate(), SVDSetUp(), SVDSolve()
244: @*/
245: PetscErrorCode SVDDestroy(SVD svd)
246: {
248:   int            i;
249: 
252:   if (--svd->refct > 0) return(0);

254:   /* if memory was published with AMS then destroy it */
255:   PetscObjectDepublish(svd);

257:   if (svd->ops->destroy) {
258:     (*svd->ops->destroy)(svd);
259:   }

261:   if (svd->OP) { MatDestroy(svd->OP); }
262:   if (svd->A) { MatDestroy(svd->A); }
263:   if (svd->AT) { MatDestroy(svd->AT); }
264:   if (svd->n) {
265:     PetscFree(svd->sigma);
266:     PetscFree(svd->errest);
267:     if (svd->U) {
268:       for (i=0;i<svd->n;i++) {
269:         VecDestroy(svd->U[i]);
270:       }
271:       PetscFree(svd->U);
272:     }
273:     for (i=0;i<svd->n;i++) {
274:       VecDestroy(svd->V[i]);
275:     }
276:     PetscFree(svd->V);
277:   }
278:   if (svd->vec_initial) { VecDestroy(svd->vec_initial); }
279:   SVDMonitorCancel(svd);
280: 
281:   IPDestroy(svd->ip);
282: 
283:   PetscHeaderDestroy(svd);
284:   return(0);
285: }

289: PetscErrorCode SVDDestroy_Default(SVD svd)
290: {

295:   PetscFree(svd->data);
296:   return(0);
297: }

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

304:    Collective on SVD

306:    Input Parameters:
307: +  svd      - the singular value solver context
308: -  type     - a known method

310:    Options Database Key:
311: .  -svd_type <method> - Sets the method; use -help for a list 
312:     of available methods 
313:     
314:    Notes:  
315:    See "slepc/include/slepcsvd.h" for available methods. The default
316:    is SVDEIGENSOLVER.

318:    Normally, it is best to use the SVDSetFromOptions() command and
319:    then set the SVD type from the options database rather than by using
320:    this routine.  Using the options database provides the user with
321:    maximum flexibility in evaluating the different available methods.
322:    The SVDSetType() routine is provided for those situations where it
323:    is necessary to set the iterative solver independently of the command
324:    line or options database. 

326:    Level: intermediate

328: .seealso: SVDType
329: @*/
330: PetscErrorCode SVDSetType(SVD svd,SVDType type)
331: {
332:   PetscErrorCode ierr,(*r)(SVD);
333:   PetscTruth match;


339:   PetscTypeCompare((PetscObject)svd,type,&match);
340:   if (match) return(0);

342:   if (svd->data) {
343:     /* destroy the old private SVD context */
344:     (*svd->ops->destroy)(svd);
345:     svd->data = 0;
346:   }

348:   PetscFListFind(SVDList,svd->comm,type,(void (**)(void)) &r);

350:   if (!r) SETERRQ1(1,"Unknown SVD type given: %s",type);

352:   svd->setupcalled = 0;
353:   PetscMemzero(svd->ops,sizeof(struct _SVDOps));
354:   (*r)(svd);

356:   PetscObjectChangeTypeName((PetscObject)svd,type);
357:   return(0);
358: }

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

365:    Not Collective

367:    Input Parameter:
368: .  svd - the singular value solver context 

370:    Output Parameter:
371: .  name - name of SVD method 

373:    Level: intermediate

375: .seealso: SVDSetType()
376: @*/
377: PetscErrorCode SVDGetType(SVD svd,SVDType *type)
378: {
381:   *type = svd->type_name;
382:   return(0);
383: }

385: /*MC
386:    SVDRegisterDynamic - Adds a method to the singular value solver package.

388:    Synopsis:
389:    SVDRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(SVD))

391:    Not Collective

393:    Input Parameters:
394: +  name_solver - name of a new user-defined solver
395: .  path - path (either absolute or relative) the library containing this solver
396: .  name_create - name of routine to create the solver context
397: -  routine_create - routine to create the solver context

399:    Notes:
400:    SVDRegisterDynamic() may be called multiple times to add several user-defined solvers.

402:    If dynamic libraries are used, then the fourth input argument (routine_create)
403:    is ignored.

405:    Sample usage:
406: .vb
407:    SVDRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
408:                "MySolverCreate",MySolverCreate);
409: .ve

411:    Then, your solver can be chosen with the procedural interface via
412: $     SVDSetType(svd,"my_solver")
413:    or at runtime via the option
414: $     -svd_type my_solver

416:    Level: advanced

418:    Environmental variables such as ${PETSC_ARCH}, ${SLEPC_DIR},
419:    and others of the form ${any_environmental_variable} occuring in pathname will be 
420:    replaced with appropriate values.

422: .seealso: SVDRegisterDestroy(), SVDRegisterAll()

424: M*/

428: /*@C
429:   SVDRegister - See SVDRegisterDynamic()

431:   Level: advanced
432: @*/
433: PetscErrorCode SVDRegister(const char *sname,const char *path,const char *name,int (*function)(SVD))
434: {
436:   char           fullname[256];

439:   PetscFListConcat(path,name,fullname);
440:   PetscFListAdd(&SVDList,sname,fullname,(void (*)(void))function);
441:   return(0);
442: }

446: /*@
447:    SVDRegisterDestroy - Frees the list of SVD methods that were
448:    registered by SVDRegisterDynamic().

450:    Not Collective

452:    Level: advanced

454: .seealso: SVDRegisterDynamic(), SVDRegisterAll()
455: @*/
456: PetscErrorCode SVDRegisterDestroy(void)
457: {

461:   PetscFListDestroy(&SVDList);
462:   SVDRegisterAll(PETSC_NULL);
463:   return(0);
464: }

468: /*@
469:    SVDSetIP - Associates an inner product object to the
470:    singular value solver. 

472:    Collective on SVD

474:    Input Parameters:
475: +  svd - singular value solver context obtained from SVDCreate()
476: -  ip  - the inner product object

478:    Note:
479:    Use SVDGetIP() to retrieve the inner product context (for example,
480:    to free it at the end of the computations).

482:    Level: advanced

484: .seealso: SVDGetIP()
485: @*/
486: PetscErrorCode SVDSetIP(SVD svd,IP ip)
487: {

494:   PetscObjectReference((PetscObject)ip);
495:   IPDestroy(svd->ip);
496:   svd->ip = ip;
497:   return(0);
498: }

502: /*@C
503:    SVDGetIP - Obtain the inner product object associated
504:    to the singular value solver object.

506:    Not Collective

508:    Input Parameters:
509: .  svd - singular value solver context obtained from SVDCreate()

511:    Output Parameter:
512: .  ip - inner product context

514:    Level: advanced

516: .seealso: SVDSetIP()
517: @*/
518: PetscErrorCode SVDGetIP(SVD svd,IP *ip)
519: {
523:   *ip = svd->ip;
524:   return(0);
525: }