Actual source code: qepbasic.c

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

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.
  9:       
 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 private/qepimpl.h

 26: PetscFList QEPList = 0;
 27: PetscCookie QEP_COOKIE = 0;
 28: PetscLogEvent QEP_SetUp = 0, QEP_Solve = 0, QEP_Dense = 0;
 29: static PetscTruth QEPPackageInitialized = PETSC_FALSE;

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

 37:   Level: developer

 39: .seealso: SlepcFinalize()
 40: @*/
 41: PetscErrorCode QEPFinalizePackage(void)
 42: {
 44:   QEPPackageInitialized = PETSC_FALSE;
 45:   QEPList               = 0;
 46:   return(0);
 47: }

 51: /*@C
 52:   QEPInitializePackage - This function initializes everything in the QEP package. It is called
 53:   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to QEPCreate()
 54:   when using static libraries.

 56:   Input Parameter:
 57:   path - The dynamic library path, or PETSC_NULL

 59:   Level: developer

 61: .seealso: SlepcInitialize()
 62: @*/
 63: PetscErrorCode QEPInitializePackage(char *path) {
 64:   char           logList[256];
 65:   char           *className;
 66:   PetscTruth     opt;

 70:   if (QEPPackageInitialized) return(0);
 71:   QEPPackageInitialized = PETSC_TRUE;
 72:   /* Register Classes */
 73:   PetscCookieRegister("Quadratic Eigenproblem Solver",&QEP_COOKIE);
 74:   /* Register Constructors */
 75:   QEPRegisterAll(path);
 76:   /* Register Events */
 77:   PetscLogEventRegister("QEPSetUp",QEP_COOKIE,&QEP_SetUp);
 78:   PetscLogEventRegister("QEPSolve",QEP_COOKIE,&QEP_Solve);
 79:   PetscLogEventRegister("QEPDense",QEP_COOKIE,&QEP_Dense);
 80:   /* Process info exclusions */
 81:   PetscOptionsGetString(PETSC_NULL,"-log_info_exclude",logList,256,&opt);
 82:   if (opt) {
 83:     PetscStrstr(logList,"qep",&className);
 84:     if (className) {
 85:       PetscInfoDeactivateClass(QEP_COOKIE);
 86:     }
 87:   }
 88:   /* Process summary exclusions */
 89:   PetscOptionsGetString(PETSC_NULL,"-log_summary_exclude",logList,256,&opt);
 90:   if (opt) {
 91:     PetscStrstr(logList,"qep",&className);
 92:     if (className) {
 93:       PetscLogEventDeactivateClass(QEP_COOKIE);
 94:     }
 95:   }
 96:   PetscRegisterFinalize(QEPFinalizePackage);
 97:   return(0);
 98: }

102: /*@C
103:    QEPView - Prints the QEP data structure.

105:    Collective on QEP

107:    Input Parameters:
108: +  qep - the quadratic eigenproblem solver context
109: -  viewer - optional visualization context

111:    Options Database Key:
112: .  -qep_view -  Calls QEPView() at end of QEPSolve()

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

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

125:    Level: beginner

127: .seealso: PetscViewerASCIIOpen()
128: @*/
129: PetscErrorCode QEPView(QEP qep,PetscViewer viewer)
130: {
132:   const char     *type;
133:   PetscTruth     isascii;

137:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(((PetscObject)qep)->comm);

141: #if defined(PETSC_USE_COMPLEX)
142: #define HERM "hermitian"
143: #else
144: #define HERM "symmetric"
145: #endif
146:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
147:   if (isascii) {
148:     PetscViewerASCIIPrintf(viewer,"QEP Object:\n");
149:     if (qep->problem_type) {
150:       switch (qep->problem_type) {
151:         case QEP_GENERAL:    type = "general quadratic eigenvalue problem"; break;
152:         case QEP_HERMITIAN:  type = HERM " quadratic eigenvalue problem"; break;
153:         case QEP_GYROSCOPIC: type = "gyroscopic quadratic eigenvalue problem"; break;
154:         default: SETERRQ(1,"Wrong value of qep->problem_type");
155:       }
156:     } else type = "not yet set";
157:     PetscViewerASCIIPrintf(viewer,"  problem type: %s\n",type);
158:     QEPGetType(qep,&type);
159:     if (type) {
160:       PetscViewerASCIIPrintf(viewer,"  method: %s\n",type);
161:     } else {
162:       PetscViewerASCIIPrintf(viewer,"  method: not yet set\n");
163:     }
164:     PetscViewerASCIIPrintf(viewer,"  selected portion of the spectrum: ");
165:     if (!qep->which) {
166:       PetscViewerASCIIPrintf(viewer,"not yet set\n");
167:     } else switch (qep->which) {
168:       case QEP_LARGEST_MAGNITUDE:
169:         PetscViewerASCIIPrintf(viewer,"largest eigenvalues in magnitude\n");
170:         break;
171:       case QEP_SMALLEST_MAGNITUDE:
172:         PetscViewerASCIIPrintf(viewer,"smallest eigenvalues in magnitude\n");
173:         break;
174:       case QEP_LARGEST_REAL:
175:         PetscViewerASCIIPrintf(viewer,"largest real parts\n");
176:         break;
177:       case QEP_SMALLEST_REAL:
178:         PetscViewerASCIIPrintf(viewer,"smallest real parts\n");
179:         break;
180:       case QEP_LARGEST_IMAGINARY:
181:         PetscViewerASCIIPrintf(viewer,"largest imaginary parts\n");
182:         break;
183:       case QEP_SMALLEST_IMAGINARY:
184:         PetscViewerASCIIPrintf(viewer,"smallest imaginary parts\n");
185:         break;
186:       default: SETERRQ(1,"Wrong value of qep->which");
187:     }
188:     if (qep->leftvecs) {
189:       PetscViewerASCIIPrintf(viewer,"  computing left eigenvectors also\n");
190:     }
191:     PetscViewerASCIIPrintf(viewer,"  number of eigenvalues (nev): %d\n",qep->nev);
192:     PetscViewerASCIIPrintf(viewer,"  number of column vectors (ncv): %d\n",qep->ncv);
193:     PetscViewerASCIIPrintf(viewer,"  maximum dimension of projected problem (mpd): %d\n",qep->mpd);
194:     PetscViewerASCIIPrintf(viewer,"  maximum number of iterations: %d\n", qep->max_it);
195:     PetscViewerASCIIPrintf(viewer,"  tolerance: %g\n",qep->tol);
196:     PetscViewerASCIIPrintf(viewer,"  scaling factor: %g\n",qep->sfactor);
197:     if (qep->nini!=0) {
198:       PetscViewerASCIIPrintf(viewer,"  dimension of user-provided initial space: %d\n",PetscAbs(qep->nini));
199:     }
200:     if (qep->ninil!=0) {
201:       PetscViewerASCIIPrintf(viewer,"  dimension of user-provided initial left space: %d\n",PetscAbs(qep->ninil));
202:     }
203:     if (qep->ops->view) {
204:       PetscViewerASCIIPushTab(viewer);
205:       (*qep->ops->view)(qep,viewer);
206:       PetscViewerASCIIPopTab(viewer);
207:     }
208:     PetscViewerASCIIPushTab(viewer);
209:     IPView(qep->ip,viewer);
210:     PetscViewerASCIIPopTab(viewer);
211:   } else {
212:     if (qep->ops->view) {
213:       (*qep->ops->view)(qep,viewer);
214:     }
215:   }
216:   return(0);
217: }

221: /*@C
222:    QEPCreate - Creates the default QEP context.

224:    Collective on MPI_Comm

226:    Input Parameter:
227: .  comm - MPI communicator

229:    Output Parameter:
230: .  qep - location to put the QEP context

232:    Note:
233:    The default QEP type is QEPLINEAR

235:    Level: beginner

237: .seealso: QEPSetUp(), QEPSolve(), QEPDestroy(), QEP
238: @*/
239: PetscErrorCode QEPCreate(MPI_Comm comm,QEP *outqep)
240: {
242:   QEP            qep;

246:   *outqep = 0;

248:   PetscHeaderCreate(qep,_p_QEP,struct _QEPOps,QEP_COOKIE,-1,"QEP",comm,QEPDestroy,QEPView);
249:   *outqep = qep;

251:   PetscMemzero(qep->ops,sizeof(struct _QEPOps));

253:   qep->max_it          = 0;
254:   qep->nev             = 1;
255:   qep->ncv             = 0;
256:   qep->mpd             = 0;
257:   qep->nini            = 0;
258:   qep->ninil           = 0;
259:   qep->tol             = 1e-7;
260:   qep->sfactor         = 0.0;
261:   qep->conv_func       = QEPDefaultConverged;
262:   qep->conv_ctx        = PETSC_NULL;
263:   qep->which           = (QEPWhich)0;
264:   qep->which_func      = PETSC_NULL;
265:   qep->which_ctx       = PETSC_NULL;
266:   qep->leftvecs        = PETSC_FALSE;
267:   qep->problem_type    = (QEPProblemType)0;
268:   qep->V               = PETSC_NULL;
269:   qep->IS              = PETSC_NULL;
270:   qep->ISL             = PETSC_NULL;
271:   qep->T               = PETSC_NULL;
272:   qep->eigr            = PETSC_NULL;
273:   qep->eigi            = PETSC_NULL;
274:   qep->errest          = PETSC_NULL;
275:   qep->data            = PETSC_NULL;
276:   qep->nconv           = 0;
277:   qep->its             = 0;
278:   qep->perm            = PETSC_NULL;
279:   qep->matvecs         = 0;
280:   qep->linits          = 0;
281:   qep->nwork           = 0;
282:   qep->work            = PETSC_NULL;
283:   qep->setupcalled     = 0;
284:   qep->reason          = QEP_CONVERGED_ITERATING;
285:   qep->numbermonitors  = 0;
286:   qep->trackall        = PETSC_FALSE;

288:   IPCreate(comm,&qep->ip);
289:   IPSetOptionsPrefix(qep->ip,((PetscObject)qep)->prefix);
290:   IPAppendOptionsPrefix(qep->ip,"qep_");
291:   PetscLogObjectParent(qep,qep->ip);
292:   PetscPublishAll(qep);
293:   return(0);
294: }
295: 
298: /*@C
299:    QEPSetType - Selects the particular solver to be used in the QEP object. 

301:    Collective on QEP

303:    Input Parameters:
304: +  qep      - the quadratic eigensolver context
305: -  type     - a known method

307:    Options Database Key:
308: .  -qep_type <method> - Sets the method; use -help for a list 
309:     of available methods 
310:     
311:    Notes:  
312:    See "slepc/include/slepcqep.h" for available methods. The default
313:    is QEPLINEAR.

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

323:    Level: intermediate

325: .seealso: QEPType
326: @*/
327: PetscErrorCode QEPSetType(QEP qep,const QEPType type)
328: {
329:   PetscErrorCode ierr,(*r)(QEP);
330:   PetscTruth match;


336:   PetscTypeCompare((PetscObject)qep,type,&match);
337:   if (match) return(0);

339:   if (qep->data) {
340:     /* destroy the old private QEP context */
341:     (*qep->ops->destroy)(qep);
342:     qep->data = 0;
343:   }

345:   PetscFListFind(QEPList,((PetscObject)qep)->comm,type,(void (**)(void))&r);

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

349:   qep->setupcalled = 0;
350:   PetscMemzero(qep->ops,sizeof(struct _QEPOps));
351:   (*r)(qep);

353:   PetscObjectChangeTypeName((PetscObject)qep,type);
354:   return(0);
355: }

359: /*@C
360:    QEPGetType - Gets the QEP type as a string from the QEP object.

362:    Not Collective

364:    Input Parameter:
365: .  qep - the eigensolver context 

367:    Output Parameter:
368: .  name - name of QEP method 

370:    Level: intermediate

372: .seealso: QEPSetType()
373: @*/
374: PetscErrorCode QEPGetType(QEP qep,const QEPType *type)
375: {
379:   *type = ((PetscObject)qep)->type_name;
380:   return(0);
381: }

383: /*MC
384:    QEPRegisterDynamic - Adds a method to the quadratic eigenproblem solver package.

386:    Synopsis:
387:    QEPRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(QEP))

389:    Not Collective

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

397:    Notes:
398:    QEPRegisterDynamic() may be called multiple times to add several user-defined solvers.

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

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

409:    Then, your solver can be chosen with the procedural interface via
410: $     QEPSetType(qep,"my_solver")
411:    or at runtime via the option
412: $     -qep_type my_solver

414:    Level: advanced

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

420: .seealso: QEPRegisterDestroy(), QEPRegisterAll()

422: M*/

426: /*@C
427:   QEPRegister - See QEPRegisterDynamic()

429:   Level: advanced
430: @*/
431: PetscErrorCode QEPRegister(const char *sname,const char *path,const char *name,PetscErrorCode (*function)(QEP))
432: {
434:   char           fullname[256];

437:   PetscFListConcat(path,name,fullname);
438:   PetscFListAdd(&QEPList,sname,fullname,(void (*)(void))function);
439:   return(0);
440: }

444: /*@
445:    QEPRegisterDestroy - Frees the list of QEP methods that were
446:    registered by QEPRegisterDynamic().

448:    Not Collective

450:    Level: advanced

452: .seealso: QEPRegisterDynamic(), QEPRegisterAll()
453: @*/
454: PetscErrorCode QEPRegisterDestroy(void)
455: {

459:   PetscFListDestroy(&QEPList);
460:   QEPRegisterAll(PETSC_NULL);
461:   return(0);
462: }

466: /*@
467:    QEPDestroy - Destroys the QEP context.

469:    Collective on QEP

471:    Input Parameter:
472: .  qep - eigensolver context obtained from QEPCreate()

474:    Level: beginner

476: .seealso: QEPCreate(), QEPSetUp(), QEPSolve()
477: @*/
478: PetscErrorCode QEPDestroy(QEP qep)
479: {
481:   PetscScalar    *pV;
482:   PetscInt       i;

486:   if (--((PetscObject)qep)->refct > 0) return(0);

488:   /* if memory was published with AMS then destroy it */
489:   PetscObjectDepublish(qep);

491:   if (qep->ops->destroy) {
492:     (*qep->ops->destroy)(qep);
493:   }
494: 
495:   PetscFree(qep->T);

497:   if (qep->eigr) {
498:     PetscFree(qep->eigr);
499:     PetscFree(qep->eigi);
500:     PetscFree(qep->perm);
501:     PetscFree(qep->errest);
502:     VecGetArray(qep->V[0],&pV);
503:     for (i=0;i<qep->ncv;i++) {
504:       VecDestroy(qep->V[i]);
505:     }
506:     PetscFree(pV);
507:     PetscFree(qep->V);
508:   }

510:   QEPMonitorCancel(qep);

512:   IPDestroy(qep->ip);
513:   if (qep->rand) {
514:     PetscRandomDestroy(qep->rand);
515:   }

517:   if (qep->M) { MatDestroy(qep->M); }
518:   if (qep->C) { MatDestroy(qep->C); }
519:   if (qep->K) { MatDestroy(qep->K); }

521:   PetscHeaderDestroy(qep);
522:   return(0);
523: }

527: /*@
528:    QEPSetIP - Associates an inner product object to the quadratic eigensolver. 

530:    Collective on QEP and IP

532:    Input Parameters:
533: +  qep - eigensolver context obtained from QEPCreate()
534: -  ip  - the inner product object

536:    Note:
537:    Use QEPGetIP() to retrieve the inner product context (for example,
538:    to free it at the end of the computations).

540:    Level: advanced

542: .seealso: QEPGetIP()
543: @*/
544: PetscErrorCode QEPSetIP(QEP qep,IP ip)
545: {

552:   PetscObjectReference((PetscObject)ip);
553:   IPDestroy(qep->ip);
554:   qep->ip = ip;
555:   return(0);
556: }

560: /*@C
561:    QEPGetIP - Obtain the inner product object associated
562:    to the quadratic eigensolver object.

564:    Not Collective

566:    Input Parameters:
567: .  qep - eigensolver context obtained from QEPCreate()

569:    Output Parameter:
570: .  ip - inner product context

572:    Level: advanced

574: .seealso: QEPSetIP()
575: @*/
576: PetscErrorCode QEPGetIP(QEP qep,IP *ip)
577: {
581:   *ip = qep->ip;
582:   return(0);
583: }