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
519 dsic.upv.es!antodo 1
/*
2
     The basic EPS routines, Create, View, etc. are here.
3
*/
4
#include "src/eps/epsimpl.h"      /*I "slepceps.h" I*/
5
 
6
PetscFList EPSList = 0;
842 dsic.upv.es!antodo 7
PetscCookie EPS_COOKIE = 0;
870 dsic.upv.es!antodo 8
PetscEvent EPS_SetUp = 0, EPS_Solve = 0, EPS_Orthogonalize = 0;
519 dsic.upv.es!antodo 9
 
10
#undef __FUNCT__  
842 dsic.upv.es!antodo 11
#define __FUNCT__ "EPSInitializePackage"
12
/*@C
13
  EPSInitializePackage - This function initializes everything in the EPS package. It is called
14
  from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to EPSCreate()
15
  when using static libraries.
16
 
17
  Input Parameter:
18
  path - The dynamic library path, or PETSC_NULL
19
 
20
  Level: developer
21
 
22
.seealso: SlepcInitialize()
23
@*/
24
PetscErrorCode EPSInitializePackage(char *path) {
25
  static PetscTruth initialized = PETSC_FALSE;
26
  char              logList[256];
27
  char             *className;
28
  PetscTruth        opt;
29
  PetscErrorCode ierr;
30
 
31
  PetscFunctionBegin;
32
  if (initialized) PetscFunctionReturn(0);
33
  initialized = PETSC_TRUE;
34
  /* Register Classes */
35
  ierr = PetscLogClassRegister(&EPS_COOKIE,"Eigenproblem Solver");CHKERRQ(ierr);
36
  /* Register Constructors */
37
  ierr = EPSRegisterAll(path);CHKERRQ(ierr);
38
  /* Register Events */
39
  ierr = PetscLogEventRegister(&EPS_SetUp,"EPSSetUp",EPS_COOKIE);CHKERRQ(ierr);
40
  ierr = PetscLogEventRegister(&EPS_Solve,"EPSSolve",EPS_COOKIE);CHKERRQ(ierr);
41
  ierr = PetscLogEventRegister(&EPS_Orthogonalize,"EPSOrthogonalize",EPS_COOKIE); CHKERRQ(ierr);
42
  /* Process info exclusions */
43
  ierr = PetscOptionsGetString(PETSC_NULL, "-log_info_exclude", logList, 256, &opt);CHKERRQ(ierr);
44
  if (opt) {
45
    ierr = PetscStrstr(logList, "eps", &className);CHKERRQ(ierr);
46
    if (className) {
1037 slepc 47
      ierr = PetscInfoDeactivateClass(EPS_COOKIE);CHKERRQ(ierr);
842 dsic.upv.es!antodo 48
    }
49
  }
50
  /* Process summary exclusions */
51
  ierr = PetscOptionsGetString(PETSC_NULL, "-log_summary_exclude", logList, 256, &opt);CHKERRQ(ierr);
52
  if (opt) {
53
    ierr = PetscStrstr(logList, "eps", &className);CHKERRQ(ierr);
54
    if (className) {
55
      ierr = PetscLogEventDeactivateClass(EPS_COOKIE);CHKERRQ(ierr);
56
    }
57
  }
58
  PetscFunctionReturn(0);
59
}
60
 
61
#undef __FUNCT__  
519 dsic.upv.es!antodo 62
#define __FUNCT__ "EPSView"
63
/*@C
64
   EPSView - Prints the EPS data structure.
65
 
66
   Collective on EPS
67
 
68
   Input Parameters:
69
+  eps - the eigenproblem solver context
70
-  viewer - optional visualization context
71
 
72
   Options Database Key:
73
.  -eps_view -  Calls EPSView() at end of EPSSolve()
74
 
75
   Note:
76
   The available visualization contexts include
77
+     PETSC_VIEWER_STDOUT_SELF - standard output (default)
78
-     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
79
         output where only the first processor opens
80
         the file.  All other processors send their
81
         data to the first processor to print.
82
 
83
   The user can open an alternative visualization context with
84
   PetscViewerASCIIOpen() - output to a specified file.
85
 
86
   Level: beginner
87
 
88
.seealso: STView(), PetscViewerASCIIOpen()
89
@*/
90
PetscErrorCode EPSView(EPS eps,PetscViewer viewer)
91
{
92
  PetscErrorCode ierr;
1004 slepc 93
  const char     *type, *which;
94
  PetscTruth     isascii;
519 dsic.upv.es!antodo 95
 
96
  PetscFunctionBegin;
97
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
98
  if (!viewer) viewer = PETSC_VIEWER_STDOUT_(eps->comm);
99
  PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2);
100
  PetscCheckSameComm(eps,1,viewer,2);
101
 
102
#if defined(PETSC_USE_COMPLEX)
103
#define HERM "hermitian"
104
#else
105
#define HERM "symmetric"
106
#endif
107
  ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
108
  if (isascii) {
109
    ierr = PetscViewerASCIIPrintf(viewer,"EPS Object:\n");CHKERRQ(ierr);
110
    switch (eps->problem_type) {
111
      case EPS_HEP:   type = HERM " eigenvalue problem"; break;
112
      case EPS_GHEP:  type = "generalized " HERM " eigenvalue problem"; break;
113
      case EPS_NHEP:  type = "non-" HERM " eigenvalue problem"; break;
114
      case EPS_GNHEP: type = "generalized non-" HERM " eigenvalue problem"; break;
115
      case 0:         type = "not yet set"; break;
116
      default: SETERRQ(1,"Wrong value of eps->problem_type");
117
    }
118
    ierr = PetscViewerASCIIPrintf(viewer,"  problem type: %s\n",type);CHKERRQ(ierr);
119
    ierr = EPSGetType(eps,&type);CHKERRQ(ierr);
120
    if (type) {
780 dsic.upv.es!jroman 121
      ierr = PetscViewerASCIIPrintf(viewer,"  method: %s",type);CHKERRQ(ierr);
122
      switch (eps->solverclass) {
123
        case EPS_ONE_SIDE:
124
          ierr = PetscViewerASCIIPrintf(viewer,"\n",type);CHKERRQ(ierr); break;
125
        case EPS_TWO_SIDE:
126
          ierr = PetscViewerASCIIPrintf(viewer," (two-sided)\n",type);CHKERRQ(ierr); break;
127
        default: SETERRQ(1,"Wrong value of eps->solverclass");
128
      }
519 dsic.upv.es!antodo 129
    } else {
130
      ierr = PetscViewerASCIIPrintf(viewer,"  method: not yet set\n");CHKERRQ(ierr);
131
    }
132
    if (eps->ops->view) {
133
      ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
134
      ierr = (*eps->ops->view)(eps,viewer);CHKERRQ(ierr);
135
      ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
136
    }
137
    switch (eps->which) {
138
      case EPS_LARGEST_MAGNITUDE:  which = "largest eigenvalues in magnitude"; break;
139
      case EPS_SMALLEST_MAGNITUDE: which = "smallest eigenvalues in magnitude"; break;
140
      case EPS_LARGEST_REAL:       which = "largest real parts"; break;
141
      case EPS_SMALLEST_REAL:      which = "smallest real parts"; break;
142
      case EPS_LARGEST_IMAGINARY:  which = "largest imaginary parts"; break;
143
      case EPS_SMALLEST_IMAGINARY: which = "smallest imaginary parts"; break;
144
      default: SETERRQ(1,"Wrong value of eps->which");
145
    }
146
    ierr = PetscViewerASCIIPrintf(viewer,"  selected portion of the spectrum: %s\n",which);CHKERRQ(ierr);
147
    ierr = PetscViewerASCIIPrintf(viewer,"  number of eigenvalues (nev): %d\n",eps->nev);CHKERRQ(ierr);
600 dsic.upv.es!jroman 148
    ierr = PetscViewerASCIIPrintf(viewer,"  number of column vectors (ncv): %d\n",eps->ncv);CHKERRQ(ierr);
519 dsic.upv.es!antodo 149
    ierr = PetscViewerASCIIPrintf(viewer,"  maximum number of iterations: %d\n", eps->max_it);
150
    ierr = PetscViewerASCIIPrintf(viewer,"  tolerance: %g\n",eps->tol);CHKERRQ(ierr);
151
    ierr = PetscViewerASCIIPrintf(viewer,"  orthogonalization method: ");CHKERRQ(ierr);
152
    switch (eps->orthog_type) {
153
      case EPS_MGS_ORTH:
154
        ierr = PetscViewerASCIIPrintf(viewer,"modified Gram-Schmidt\n");CHKERRQ(ierr);
155
        break;
156
      case EPS_CGS_ORTH:
157
        ierr = PetscViewerASCIIPrintf(viewer,"classical Gram-Schmidt\n");CHKERRQ(ierr);
158
        break;
159
      default: SETERRQ(1,"Wrong value of eps->orth_type");
160
    }
161
    ierr = PetscViewerASCIIPrintf(viewer,"  orthogonalization refinement: ");CHKERRQ(ierr);
162
    switch (eps->orthog_ref) {
163
      case EPS_ORTH_REFINE_NEVER:
164
        ierr = PetscViewerASCIIPrintf(viewer,"never\n");CHKERRQ(ierr);
165
        break;
166
      case EPS_ORTH_REFINE_IFNEEDED:
167
        ierr = PetscViewerASCIIPrintf(viewer,"if needed (eta: %f)\n",eps->orthog_eta);CHKERRQ(ierr);
168
        break;
169
      case EPS_ORTH_REFINE_ALWAYS:
170
        ierr = PetscViewerASCIIPrintf(viewer,"always\n");CHKERRQ(ierr);
171
        break;
172
      default: SETERRQ(1,"Wrong value of eps->orth_ref");
173
    }
174
    ierr = PetscViewerASCIIPrintf(viewer,"  dimension of user-provided deflation space: %d\n",eps->nds);CHKERRQ(ierr);
175
    ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
176
    ierr = STView(eps->OP,viewer); CHKERRQ(ierr);
177
    ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
178
  } else {
179
    if (eps->ops->view) {
180
      ierr = (*eps->ops->view)(eps,viewer);CHKERRQ(ierr);
181
    }
182
    ierr = STView(eps->OP,viewer); CHKERRQ(ierr);
183
  }
184
  PetscFunctionReturn(0);
185
}
186
 
187
#undef __FUNCT__  
188
#define __FUNCT__ "EPSPublish_Petsc"
189
static PetscErrorCode EPSPublish_Petsc(PetscObject object)
190
{
191
  PetscFunctionBegin;
192
  PetscFunctionReturn(0);
193
}
194
 
195
#undef __FUNCT__  
196
#define __FUNCT__ "EPSCreate"
197
/*@C
198
   EPSCreate - Creates the default EPS context.
199
 
200
   Collective on MPI_Comm
201
 
202
   Input Parameter:
203
.  comm - MPI communicator
204
 
205
   Output Parameter:
206
.  eps - location to put the EPS context
207
 
208
   Note:
1198 slepc 209
   The default EPS type is EPSKRYLOVSCHUR
519 dsic.upv.es!antodo 210
 
211
   Level: beginner
212
 
213
.seealso: EPSSetUp(), EPSSolve(), EPSDestroy(), EPS
214
@*/
215
PetscErrorCode EPSCreate(MPI_Comm comm,EPS *outeps)
216
{
217
  PetscErrorCode ierr;
218
  EPS            eps;
219
 
220
  PetscFunctionBegin;
221
  PetscValidPointer(outeps,2);
222
  *outeps = 0;
223
 
224
  PetscHeaderCreate(eps,_p_EPS,struct _EPSOps,EPS_COOKIE,-1,"EPS",comm,EPSDestroy,EPSView);
225
  PetscLogObjectCreate(eps);
226
  *outeps = eps;
227
 
228
  eps->bops->publish   = EPSPublish_Petsc;
229
  ierr = PetscMemzero(eps->ops,sizeof(struct _EPSOps));CHKERRQ(ierr);
230
 
231
  eps->type            = -1;
232
  eps->max_it          = 0;
233
  eps->nev             = 1;
234
  eps->ncv             = 0;
235
  eps->allocated_ncv   = 0;
236
  eps->nds             = 0;
237
  eps->tol             = 0.0;
238
  eps->which           = EPS_LARGEST_MAGNITUDE;
239
  eps->evecsavailable  = PETSC_FALSE;
240
  eps->problem_type    = (EPSProblemType)0;
780 dsic.upv.es!jroman 241
  eps->solverclass     = (EPSClass)0;
519 dsic.upv.es!antodo 242
 
243
  eps->vec_initial     = 0;
780 dsic.upv.es!jroman 244
  eps->vec_initial_left= 0;
519 dsic.upv.es!antodo 245
  eps->V               = 0;
246
  eps->AV              = 0;
780 dsic.upv.es!jroman 247
  eps->W               = 0;
248
  eps->AW              = 0;
519 dsic.upv.es!antodo 249
  eps->T               = 0;
250
  eps->DS              = 0;
251
  eps->ds_ortho        = PETSC_TRUE;
252
  eps->eigr            = 0;
253
  eps->eigi            = 0;
254
  eps->errest          = 0;
780 dsic.upv.es!jroman 255
  eps->errest_left     = 0;
519 dsic.upv.es!antodo 256
  eps->OP              = 0;
257
  eps->data            = 0;
258
  eps->nconv           = 0;
259
  eps->its             = 0;
260
  eps->perm            = PETSC_NULL;
261
 
262
  eps->nwork           = 0;
263
  eps->work            = 0;
264
  eps->isgeneralized   = PETSC_FALSE;
265
  eps->ishermitian     = PETSC_FALSE;
266
  eps->setupcalled     = 0;
267
  eps->reason          = EPS_CONVERGED_ITERATING;
268
 
269
  eps->numbermonitors  = 0;
270
 
271
  eps->orthog_type     = EPS_CGS_ORTH;
272
  eps->orthog_ref      = EPS_ORTH_REFINE_IFNEEDED;
691 dsic.upv.es!antodo 273
  eps->orthog_eta      = PETSC_DEFAULT;
519 dsic.upv.es!antodo 274
 
275
  ierr = STCreate(comm,&eps->OP); CHKERRQ(ierr);
276
  PetscLogObjectParent(eps,eps->OP);
277
  ierr = PetscPublishAll(eps);CHKERRQ(ierr);
278
  PetscFunctionReturn(0);
279
}
280
 
281
#undef __FUNCT__  
282
#define __FUNCT__ "EPSSetType"
283
/*@C
284
   EPSSetType - Selects the particular solver to be used in the EPS object.
285
 
286
   Collective on EPS
287
 
288
   Input Parameters:
289
+  eps      - the eigensolver context
290
-  type     - a known method
291
 
292
   Options Database Key:
293
.  -eps_type <method> - Sets the method; use -help for a list
294
    of available methods
295
 
296
   Notes:  
545 dsic.upv.es!jroman 297
   See "slepc/include/slepceps.h" for available methods. The default
1198 slepc 298
   is EPSKRYLOVSCHUR.
519 dsic.upv.es!antodo 299
 
300
   Normally, it is best to use the EPSSetFromOptions() command and
301
   then set the EPS type from the options database rather than by using
302
   this routine.  Using the options database provides the user with
303
   maximum flexibility in evaluating the different available methods.
304
   The EPSSetType() routine is provided for those situations where it
305
   is necessary to set the iterative solver independently of the command
306
   line or options database.
307
 
308
   Level: intermediate
309
 
310
.seealso: STSetType(), EPSType
311
@*/
312
PetscErrorCode EPSSetType(EPS eps,EPSType type)
313
{
314
  PetscErrorCode ierr,(*r)(EPS);
315
  PetscTruth match;
316
 
317
  PetscFunctionBegin;
318
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
319
  PetscValidCharPointer(type,2);
320
 
321
  ierr = PetscTypeCompare((PetscObject)eps,type,&match);CHKERRQ(ierr);
322
  if (match) PetscFunctionReturn(0);
323
 
324
  if (eps->data) {
325
    /* destroy the old private EPS context */
326
    ierr = (*eps->ops->destroy)(eps); CHKERRQ(ierr);
327
    eps->data = 0;
328
  }
329
 
330
  ierr = PetscFListFind(eps->comm,EPSList,type,(void (**)(void)) &r);CHKERRQ(ierr);
331
 
332
  if (!r) SETERRQ1(1,"Unknown EPS type given: %s",type);
333
 
334
  eps->setupcalled = 0;
335
  ierr = PetscMemzero(eps->ops,sizeof(struct _EPSOps));CHKERRQ(ierr);
336
  ierr = (*r)(eps); CHKERRQ(ierr);
337
 
338
  ierr = PetscObjectChangeTypeName((PetscObject)eps,type);CHKERRQ(ierr);
339
  PetscFunctionReturn(0);
340
}
341
 
342
#undef __FUNCT__  
343
#define __FUNCT__ "EPSGetType"
344
/*@C
345
   EPSGetType - Gets the EPS type as a string from the EPS object.
346
 
347
   Not Collective
348
 
349
   Input Parameter:
350
.  eps - the eigensolver context
351
 
352
   Output Parameter:
353
.  name - name of EPS method
354
 
355
   Level: intermediate
356
 
357
.seealso: EPSSetType()
358
@*/
359
PetscErrorCode EPSGetType(EPS eps,EPSType *type)
360
{
361
  PetscFunctionBegin;
362
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
363
  *type = eps->type_name;
364
  PetscFunctionReturn(0);
365
}
366
 
367
/*MC
368
   EPSRegisterDynamic - Adds a method to the eigenproblem solver package.
369
 
370
   Synopsis:
371
   EPSRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(EPS))
372
 
373
   Not Collective
374
 
375
   Input Parameters:
376
+  name_solver - name of a new user-defined solver
377
.  path - path (either absolute or relative) the library containing this solver
378
.  name_create - name of routine to create the solver context
379
-  routine_create - routine to create the solver context
380
 
381
   Notes:
382
   EPSRegisterDynamic() may be called multiple times to add several user-defined solvers.
383
 
384
   If dynamic libraries are used, then the fourth input argument (routine_create)
385
   is ignored.
386
 
387
   Sample usage:
388
.vb
389
   EPSRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
390
               "MySolverCreate",MySolverCreate);
391
.ve
392
 
393
   Then, your solver can be chosen with the procedural interface via
394
$     EPSSetType(eps,"my_solver")
395
   or at runtime via the option
396
$     -eps_type my_solver
397
 
398
   Level: advanced
399
 
400
   Environmental variables such as ${PETSC_ARCH}, ${SLEPC_DIR}, ${BOPT},
401
   and others of the form ${any_environmental_variable} occuring in pathname will be
402
   replaced with appropriate values.
403
 
842 dsic.upv.es!antodo 404
.seealso: EPSRegisterAll()
519 dsic.upv.es!antodo 405
 
406
M*/
407
 
408
#undef __FUNCT__  
409
#define __FUNCT__ "EPSRegister"
1004 slepc 410
PetscErrorCode EPSRegister(const char *sname,const char *path,const char *name,int (*function)(EPS))
519 dsic.upv.es!antodo 411
{
412
  PetscErrorCode ierr;
413
  char           fullname[256];
414
 
415
  PetscFunctionBegin;
416
  ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
417
  ierr = PetscFListAdd(&EPSList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
418
  PetscFunctionReturn(0);
419
}
420
 
421
#undef __FUNCT__  
422
#define __FUNCT__ "EPSDestroy"
423
/*@
424
   EPSDestroy - Destroys the EPS context.
425
 
426
   Collective on EPS
427
 
428
   Input Parameter:
429
.  eps - eigensolver context obtained from EPSCreate()
430
 
431
   Level: beginner
432
 
433
.seealso: EPSCreate(), EPSSetUp(), EPSSolve()
434
@*/
435
PetscErrorCode EPSDestroy(EPS eps)
436
{
437
  PetscErrorCode ierr;
438
 
439
  PetscFunctionBegin;
440
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
441
  if (--eps->refct > 0) PetscFunctionReturn(0);
442
 
443
  /* if memory was published with AMS then destroy it */
444
  ierr = PetscObjectDepublish(eps);CHKERRQ(ierr);
445
 
446
  ierr = STDestroy(eps->OP);CHKERRQ(ierr);
447
 
448
  if (eps->ops->destroy) {
449
    ierr = (*eps->ops->destroy)(eps); CHKERRQ(ierr);
450
  }
451
 
1040 slepc 452
  ierr = PetscFree(eps->T);CHKERRQ(ierr);
453
  ierr = PetscFree(eps->Tl);CHKERRQ(ierr);
454
  ierr = PetscFree(eps->perm);CHKERRQ(ierr);
519 dsic.upv.es!antodo 455
 
456
  if (eps->vec_initial) {
457
    ierr = VecDestroy(eps->vec_initial);CHKERRQ(ierr);
458
  }
459
 
780 dsic.upv.es!jroman 460
  if (eps->vec_initial_left) {
461
    ierr = VecDestroy(eps->vec_initial_left);CHKERRQ(ierr);
462
  }
463
 
519 dsic.upv.es!antodo 464
  if (eps->nds > 0) {
465
    ierr = VecDestroyVecs(eps->DS, eps->nds);
466
  }
467
 
1040 slepc 468
  ierr = PetscFree(eps->DSV);CHKERRQ(ierr);
519 dsic.upv.es!antodo 469
 
470
  PetscLogObjectDestroy(eps);
471
  PetscHeaderDestroy(eps);
472
  PetscFunctionReturn(0);
473
}
474
 
475
#undef __FUNCT__  
476
#define __FUNCT__ "EPSSetST"
477
/*@
478
   EPSSetST - Associates a spectral transformation object to the
479
   eigensolver.
480
 
481
   Collective on EPS
482
 
483
   Input Parameters:
484
+  eps - eigensolver context obtained from EPSCreate()
485
-  st   - the spectral transformation object
486
 
487
   Note:
488
   Use EPSGetST() to retrieve the spectral transformation context (for example,
489
   to free it at the end of the computations).
490
 
491
   Level: advanced
492
 
493
.seealso: EPSGetST()
494
@*/
495
PetscErrorCode EPSSetST(EPS eps,ST st)
496
{
497
  PetscErrorCode ierr;
498
 
499
  PetscFunctionBegin;
500
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
501
  PetscValidHeaderSpecific(st,ST_COOKIE,2);
502
  PetscCheckSameComm(eps,1,st,2);
503
  ierr = STDestroy(eps->OP); CHKERRQ(ierr);
504
  eps->OP = st;
863 dsic.upv.es!jroman 505
  ierr = PetscObjectReference((PetscObject)eps->OP);CHKERRQ(ierr);
519 dsic.upv.es!antodo 506
  PetscFunctionReturn(0);
507
}
508
 
509
#undef __FUNCT__  
510
#define __FUNCT__ "EPSGetST"
707 dsic.upv.es!antodo 511
/*@C
519 dsic.upv.es!antodo 512
   EPSGetST - Obtain the spectral transformation (ST) object associated
513
   to the eigensolver object.
514
 
515
   Not Collective
516
 
517
   Input Parameters:
518
.  eps - eigensolver context obtained from EPSCreate()
519
 
520
   Output Parameter:
521
.  st - spectral transformation context
522
 
523
   Level: beginner
524
 
525
.seealso: EPSSetST()
526
@*/
527
PetscErrorCode EPSGetST(EPS eps, ST *st)
528
{
529
  PetscFunctionBegin;
530
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
531
  *st = eps->OP;
532
  PetscFunctionReturn(0);
533
}
534
 
535
#undef __FUNCT__  
536
#define __FUNCT__ "EPSIsGeneralized"
537
/*@
538
   EPSIsGeneralized - Ask if the EPS object corresponds to a generalized
539
   eigenvalue problem.
540
 
541
   Not collective
542
 
543
   Input Parameter:
544
.  eps - the eigenproblem solver context
545
 
546
   Output Parameter:
547
.  is - the answer
548
 
549
   Level: intermediate
550
 
551
@*/
552
PetscErrorCode EPSIsGeneralized(EPS eps,PetscTruth* is)
553
{
554
  PetscErrorCode ierr;
555
  Mat            B;
556
 
557
  PetscFunctionBegin;
558
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
559
  ierr = STGetOperators(eps->OP,PETSC_NULL,&B);CHKERRQ(ierr);
560
  if( B ) *is = PETSC_TRUE;
561
  else *is = PETSC_FALSE;
562
  if( eps->setupcalled ) {
563
    if( eps->isgeneralized != *is ) {
564
      SETERRQ(0,"Warning: Inconsistent EPS state");
565
    }
566
  }
567
  PetscFunctionReturn(0);
568
}
569
 
570
#undef __FUNCT__  
571
#define __FUNCT__ "EPSIsHermitian"
572
/*@
573
   EPSIsHermitian - Ask if the EPS object corresponds to a Hermitian
574
   eigenvalue problem.
575
 
576
   Not collective
577
 
578
   Input Parameter:
579
.  eps - the eigenproblem solver context
580
 
581
   Output Parameter:
582
.  is - the answer
583
 
584
   Level: intermediate
585
 
586
@*/
587
PetscErrorCode EPSIsHermitian(EPS eps,PetscTruth* is)
588
{
589
  PetscFunctionBegin;
590
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
591
  if( eps->ishermitian ) *is = PETSC_TRUE;
592
  else *is = PETSC_FALSE;
593
  PetscFunctionReturn(0);
594
}