/*
EPS routines related to problem setup.
*/
#include "src/eps/epsimpl.h" /*I "slepceps.h" I*/
#undef __FUNCT__
#define __FUNCT__ "EPSSetUp"
/*@
EPSSetUp - Sets up all the internal data structures necessary for the
execution of the eigensolver. Then calls STSetUp() for any set-up
operations associated to the ST object.
Collective on EPS
Input Parameter:
. eps - eigenproblem solver context
Level: advanced
Notes:
This function need not be called explicitly in most cases, since EPSSolve()
calls it. It can be useful when one wants to measure the set-up time
separately from the solve time.
This function sets a random initial vector if none has been provided.
.seealso: EPSCreate(), EPSSolve(), EPSDestroy(), STSetUp()
@*/
PetscErrorCode EPSSetUp(EPS eps)
{
PetscErrorCode ierr;
int i;
Vec v0,w0;
Mat A,B;
PetscInt N;
PetscFunctionBegin;
PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
if (eps->setupcalled) PetscFunctionReturn(0);
ierr = PetscLogEventBegin(EPS_SetUp,eps,0,0,0);CHKERRQ(ierr);
/* Set default solver type */
if (!eps->type_name) {
ierr = EPSSetType(eps,EPSKRYLOVSCHUR);CHKERRQ(ierr);
}
/* Set default eta for orthogonalization */
if (eps->orthog_eta == PETSC_DEFAULT)
eps->orthog_eta = 0.7071;
ierr = STGetOperators(eps->OP,&A,&B);CHKERRQ(ierr);
/* Set default problem type */
if (!eps->problem_type) {
if (B==PETSC_NULL) {
ierr = EPSSetProblemType(eps,EPS_NHEP);CHKERRQ(ierr);
}
else {
ierr = EPSSetProblemType(eps,EPS_GNHEP);CHKERRQ(ierr);
}
} else if ((B && !eps->isgeneralized) || (!B && eps->isgeneralized)) {
SETERRQ(0,"Warning: Inconsistent EPS state");
}
/* Create random initial vectors if not set */
/* right */
ierr = EPSGetInitialVector(eps,&v0);CHKERRQ(ierr);
if (!v0) {
ierr = MatGetVecs(A,&v0,PETSC_NULL);CHKERRQ(ierr);
ierr = SlepcVecSetRandom(v0);CHKERRQ(ierr);
eps->vec_initial = v0;
}
/* left */
ierr = EPSGetLeftInitialVector(eps,&w0);CHKERRQ(ierr);
if (!w0) {
ierr = MatGetVecs(A,PETSC_NULL,&w0);CHKERRQ(ierr);
ierr = SlepcVecSetRandom(w0);CHKERRQ(ierr);
eps->vec_initial_left = w0;
}
ierr = VecGetSize(eps->vec_initial,&N);CHKERRQ(ierr);
if (eps->nev > N) eps->nev = N;
if (eps->ncv > N) eps->ncv = N;
if (!eps->tol) eps->tol = 1.e-7;
ierr = (*eps->ops->setup)(eps);CHKERRQ(ierr);
ierr = STSetUp(eps->OP); CHKERRQ(ierr);
/* DSV is equal to the columns of DS followed by the ones in V */
ierr = PetscFree(eps->DSV);CHKERRQ(ierr);
ierr = PetscMalloc((eps->ncv+eps->nds)*sizeof(Vec),&eps->DSV);CHKERRQ(ierr);
for (i = 0; i < eps->nds; i++) eps->DSV[i] = eps->DS[i];
for (i = 0; i < eps->ncv; i++) eps->DSV[i+eps->nds] = eps->V[i];
if (eps->nds>0) {
if (!eps->ds_ortho) {
/* orthonormalize vectors in DS if necessary */
ierr = EPSQRDecomposition(eps,eps->DS,0,eps->nds,PETSC_NULL,0);CHKERRQ(ierr);
}
ierr = EPSOrthogonalize(eps,eps->nds,PETSC_NULL,eps->DS,eps->vec_initial,PETSC_NULL,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
}
ierr = STCheckNullSpace(eps->OP,eps->nds,eps->DS);CHKERRQ(ierr);
ierr = PetscLogEventEnd(EPS_SetUp,eps,0,0,0);CHKERRQ(ierr);
eps->setupcalled = 1;
PetscFunctionReturn(0);
}
#undef __FUNCT__
#define __FUNCT__ "EPSSetInitialVector"
/*@
EPSSetInitialVector - Sets the initial vector from which the
eigensolver starts to iterate.
Collective on EPS and Vec
Input Parameters:
+ eps - the eigensolver context
- vec - the vector
Level: intermediate
.seealso: EPSGetInitialVector(), EPSSetLeftInitialVector()
@*/
PetscErrorCode EPSSetInitialVector(EPS eps,Vec vec)
{
PetscErrorCode ierr;
PetscFunctionBegin;
PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
PetscValidHeaderSpecific(vec,VEC_COOKIE,2);
PetscCheckSameComm(eps,1,vec,2);
ierr = PetscObjectReference((PetscObject)eps->vec_initial);CHKERRQ(ierr);
if (eps->vec_initial) {
ierr = VecDestroy(eps->vec_initial); CHKERRQ(ierr);
}
eps->vec_initial = vec;
PetscFunctionReturn(0);
}
#undef __FUNCT__
#define __FUNCT__ "EPSGetInitialVector"
/*@
EPSGetInitialVector - Gets the initial vector associated with the
eigensolver; if the vector was not set it will return a 0 pointer or
a vector randomly generated by EPSSetUp().
Not collective, but vector is shared by all processors that share the EPS
Input Parameter:
. eps - the eigensolver context
Output Parameter:
. vec - the vector
Level: intermediate
.seealso: EPSSetInitialVector(), EPSGetLeftInitialVector()
@*/
PetscErrorCode EPSGetInitialVector(EPS eps,Vec *vec)
{
PetscFunctionBegin;
PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
PetscValidPointer(vec,2);
*vec = eps->vec_initial;
PetscFunctionReturn(0);
}
#undef __FUNCT__
#define __FUNCT__ "EPSSetLeftInitialVector"
/*@
EPSSetLeftInitialVector - Sets the initial vector from which the eigensolver
starts to iterate, corresponding to the left recurrence (two-sided solvers).
Collective on EPS and Vec
Input Parameters:
+ eps - the eigensolver context
- vec - the vector
Level: intermediate
.seealso: EPSGetLeftInitialVector(), EPSSetInitialVector()
@*/
PetscErrorCode EPSSetLeftInitialVector(EPS eps,Vec vec)
{
PetscErrorCode ierr;
PetscFunctionBegin;
PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
PetscValidHeaderSpecific(vec,VEC_COOKIE,2);
PetscCheckSameComm(eps,1,vec,2);
ierr = PetscObjectReference((PetscObject)eps->vec_initial_left);CHKERRQ(ierr);
if (eps->vec_initial_left) {
ierr = VecDestroy(eps->vec_initial_left); CHKERRQ(ierr);
}
eps->vec_initial_left = vec;
PetscFunctionReturn(0);
}
#undef __FUNCT__
#define __FUNCT__ "EPSGetLeftInitialVector"
/*@
EPSGetLeftInitialVector - Gets the left initial vector associated with the
eigensolver; if the vector was not set it will return a 0 pointer or
a vector randomly generated by EPSSetUp().
Not collective, but vector is shared by all processors that share the EPS
Input Parameter:
. eps - the eigensolver context
Output Parameter:
. vec - the vector
Level: intermediate
.seealso: EPSSetLeftInitialVector(), EPSGetLeftInitialVector()
@*/
PetscErrorCode EPSGetLeftInitialVector(EPS eps,Vec *vec)
{
PetscFunctionBegin;
PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
PetscValidPointer(vec,2);
*vec = eps->vec_initial_left;
PetscFunctionReturn(0);
}
#undef __FUNCT__
#define __FUNCT__ "EPSSetOperators"
/*@
EPSSetOperators - Sets the matrices associated with the eigenvalue problem.
Collective on EPS and Mat
Input Parameters:
+ eps - the eigenproblem solver context
. A - the matrix associated with the eigensystem
- B - the second matrix in the case of generalized eigenproblems
Notes:
To specify a standard eigenproblem, use PETSC_NULL for parameter B.
Level: beginner
.seealso: EPSSolve(), EPSGetST(), STGetOperators()
@*/
PetscErrorCode EPSSetOperators(EPS eps,Mat A,Mat B)
{
PetscErrorCode ierr;
PetscInt m,n;
PetscFunctionBegin;
PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
PetscValidHeaderSpecific(A,MAT_COOKIE,2);
if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3);
PetscCheckSameComm(eps,1,A,2);
if (B) PetscCheckSameComm(eps,1,B,3);
/* Check for square matrices */
ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr);
if (m!=n) { SETERRQ(1,"A is a non-square matrix"); }
if (B) {
ierr = MatGetSize(B,&m,&n);CHKERRQ(ierr);
if (m!=n) { SETERRQ(1,"B is a non-square matrix"); }
}
ierr = STSetOperators(eps->OP,A,B);CHKERRQ(ierr);
eps->setupcalled = 0; /* so that next solve call will call setup */
/* Destroy randomly generated initial vectors */
if (eps->vec_initial) {
ierr = VecDestroy(eps->vec_initial);CHKERRQ(ierr);
eps->vec_initial = PETSC_NULL;
}
if (eps->vec_initial_left) {
ierr = VecDestroy(eps->vec_initial_left);CHKERRQ(ierr);
eps->vec_initial_left = PETSC_NULL;
}
PetscFunctionReturn(0);
}
#undef __FUNCT__
#define __FUNCT__ "EPSAttachDeflationSpace"
/*@
EPSAttachDeflationSpace - Add vectors to the basis of the deflation space.
Not Collective
Input Parameter:
+ eps - the eigenproblem solver context
. n - number of vectors to add
. ds - set of basis vectors of the deflation space
- ortho - PETSC_TRUE if basis vectors of deflation space are orthonormal
Notes:
When a deflation space is given, the eigensolver seeks the eigensolution
in the restriction of the problem to the orthogonal complement of this
space. This can be used for instance in the case that an invariant
subspace is known beforehand (such as the nullspace of the matrix).
The basis vectors can be provided all at once or incrementally with
several calls to EPSAttachDeflationSpace().
Use a value of PETSC_TRUE for parameter ortho if all the vectors passed
in are known to be mutually orthonormal.
Level: intermediate
.seealso: EPSRemoveDeflationSpace()
@*/
PetscErrorCode EPSAttachDeflationSpace(EPS eps,int n,Vec *ds,PetscTruth ortho)
{
PetscErrorCode ierr;
int i;
Vec *tvec;
PetscFunctionBegin;
PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
tvec = eps->DS;
if (n+eps->nds > 0) {
ierr = PetscMalloc((n+eps->nds)*sizeof(Vec), &eps->DS);CHKERRQ(ierr);
}
if (eps->nds > 0) {
for (i=0; i<eps->nds; i++) eps->DS[i] = tvec[i];
ierr = PetscFree(tvec);CHKERRQ(ierr);
}
for (i=0; i<n; i++) {
ierr = VecDuplicate(ds[i],&eps->DS[i + eps->nds]);CHKERRQ(ierr);
ierr = VecCopy(ds[i],eps->DS[i + eps->nds]);CHKERRQ(ierr);
}
eps->nds += n;
if (!ortho) eps->ds_ortho = PETSC_FALSE;
eps->setupcalled = 0;
PetscFunctionReturn(0);
}
#undef __FUNCT__
#define __FUNCT__ "EPSRemoveDeflationSpace"
/*@
EPSRemoveDeflationSpace - Removes the deflation space.
Not Collective
Input Parameter:
. eps - the eigenproblem solver context
Level: intermediate
.seealso: EPSAttachDeflationSpace()
@*/
PetscErrorCode EPSRemoveDeflationSpace(EPS eps)
{
PetscErrorCode ierr;
PetscFunctionBegin;
PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
if (eps->nds > 0) {
ierr = VecDestroyVecs(eps->DS, eps->nds);CHKERRQ(ierr);
}
eps->ds_ortho = PETSC_TRUE;
eps->setupcalled = 0;
PetscFunctionReturn(0);
}