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