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;
2331 jroman 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);
2367 jroman 150
    if (qep->ops->view) {
151
      ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
152
      ierr = (*qep->ops->view)(qep,viewer);CHKERRQ(ierr);
153
      ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
154
    }
2231 jroman 155
    if (qep->problem_type) {
156
      switch (qep->problem_type) {
157
        case QEP_GENERAL:    type = "general quadratic eigenvalue problem"; break;
158
        case QEP_HERMITIAN:  type = HERM " quadratic eigenvalue problem"; break;
159
        case QEP_GYROSCOPIC: type = "gyroscopic quadratic eigenvalue problem"; break;
160
        default: SETERRQ(((PetscObject)qep)->comm,1,"Wrong value of qep->problem_type");
161
      }
162
    } else type = "not yet set";
1907 jroman 163
    ierr = PetscViewerASCIIPrintf(viewer,"  problem type: %s\n",type);CHKERRQ(ierr);
1887 jroman 164
    ierr = PetscViewerASCIIPrintf(viewer,"  selected portion of the spectrum: ");CHKERRQ(ierr);
1942 jroman 165
    if (!qep->which) {
166
      ierr = PetscViewerASCIIPrintf(viewer,"not yet set\n");CHKERRQ(ierr);
167
    } else switch (qep->which) {
1887 jroman 168
      case QEP_LARGEST_MAGNITUDE:
169
        ierr = PetscViewerASCIIPrintf(viewer,"largest eigenvalues in magnitude\n");CHKERRQ(ierr);
170
        break;
171
      case QEP_SMALLEST_MAGNITUDE:
172
        ierr = PetscViewerASCIIPrintf(viewer,"smallest eigenvalues in magnitude\n");CHKERRQ(ierr);
173
        break;
174
      case QEP_LARGEST_REAL:
175
        ierr = PetscViewerASCIIPrintf(viewer,"largest real parts\n");CHKERRQ(ierr);
176
        break;
177
      case QEP_SMALLEST_REAL:
178
        ierr = PetscViewerASCIIPrintf(viewer,"smallest real parts\n");CHKERRQ(ierr);
179
        break;
180
      case QEP_LARGEST_IMAGINARY:
181
        ierr = PetscViewerASCIIPrintf(viewer,"largest imaginary parts\n");CHKERRQ(ierr);
182
        break;
183
      case QEP_SMALLEST_IMAGINARY:
184
        ierr = PetscViewerASCIIPrintf(viewer,"smallest imaginary parts\n");CHKERRQ(ierr);
185
        break;
2214 jroman 186
      default: SETERRQ(((PetscObject)qep)->comm,1,"Wrong value of qep->which");
1887 jroman 187
    }    
1944 jroman 188
    if (qep->leftvecs) {
189
      ierr = PetscViewerASCIIPrintf(viewer,"  computing left eigenvectors also\n");CHKERRQ(ierr);
190
    }
1887 jroman 191
    ierr = PetscViewerASCIIPrintf(viewer,"  number of eigenvalues (nev): %d\n",qep->nev);CHKERRQ(ierr);
192
    ierr = PetscViewerASCIIPrintf(viewer,"  number of column vectors (ncv): %d\n",qep->ncv);CHKERRQ(ierr);
193
    ierr = PetscViewerASCIIPrintf(viewer,"  maximum dimension of projected problem (mpd): %d\n",qep->mpd);CHKERRQ(ierr);
2331 jroman 194
    ierr = PetscViewerASCIIPrintf(viewer,"  maximum number of iterations: %d\n",qep->max_it);
1887 jroman 195
    ierr = PetscViewerASCIIPrintf(viewer,"  tolerance: %g\n",qep->tol);CHKERRQ(ierr);
1918 jroman 196
    ierr = PetscViewerASCIIPrintf(viewer,"  scaling factor: %g\n",qep->sfactor);CHKERRQ(ierr);
1952 jroman 197
    if (qep->nini!=0) {
198
      ierr = PetscViewerASCIIPrintf(viewer,"  dimension of user-provided initial space: %d\n",PetscAbs(qep->nini));CHKERRQ(ierr);
199
    }
200
    if (qep->ninil!=0) {
201
      ierr = PetscViewerASCIIPrintf(viewer,"  dimension of user-provided initial left space: %d\n",PetscAbs(qep->ninil));CHKERRQ(ierr);
202
    }
1887 jroman 203
  } else {
204
    if (qep->ops->view) {
205
      ierr = (*qep->ops->view)(qep,viewer);CHKERRQ(ierr);
206
    }
207
  }
2367 jroman 208
  ierr = IPView(qep->ip,viewer);CHKERRQ(ierr);
1887 jroman 209
  PetscFunctionReturn(0);
210
}
211
 
212
#undef __FUNCT__  
213
#define __FUNCT__ "QEPCreate"
214
/*@C
215
   QEPCreate - Creates the default QEP context.
216
 
217
   Collective on MPI_Comm
218
 
219
   Input Parameter:
220
.  comm - MPI communicator
221
 
222
   Output Parameter:
223
.  qep - location to put the QEP context
224
 
225
   Note:
226
   The default QEP type is QEPLINEAR
227
 
228
   Level: beginner
229
 
230
.seealso: QEPSetUp(), QEPSolve(), QEPDestroy(), QEP
231
@*/
232
PetscErrorCode QEPCreate(MPI_Comm comm,QEP *outqep)
233
{
234
  PetscErrorCode ierr;
235
  QEP            qep;
236
 
237
  PetscFunctionBegin;
238
  PetscValidPointer(outqep,2);
239
  *outqep = 0;
2213 jroman 240
  ierr = PetscHeaderCreate(qep,_p_QEP,struct _QEPOps,QEP_CLASSID,-1,"QEP",comm,QEPDestroy,QEPView);CHKERRQ(ierr);
1887 jroman 241
 
2309 jroman 242
  qep->M               = 0;
243
  qep->C               = 0;
244
  qep->K               = 0;
1887 jroman 245
  qep->max_it          = 0;
246
  qep->nev             = 1;
247
  qep->ncv             = 0;
248
  qep->mpd             = 0;
1952 jroman 249
  qep->nini            = 0;
250
  qep->ninil           = 0;
2369 jroman 251
  qep->allocated_ncv   = 0;
1887 jroman 252
  qep->tol             = 1e-7;
1918 jroman 253
  qep->sfactor         = 0.0;
1887 jroman 254
  qep->conv_func       = QEPDefaultConverged;
255
  qep->conv_ctx        = PETSC_NULL;
1942 jroman 256
  qep->which           = (QEPWhich)0;
1887 jroman 257
  qep->which_func      = PETSC_NULL;
258
  qep->which_ctx       = PETSC_NULL;
1944 jroman 259
  qep->leftvecs        = PETSC_FALSE;
1906 jroman 260
  qep->problem_type    = (QEPProblemType)0;
1887 jroman 261
  qep->V               = PETSC_NULL;
2309 jroman 262
  qep->W               = PETSC_NULL;
1952 jroman 263
  qep->IS              = PETSC_NULL;
264
  qep->ISL             = PETSC_NULL;
1887 jroman 265
  qep->T               = PETSC_NULL;
266
  qep->eigr            = PETSC_NULL;
267
  qep->eigi            = PETSC_NULL;
268
  qep->errest          = PETSC_NULL;
269
  qep->data            = PETSC_NULL;
270
  qep->nconv           = 0;
271
  qep->its             = 0;
272
  qep->perm            = PETSC_NULL;
1896 jroman 273
  qep->matvecs         = 0;
274
  qep->linits          = 0;
1887 jroman 275
  qep->nwork           = 0;
276
  qep->work            = PETSC_NULL;
277
  qep->setupcalled     = 0;
278
  qep->reason          = QEP_CONVERGED_ITERATING;
279
  qep->numbermonitors  = 0;
2054 eromero 280
  qep->trackall        = PETSC_FALSE;
2309 jroman 281
  qep->rand            = 0;
1887 jroman 282
 
2361 jroman 283
  ierr = PetscRandomCreate(comm,&qep->rand);CHKERRQ(ierr);
284
  ierr = PetscLogObjectParent(qep,qep->rand);CHKERRQ(ierr);
2330 jroman 285
  ierr = IPCreate(comm,&qep->ip);CHKERRQ(ierr);
2327 jroman 286
  ierr = PetscLogObjectParent(qep,qep->ip);CHKERRQ(ierr);
2361 jroman 287
  *outqep = qep;
1887 jroman 288
  PetscFunctionReturn(0);
289
}
290
 
291
#undef __FUNCT__  
292
#define __FUNCT__ "QEPSetType"
293
/*@C
294
   QEPSetType - Selects the particular solver to be used in the QEP object.
295
 
2328 jroman 296
   Logically Collective on QEP
1887 jroman 297
 
298
   Input Parameters:
299
+  qep      - the quadratic eigensolver context
300
-  type     - a known method
301
 
302
   Options Database Key:
303
.  -qep_type <method> - Sets the method; use -help for a list
304
    of available methods
305
 
306
   Notes:  
307
   See "slepc/include/slepcqep.h" for available methods. The default
308
   is QEPLINEAR.
309
 
310
   Normally, it is best to use the QEPSetFromOptions() command and
311
   then set the QEP type from the options database rather than by using
312
   this routine.  Using the options database provides the user with
313
   maximum flexibility in evaluating the different available methods.
314
   The QEPSetType() routine is provided for those situations where it
315
   is necessary to set the iterative solver independently of the command
316
   line or options database.
317
 
318
   Level: intermediate
319
 
320
.seealso: QEPType
321
@*/
322
PetscErrorCode QEPSetType(QEP qep,const QEPType type)
323
{
324
  PetscErrorCode ierr,(*r)(QEP);
2216 jroman 325
  PetscBool      match;
1887 jroman 326
 
327
  PetscFunctionBegin;
2213 jroman 328
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1887 jroman 329
  PetscValidCharPointer(type,2);
2361 jroman 330
 
1887 jroman 331
  ierr = PetscTypeCompare((PetscObject)qep,type,&match);CHKERRQ(ierr);
332
  if (match) PetscFunctionReturn(0);
333
 
2300 jroman 334
  ierr = PetscFListFind(QEPList,((PetscObject)qep)->comm,type,PETSC_TRUE,(void (**)(void))&r);CHKERRQ(ierr);
2361 jroman 335
  if (!r) SETERRQ1(((PetscObject)qep)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown QEP type given: %s",type);
1887 jroman 336
 
2361 jroman 337
  if (qep->ops->destroy) { ierr = (*qep->ops->destroy)(qep);CHKERRQ(ierr); }
338
  ierr = PetscMemzero(qep->ops,sizeof(struct _QEPOps));CHKERRQ(ierr);
1887 jroman 339
 
340
  qep->setupcalled = 0;
2361 jroman 341
  ierr = PetscObjectChangeTypeName((PetscObject)qep,type);CHKERRQ(ierr);
1887 jroman 342
  ierr = (*r)(qep);CHKERRQ(ierr);
343
  PetscFunctionReturn(0);
344
}
345
 
346
#undef __FUNCT__  
347
#define __FUNCT__ "QEPGetType"
348
/*@C
349
   QEPGetType - Gets the QEP type as a string from the QEP object.
350
 
351
   Not Collective
352
 
353
   Input Parameter:
354
.  qep - the eigensolver context
355
 
356
   Output Parameter:
357
.  name - name of QEP method
358
 
359
   Level: intermediate
360
 
361
.seealso: QEPSetType()
362
@*/
363
PetscErrorCode QEPGetType(QEP qep,const QEPType *type)
364
{
365
  PetscFunctionBegin;
2213 jroman 366
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1887 jroman 367
  PetscValidPointer(type,2);
368
  *type = ((PetscObject)qep)->type_name;
369
  PetscFunctionReturn(0);
370
}
371
 
372
/*MC
373
   QEPRegisterDynamic - Adds a method to the quadratic eigenproblem solver package.
374
 
375
   Synopsis:
376
   QEPRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(QEP))
377
 
378
   Not Collective
379
 
380
   Input Parameters:
381
+  name_solver - name of a new user-defined solver
382
.  path - path (either absolute or relative) the library containing this solver
383
.  name_create - name of routine to create the solver context
384
-  routine_create - routine to create the solver context
385
 
386
   Notes:
387
   QEPRegisterDynamic() may be called multiple times to add several user-defined solvers.
388
 
389
   If dynamic libraries are used, then the fourth input argument (routine_create)
390
   is ignored.
391
 
392
   Sample usage:
393
.vb
394
   QEPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
395
               "MySolverCreate",MySolverCreate);
396
.ve
397
 
398
   Then, your solver can be chosen with the procedural interface via
399
$     QEPSetType(qep,"my_solver")
400
   or at runtime via the option
401
$     -qep_type my_solver
402
 
403
   Level: advanced
404
 
405
   Environmental variables such as ${PETSC_ARCH}, ${SLEPC_DIR},
406
   and others of the form ${any_environmental_variable} occuring in pathname will be
407
   replaced with appropriate values.
408
 
409
.seealso: QEPRegisterDestroy(), QEPRegisterAll()
410
 
411
M*/
412
 
413
#undef __FUNCT__  
414
#define __FUNCT__ "QEPRegister"
415
/*@C
416
  QEPRegister - See QEPRegisterDynamic()
417
 
418
  Level: advanced
419
@*/
420
PetscErrorCode QEPRegister(const char *sname,const char *path,const char *name,PetscErrorCode (*function)(QEP))
421
{
422
  PetscErrorCode ierr;
2374 jroman 423
  char           fullname[PETSC_MAX_PATH_LEN];
1887 jroman 424
 
425
  PetscFunctionBegin;
426
  ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
427
  ierr = PetscFListAdd(&QEPList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
428
  PetscFunctionReturn(0);
429
}
430
 
431
#undef __FUNCT__  
432
#define __FUNCT__ "QEPRegisterDestroy"
433
/*@
434
   QEPRegisterDestroy - Frees the list of QEP methods that were
435
   registered by QEPRegisterDynamic().
436
 
437
   Not Collective
438
 
439
   Level: advanced
440
 
441
.seealso: QEPRegisterDynamic(), QEPRegisterAll()
442
@*/
443
PetscErrorCode QEPRegisterDestroy(void)
444
{
445
  PetscErrorCode ierr;
446
 
447
  PetscFunctionBegin;
448
  ierr = PetscFListDestroy(&QEPList);CHKERRQ(ierr);
449
  ierr = QEPRegisterAll(PETSC_NULL);CHKERRQ(ierr);
450
  PetscFunctionReturn(0);
451
}
452
 
453
#undef __FUNCT__  
2361 jroman 454
#define __FUNCT__ "QEPReset"
455
/*@C
456
   QEPReset - Resets the QEP context to the setupcalled=0 state and removes any
457
   allocated objects.
458
 
459
   Collective on QEP
460
 
461
   Input Parameter:
462
.  qep - eigensolver context obtained from QEPCreate()
463
 
464
   Level: advanced
465
 
466
.seealso: QEPDestroy()
467
@*/
468
PetscErrorCode QEPReset(QEP qep)
469
{
470
  PetscErrorCode ierr;
471
 
472
  PetscFunctionBegin;
473
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
474
  if (qep->ops->reset) { ierr = (qep->ops->reset)(qep);CHKERRQ(ierr); }
475
  if (qep->ip) { ierr = IPReset(qep->ip);CHKERRQ(ierr); }
476
  ierr = MatDestroy(&qep->M);CHKERRQ(ierr);
477
  ierr = MatDestroy(&qep->C);CHKERRQ(ierr);
478
  ierr = MatDestroy(&qep->K);CHKERRQ(ierr);
2369 jroman 479
  ierr = QEPFreeSolution(qep);CHKERRQ(ierr);
2361 jroman 480
  qep->matvecs     = 0;
481
  qep->linits      = 0;
482
  qep->setupcalled = 0;  
483
  PetscFunctionReturn(0);
484
}
485
 
486
#undef __FUNCT__  
1887 jroman 487
#define __FUNCT__ "QEPDestroy"
2312 jroman 488
/*@C
1887 jroman 489
   QEPDestroy - Destroys the QEP context.
490
 
491
   Collective on QEP
492
 
493
   Input Parameter:
494
.  qep - eigensolver context obtained from QEPCreate()
495
 
496
   Level: beginner
497
 
498
.seealso: QEPCreate(), QEPSetUp(), QEPSolve()
499
@*/
2312 jroman 500
PetscErrorCode QEPDestroy(QEP *qep)
1887 jroman 501
{
502
  PetscErrorCode ierr;
503
 
504
  PetscFunctionBegin;
2312 jroman 505
  if (!*qep) PetscFunctionReturn(0);
506
  PetscValidHeaderSpecific(*qep,QEP_CLASSID,1);
507
  if (--((PetscObject)(*qep))->refct > 0) { *qep = 0; PetscFunctionReturn(0); }
2361 jroman 508
  ierr = QEPReset(*qep);CHKERRQ(ierr);
2312 jroman 509
  ierr = PetscObjectDepublish(*qep);CHKERRQ(ierr);
2361 jroman 510
  if ((*qep)->ops->destroy) { ierr = (*(*qep)->ops->destroy)(*qep);CHKERRQ(ierr); }
2312 jroman 511
  ierr = IPDestroy(&(*qep)->ip);CHKERRQ(ierr);
512
  ierr = PetscRandomDestroy(&(*qep)->rand);CHKERRQ(ierr);
2361 jroman 513
  ierr = QEPMonitorCancel(*qep);CHKERRQ(ierr);
2312 jroman 514
  ierr = PetscHeaderDestroy(qep);CHKERRQ(ierr);
1887 jroman 515
  PetscFunctionReturn(0);
516
}
517
 
518
#undef __FUNCT__  
519
#define __FUNCT__ "QEPSetIP"
520
/*@
521
   QEPSetIP - Associates an inner product object to the quadratic eigensolver.
522
 
2328 jroman 523
   Collective on QEP
1887 jroman 524
 
525
   Input Parameters:
526
+  qep - eigensolver context obtained from QEPCreate()
527
-  ip  - the inner product object
528
 
529
   Note:
530
   Use QEPGetIP() to retrieve the inner product context (for example,
531
   to free it at the end of the computations).
532
 
533
   Level: advanced
534
 
535
.seealso: QEPGetIP()
536
@*/
537
PetscErrorCode QEPSetIP(QEP qep,IP ip)
538
{
539
  PetscErrorCode ierr;
540
 
541
  PetscFunctionBegin;
2213 jroman 542
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
543
  PetscValidHeaderSpecific(ip,IP_CLASSID,2);
1887 jroman 544
  PetscCheckSameComm(qep,1,ip,2);
545
  ierr = PetscObjectReference((PetscObject)ip);CHKERRQ(ierr);
2330 jroman 546
  ierr = IPDestroy(&qep->ip);CHKERRQ(ierr);
1887 jroman 547
  qep->ip = ip;
2327 jroman 548
  ierr = PetscLogObjectParent(qep,qep->ip);CHKERRQ(ierr);
1887 jroman 549
  PetscFunctionReturn(0);
550
}
551
 
552
#undef __FUNCT__  
553
#define __FUNCT__ "QEPGetIP"
554
/*@C
555
   QEPGetIP - Obtain the inner product object associated
556
   to the quadratic eigensolver object.
557
 
558
   Not Collective
559
 
560
   Input Parameters:
561
.  qep - eigensolver context obtained from QEPCreate()
562
 
563
   Output Parameter:
564
.  ip - inner product context
565
 
566
   Level: advanced
567
 
568
.seealso: QEPSetIP()
569
@*/
570
PetscErrorCode QEPGetIP(QEP qep,IP *ip)
571
{
572
  PetscFunctionBegin;
2213 jroman 573
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1887 jroman 574
  PetscValidPointer(ip,2);
575
  *ip = qep->ip;
576
  PetscFunctionReturn(0);
577
}
578