Subversion Repositories slepc-dev

Rev

Go to most recent revision | Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
1887 jroman 1
/*
2
     The basic QEP routines, Create, View, etc. are here.
3
 
4
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5
   SLEPc - Scalable Library for Eigenvalue Problem Computations
2116 eromero 6
   Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain
1887 jroman 7
 
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.
13
 
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.
18
 
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
*/
23
 
2284 jroman 24
#include <private/qepimpl.h>      /*I "slepcqep.h" I*/
1887 jroman 25
 
2317 jroman 26
PetscFList       QEPList = 0;
27
PetscClassId     QEP_CLASSID = 0;
28
PetscLogEvent    QEP_SetUp = 0, QEP_Solve = 0, QEP_Dense = 0;
2216 jroman 29
static PetscBool QEPPackageInitialized = PETSC_FALSE;
1887 jroman 30
 
31
#undef __FUNCT__  
32
#define __FUNCT__ "QEPFinalizePackage"
33
/*@C
2317 jroman 34
   QEPFinalizePackage - This function destroys everything in the Slepc interface
35
   to the QEP package. It is called from SlepcFinalize().
1887 jroman 36
 
2317 jroman 37
   Level: developer
1887 jroman 38
 
39
.seealso: SlepcFinalize()
40
@*/
41
PetscErrorCode QEPFinalizePackage(void)
42
{
43
  PetscFunctionBegin;
44
  QEPPackageInitialized = PETSC_FALSE;
45
  QEPList               = 0;
46
  PetscFunctionReturn(0);
47
}
48
 
49
#undef __FUNCT__  
50
#define __FUNCT__ "QEPInitializePackage"
51
/*@C
2317 jroman 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.
1887 jroman 55
 
2317 jroman 56
   Input Parameter:
57
.  path - The dynamic library path, or PETSC_NULL
1887 jroman 58
 
2317 jroman 59
   Level: developer
1887 jroman 60
 
61
.seealso: SlepcInitialize()
62
@*/
2317 jroman 63
PetscErrorCode QEPInitializePackage(const char *path)
64
{
1887 jroman 65
  char           logList[256];
66
  char           *className;
2216 jroman 67
  PetscBool      opt;
1887 jroman 68
  PetscErrorCode ierr;
69
 
70
  PetscFunctionBegin;
71
  if (QEPPackageInitialized) PetscFunctionReturn(0);
72
  QEPPackageInitialized = PETSC_TRUE;
73
  /* Register Classes */
2213 jroman 74
  ierr = PetscClassIdRegister("Quadratic Eigenproblem Solver",&QEP_CLASSID);CHKERRQ(ierr);
1887 jroman 75
  /* Register Constructors */
76
  ierr = QEPRegisterAll(path);CHKERRQ(ierr);
77
  /* Register Events */
2213 jroman 78
  ierr = PetscLogEventRegister("QEPSetUp",QEP_CLASSID,&QEP_SetUp);CHKERRQ(ierr);
79
  ierr = PetscLogEventRegister("QEPSolve",QEP_CLASSID,&QEP_Solve);CHKERRQ(ierr);
80
  ierr = PetscLogEventRegister("QEPDense",QEP_CLASSID,&QEP_Dense);CHKERRQ(ierr);
1887 jroman 81
  /* Process info exclusions */
2213 jroman 82
  ierr = PetscOptionsGetString(PETSC_NULL,"-info_exclude",logList,256,&opt);CHKERRQ(ierr);
1887 jroman 83
  if (opt) {
84
    ierr = PetscStrstr(logList,"qep",&className);CHKERRQ(ierr);
85
    if (className) {
2213 jroman 86
      ierr = PetscInfoDeactivateClass(QEP_CLASSID);CHKERRQ(ierr);
1887 jroman 87
    }
88
  }
89
  /* Process summary exclusions */
90
  ierr = PetscOptionsGetString(PETSC_NULL,"-log_summary_exclude",logList,256,&opt);CHKERRQ(ierr);
91
  if (opt) {
92
    ierr = PetscStrstr(logList,"qep",&className);CHKERRQ(ierr);
93
    if (className) {
2213 jroman 94
      ierr = PetscLogEventDeactivateClass(QEP_CLASSID);CHKERRQ(ierr);
1887 jroman 95
    }
96
  }
97
  ierr = PetscRegisterFinalize(QEPFinalizePackage);CHKERRQ(ierr);
98
  PetscFunctionReturn(0);
99
}
100
 
101
#undef __FUNCT__  
102
#define __FUNCT__ "QEPView"
103
/*@C
104
   QEPView - Prints the QEP data structure.
105
 
106
   Collective on QEP
107
 
108
   Input Parameters:
109
+  qep - the quadratic eigenproblem solver context
110
-  viewer - optional visualization context
111
 
112
   Options Database Key:
113
.  -qep_view -  Calls QEPView() at end of QEPSolve()
114
 
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.
122
 
123
   The user can open an alternative visualization context with
124
   PetscViewerASCIIOpen() - output to a specified file.
125
 
126
   Level: beginner
127
 
128
.seealso: PetscViewerASCIIOpen()
129
@*/
130
PetscErrorCode QEPView(QEP qep,PetscViewer viewer)
131
{
132
  PetscErrorCode ierr;
133
  const char     *type;
2216 jroman 134
  PetscBool      isascii;
1887 jroman 135
 
136
  PetscFunctionBegin;
2213 jroman 137
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1887 jroman 138
  if (!viewer) viewer = PETSC_VIEWER_STDOUT_(((PetscObject)qep)->comm);
2213 jroman 139
  PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1887 jroman 140
  PetscCheckSameComm(qep,1,viewer,2);
141
 
1907 jroman 142
#if defined(PETSC_USE_COMPLEX)
143
#define HERM "hermitian"
144
#else
145
#define HERM "symmetric"
146
#endif
2215 jroman 147
  ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
1887 jroman 148
  if (isascii) {
2224 jroman 149
    ierr = PetscObjectPrintClassNamePrefixType((PetscObject)qep,viewer,"QEP Object");CHKERRQ(ierr);
2231 jroman 150
    if (qep->problem_type) {
151
      switch (qep->problem_type) {
152
        case QEP_GENERAL:    type = "general quadratic eigenvalue problem"; break;
153
        case QEP_HERMITIAN:  type = HERM " quadratic eigenvalue problem"; break;
154
        case QEP_GYROSCOPIC: type = "gyroscopic quadratic eigenvalue problem"; break;
155
        default: SETERRQ(((PetscObject)qep)->comm,1,"Wrong value of qep->problem_type");
156
      }
157
    } else type = "not yet set";
1907 jroman 158
    ierr = PetscViewerASCIIPrintf(viewer,"  problem type: %s\n",type);CHKERRQ(ierr);
1887 jroman 159
    ierr = PetscViewerASCIIPrintf(viewer,"  selected portion of the spectrum: ");CHKERRQ(ierr);
1942 jroman 160
    if (!qep->which) {
161
      ierr = PetscViewerASCIIPrintf(viewer,"not yet set\n");CHKERRQ(ierr);
162
    } else switch (qep->which) {
1887 jroman 163
      case QEP_LARGEST_MAGNITUDE:
164
        ierr = PetscViewerASCIIPrintf(viewer,"largest eigenvalues in magnitude\n");CHKERRQ(ierr);
165
        break;
166
      case QEP_SMALLEST_MAGNITUDE:
167
        ierr = PetscViewerASCIIPrintf(viewer,"smallest eigenvalues in magnitude\n");CHKERRQ(ierr);
168
        break;
169
      case QEP_LARGEST_REAL:
170
        ierr = PetscViewerASCIIPrintf(viewer,"largest real parts\n");CHKERRQ(ierr);
171
        break;
172
      case QEP_SMALLEST_REAL:
173
        ierr = PetscViewerASCIIPrintf(viewer,"smallest real parts\n");CHKERRQ(ierr);
174
        break;
175
      case QEP_LARGEST_IMAGINARY:
176
        ierr = PetscViewerASCIIPrintf(viewer,"largest imaginary parts\n");CHKERRQ(ierr);
177
        break;
178
      case QEP_SMALLEST_IMAGINARY:
179
        ierr = PetscViewerASCIIPrintf(viewer,"smallest imaginary parts\n");CHKERRQ(ierr);
180
        break;
2214 jroman 181
      default: SETERRQ(((PetscObject)qep)->comm,1,"Wrong value of qep->which");
1887 jroman 182
    }    
1944 jroman 183
    if (qep->leftvecs) {
184
      ierr = PetscViewerASCIIPrintf(viewer,"  computing left eigenvectors also\n");CHKERRQ(ierr);
185
    }
1887 jroman 186
    ierr = PetscViewerASCIIPrintf(viewer,"  number of eigenvalues (nev): %d\n",qep->nev);CHKERRQ(ierr);
187
    ierr = PetscViewerASCIIPrintf(viewer,"  number of column vectors (ncv): %d\n",qep->ncv);CHKERRQ(ierr);
188
    ierr = PetscViewerASCIIPrintf(viewer,"  maximum dimension of projected problem (mpd): %d\n",qep->mpd);CHKERRQ(ierr);
189
    ierr = PetscViewerASCIIPrintf(viewer,"  maximum number of iterations: %d\n", qep->max_it);
190
    ierr = PetscViewerASCIIPrintf(viewer,"  tolerance: %g\n",qep->tol);CHKERRQ(ierr);
1918 jroman 191
    ierr = PetscViewerASCIIPrintf(viewer,"  scaling factor: %g\n",qep->sfactor);CHKERRQ(ierr);
1952 jroman 192
    if (qep->nini!=0) {
193
      ierr = PetscViewerASCIIPrintf(viewer,"  dimension of user-provided initial space: %d\n",PetscAbs(qep->nini));CHKERRQ(ierr);
194
    }
195
    if (qep->ninil!=0) {
196
      ierr = PetscViewerASCIIPrintf(viewer,"  dimension of user-provided initial left space: %d\n",PetscAbs(qep->ninil));CHKERRQ(ierr);
197
    }
1887 jroman 198
    if (qep->ops->view) {
199
      ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
200
      ierr = (*qep->ops->view)(qep,viewer);CHKERRQ(ierr);
201
      ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
202
    }
203
    ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
204
    ierr = IPView(qep->ip,viewer); CHKERRQ(ierr);
205
    ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
206
  } else {
207
    if (qep->ops->view) {
208
      ierr = (*qep->ops->view)(qep,viewer);CHKERRQ(ierr);
209
    }
210
  }
211
  PetscFunctionReturn(0);
212
}
213
 
214
#undef __FUNCT__  
215
#define __FUNCT__ "QEPCreate"
216
/*@C
217
   QEPCreate - Creates the default QEP context.
218
 
219
   Collective on MPI_Comm
220
 
221
   Input Parameter:
222
.  comm - MPI communicator
223
 
224
   Output Parameter:
225
.  qep - location to put the QEP context
226
 
227
   Note:
228
   The default QEP type is QEPLINEAR
229
 
230
   Level: beginner
231
 
232
.seealso: QEPSetUp(), QEPSolve(), QEPDestroy(), QEP
233
@*/
234
PetscErrorCode QEPCreate(MPI_Comm comm,QEP *outqep)
235
{
236
  PetscErrorCode ierr;
237
  QEP            qep;
238
 
239
  PetscFunctionBegin;
240
  PetscValidPointer(outqep,2);
241
  *outqep = 0;
242
 
2213 jroman 243
  ierr = PetscHeaderCreate(qep,_p_QEP,struct _QEPOps,QEP_CLASSID,-1,"QEP",comm,QEPDestroy,QEPView);CHKERRQ(ierr);
1887 jroman 244
  *outqep = qep;
245
 
246
  ierr = PetscMemzero(qep->ops,sizeof(struct _QEPOps));CHKERRQ(ierr);
247
 
2309 jroman 248
  qep->M               = 0;
249
  qep->C               = 0;
250
  qep->K               = 0;
1887 jroman 251
  qep->max_it          = 0;
252
  qep->nev             = 1;
253
  qep->ncv             = 0;
254
  qep->mpd             = 0;
1952 jroman 255
  qep->nini            = 0;
256
  qep->ninil           = 0;
1887 jroman 257
  qep->tol             = 1e-7;
1918 jroman 258
  qep->sfactor         = 0.0;
1887 jroman 259
  qep->conv_func       = QEPDefaultConverged;
260
  qep->conv_ctx        = PETSC_NULL;
1942 jroman 261
  qep->which           = (QEPWhich)0;
1887 jroman 262
  qep->which_func      = PETSC_NULL;
263
  qep->which_ctx       = PETSC_NULL;
1944 jroman 264
  qep->leftvecs        = PETSC_FALSE;
1906 jroman 265
  qep->problem_type    = (QEPProblemType)0;
1887 jroman 266
  qep->V               = PETSC_NULL;
2309 jroman 267
  qep->W               = PETSC_NULL;
1952 jroman 268
  qep->IS              = PETSC_NULL;
269
  qep->ISL             = PETSC_NULL;
1887 jroman 270
  qep->T               = PETSC_NULL;
271
  qep->eigr            = PETSC_NULL;
272
  qep->eigi            = PETSC_NULL;
273
  qep->errest          = PETSC_NULL;
274
  qep->data            = PETSC_NULL;
275
  qep->nconv           = 0;
276
  qep->its             = 0;
277
  qep->perm            = PETSC_NULL;
1896 jroman 278
  qep->matvecs         = 0;
279
  qep->linits          = 0;
1887 jroman 280
  qep->nwork           = 0;
281
  qep->work            = PETSC_NULL;
282
  qep->setupcalled     = 0;
283
  qep->reason          = QEP_CONVERGED_ITERATING;
284
  qep->numbermonitors  = 0;
2054 eromero 285
  qep->trackall        = PETSC_FALSE;
2309 jroman 286
  qep->rand            = 0;
1887 jroman 287
 
288
  ierr = IPCreate(comm,&qep->ip); CHKERRQ(ierr);
289
  ierr = IPSetOptionsPrefix(qep->ip,((PetscObject)qep)->prefix);
290
  ierr = IPAppendOptionsPrefix(qep->ip,"qep_");
291
  PetscLogObjectParent(qep,qep->ip);
292
  PetscFunctionReturn(0);
293
}
294
 
295
#undef __FUNCT__  
296
#define __FUNCT__ "QEPSetType"
297
/*@C
298
   QEPSetType - Selects the particular solver to be used in the QEP object.
299
 
300
   Collective on QEP
301
 
302
   Input Parameters:
303
+  qep      - the quadratic eigensolver context
304
-  type     - a known method
305
 
306
   Options Database Key:
307
.  -qep_type <method> - Sets the method; use -help for a list
308
    of available methods
309
 
310
   Notes:  
311
   See "slepc/include/slepcqep.h" for available methods. The default
312
   is QEPLINEAR.
313
 
314
   Normally, it is best to use the QEPSetFromOptions() command and
315
   then set the QEP type from the options database rather than by using
316
   this routine.  Using the options database provides the user with
317
   maximum flexibility in evaluating the different available methods.
318
   The QEPSetType() routine is provided for those situations where it
319
   is necessary to set the iterative solver independently of the command
320
   line or options database.
321
 
322
   Level: intermediate
323
 
324
.seealso: QEPType
325
@*/
326
PetscErrorCode QEPSetType(QEP qep,const QEPType type)
327
{
328
  PetscErrorCode ierr,(*r)(QEP);
2216 jroman 329
  PetscBool      match;
1887 jroman 330
 
331
  PetscFunctionBegin;
2213 jroman 332
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1887 jroman 333
  PetscValidCharPointer(type,2);
334
  ierr = PetscTypeCompare((PetscObject)qep,type,&match);CHKERRQ(ierr);
335
  if (match) PetscFunctionReturn(0);
336
 
337
  if (qep->data) {
338
    /* destroy the old private QEP context */
339
    ierr = (*qep->ops->destroy)(qep); CHKERRQ(ierr);
340
    qep->data = 0;
341
  }
342
 
2300 jroman 343
  ierr = PetscFListFind(QEPList,((PetscObject)qep)->comm,type,PETSC_TRUE,(void (**)(void))&r);CHKERRQ(ierr);
1887 jroman 344
 
2214 jroman 345
  if (!r) SETERRQ1(((PetscObject)qep)->comm,1,"Unknown QEP type given: %s",type);
1887 jroman 346
 
347
  qep->setupcalled = 0;
348
  ierr = PetscMemzero(qep->ops,sizeof(struct _QEPOps));CHKERRQ(ierr);
349
  ierr = (*r)(qep);CHKERRQ(ierr);
350
 
351
  ierr = PetscObjectChangeTypeName((PetscObject)qep,type);CHKERRQ(ierr);
352
  PetscFunctionReturn(0);
353
}
354
 
355
#undef __FUNCT__  
356
#define __FUNCT__ "QEPGetType"
357
/*@C
358
   QEPGetType - Gets the QEP type as a string from the QEP object.
359
 
360
   Not Collective
361
 
362
   Input Parameter:
363
.  qep - the eigensolver context
364
 
365
   Output Parameter:
366
.  name - name of QEP method
367
 
368
   Level: intermediate
369
 
370
.seealso: QEPSetType()
371
@*/
372
PetscErrorCode QEPGetType(QEP qep,const QEPType *type)
373
{
374
  PetscFunctionBegin;
2213 jroman 375
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1887 jroman 376
  PetscValidPointer(type,2);
377
  *type = ((PetscObject)qep)->type_name;
378
  PetscFunctionReturn(0);
379
}
380
 
381
/*MC
382
   QEPRegisterDynamic - Adds a method to the quadratic eigenproblem solver package.
383
 
384
   Synopsis:
385
   QEPRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(QEP))
386
 
387
   Not Collective
388
 
389
   Input Parameters:
390
+  name_solver - name of a new user-defined solver
391
.  path - path (either absolute or relative) the library containing this solver
392
.  name_create - name of routine to create the solver context
393
-  routine_create - routine to create the solver context
394
 
395
   Notes:
396
   QEPRegisterDynamic() may be called multiple times to add several user-defined solvers.
397
 
398
   If dynamic libraries are used, then the fourth input argument (routine_create)
399
   is ignored.
400
 
401
   Sample usage:
402
.vb
403
   QEPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
404
               "MySolverCreate",MySolverCreate);
405
.ve
406
 
407
   Then, your solver can be chosen with the procedural interface via
408
$     QEPSetType(qep,"my_solver")
409
   or at runtime via the option
410
$     -qep_type my_solver
411
 
412
   Level: advanced
413
 
414
   Environmental variables such as ${PETSC_ARCH}, ${SLEPC_DIR},
415
   and others of the form ${any_environmental_variable} occuring in pathname will be
416
   replaced with appropriate values.
417
 
418
.seealso: QEPRegisterDestroy(), QEPRegisterAll()
419
 
420
M*/
421
 
422
#undef __FUNCT__  
423
#define __FUNCT__ "QEPRegister"
424
/*@C
425
  QEPRegister - See QEPRegisterDynamic()
426
 
427
  Level: advanced
428
@*/
429
PetscErrorCode QEPRegister(const char *sname,const char *path,const char *name,PetscErrorCode (*function)(QEP))
430
{
431
  PetscErrorCode ierr;
432
  char           fullname[256];
433
 
434
  PetscFunctionBegin;
435
  ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
436
  ierr = PetscFListAdd(&QEPList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
437
  PetscFunctionReturn(0);
438
}
439
 
440
#undef __FUNCT__  
441
#define __FUNCT__ "QEPRegisterDestroy"
442
/*@
443
   QEPRegisterDestroy - Frees the list of QEP methods that were
444
   registered by QEPRegisterDynamic().
445
 
446
   Not Collective
447
 
448
   Level: advanced
449
 
450
.seealso: QEPRegisterDynamic(), QEPRegisterAll()
451
@*/
452
PetscErrorCode QEPRegisterDestroy(void)
453
{
454
  PetscErrorCode ierr;
455
 
456
  PetscFunctionBegin;
457
  ierr = PetscFListDestroy(&QEPList);CHKERRQ(ierr);
458
  ierr = QEPRegisterAll(PETSC_NULL);CHKERRQ(ierr);
459
  PetscFunctionReturn(0);
460
}
461
 
462
#undef __FUNCT__  
463
#define __FUNCT__ "QEPDestroy"
2312 jroman 464
/*@C
1887 jroman 465
   QEPDestroy - Destroys the QEP context.
466
 
467
   Collective on QEP
468
 
469
   Input Parameter:
470
.  qep - eigensolver context obtained from QEPCreate()
471
 
472
   Level: beginner
473
 
474
.seealso: QEPCreate(), QEPSetUp(), QEPSolve()
475
@*/
2312 jroman 476
PetscErrorCode QEPDestroy(QEP *qep)
1887 jroman 477
{
478
  PetscErrorCode ierr;
1903 jroman 479
  PetscScalar    *pV;
480
  PetscInt       i;
1887 jroman 481
 
482
  PetscFunctionBegin;
2312 jroman 483
  if (!*qep) PetscFunctionReturn(0);
484
  PetscValidHeaderSpecific(*qep,QEP_CLASSID,1);
485
  if (--((PetscObject)(*qep))->refct > 0) { *qep = 0; PetscFunctionReturn(0); }
486
  ierr = PetscObjectDepublish(*qep);CHKERRQ(ierr);
487
  if ((*qep)->ops->destroy) {
488
    ierr = (*(*qep)->ops->destroy)(*qep); CHKERRQ(ierr);
1887 jroman 489
  }
2312 jroman 490
  ierr = PetscFree((*qep)->T);CHKERRQ(ierr);
491
  if ((*qep)->eigr) {
492
    ierr = PetscFree((*qep)->eigr);CHKERRQ(ierr);
493
    ierr = PetscFree((*qep)->eigi);CHKERRQ(ierr);
494
    ierr = PetscFree((*qep)->perm);CHKERRQ(ierr);
495
    ierr = PetscFree((*qep)->errest);CHKERRQ(ierr);
496
    ierr = VecGetArray((*qep)->V[0],&pV);CHKERRQ(ierr);
497
    for (i=0;i<(*qep)->ncv;i++) {
498
      ierr = VecDestroy(&(*qep)->V[i]);CHKERRQ(ierr);
1903 jroman 499
    }
500
    ierr = PetscFree(pV);CHKERRQ(ierr);
2312 jroman 501
    ierr = PetscFree((*qep)->V);CHKERRQ(ierr);
1903 jroman 502
  }
2312 jroman 503
  ierr = QEPMonitorCancel(*qep);CHKERRQ(ierr);
504
  ierr = IPDestroy(&(*qep)->ip);CHKERRQ(ierr);
505
  ierr = PetscRandomDestroy(&(*qep)->rand);CHKERRQ(ierr);
506
  ierr = MatDestroy(&(*qep)->M);CHKERRQ(ierr);
507
  ierr = MatDestroy(&(*qep)->C);CHKERRQ(ierr);
508
  ierr = MatDestroy(&(*qep)->K);CHKERRQ(ierr);
509
  ierr = PetscHeaderDestroy(qep);CHKERRQ(ierr);
1887 jroman 510
  PetscFunctionReturn(0);
511
}
512
 
513
#undef __FUNCT__  
514
#define __FUNCT__ "QEPSetIP"
515
/*@
516
   QEPSetIP - Associates an inner product object to the quadratic eigensolver.
517
 
518
   Collective on QEP and IP
519
 
520
   Input Parameters:
521
+  qep - eigensolver context obtained from QEPCreate()
522
-  ip  - the inner product object
523
 
524
   Note:
525
   Use QEPGetIP() to retrieve the inner product context (for example,
526
   to free it at the end of the computations).
527
 
528
   Level: advanced
529
 
530
.seealso: QEPGetIP()
531
@*/
532
PetscErrorCode QEPSetIP(QEP qep,IP ip)
533
{
534
  PetscErrorCode ierr;
535
 
536
  PetscFunctionBegin;
2213 jroman 537
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
538
  PetscValidHeaderSpecific(ip,IP_CLASSID,2);
1887 jroman 539
  PetscCheckSameComm(qep,1,ip,2);
540
  ierr = PetscObjectReference((PetscObject)ip);CHKERRQ(ierr);
2312 jroman 541
  ierr = IPDestroy(&qep->ip); CHKERRQ(ierr);
1887 jroman 542
  qep->ip = ip;
543
  PetscFunctionReturn(0);
544
}
545
 
546
#undef __FUNCT__  
547
#define __FUNCT__ "QEPGetIP"
548
/*@C
549
   QEPGetIP - Obtain the inner product object associated
550
   to the quadratic eigensolver object.
551
 
552
   Not Collective
553
 
554
   Input Parameters:
555
.  qep - eigensolver context obtained from QEPCreate()
556
 
557
   Output Parameter:
558
.  ip - inner product context
559
 
560
   Level: advanced
561
 
562
.seealso: QEPSetIP()
563
@*/
564
PetscErrorCode QEPGetIP(QEP qep,IP *ip)
565
{
566
  PetscFunctionBegin;
2213 jroman 567
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1887 jroman 568
  PetscValidPointer(ip,2);
569
  *ip = qep->ip;
570
  PetscFunctionReturn(0);
571
}
572