Hello David,
Yes, that is very helpful.  Thanks for getting back to me.  That is very helpful!
I am tried to copy what you've done, see below.   But I do have a few questions.
1) Do you have a topetsc function defined like what I have at the top?
2) When you define the petsc matrix, petsc_a, how did you define that?
3) I guess you are solving a regular eigenvalue problem, or is it a generalized eigenvalue problem? If it's in Finite elements i presume it's generalized?
In the case of a generalized eigenvalue problem I wonder if I should have
es.setOperators(A,M)
where A and M are petsc matrices.
4) If yes, then when you have
 
vr, vi = petsc_A.getVec
I guess this is solving the eigenvalue problem but should you give it an A and M?
Sorry to bother you but I feel like I'm close, just missing an important detail.
Cheers, Francis
def topetsc(A):
 
    A.force_evaluation()
    return A.petscmat
# Set up Eigenvalue Problem: QG Basin modes
a =  beta*phi*psi.dx(0)*dx
 
m = -inner(grad(psi), grad(phi))*dx - F*psi*phi*dx
# Assemble Weak Form into a PETSc Matrix
# Impose Boundary Conditions on Mass Matrix
A = topetsc(assemble(a))
M = topetsc(assemble(m, bcs=[bc]))
# Create solver
es = SLEPc.EPS().create(comm=COMM_WORLD)
# Set type of eigenvalues
es.setWhichEigenpairs(es.Which.SMALLEST_IMAGINARY)
# Set number of eigenvalues to compute
num_eigenvalues = 5
es.setDimensions(num_eigenvalues)
#setting the matrix (petsc_a is a PETSc matrix)
es.setOperators(petsc_a)
# Fin eigenvalues an eigenvectors
vr, vi = petsc_a.getVecs()
for i in range(num_eigenvalues):
    lam = es.getEigenpair(i, vr, vi)
    eigenvalue[i] = lam.real
    eigenfunction[i].vector()[:] = vr
Hi Francis, I’m working on a project involving a SLEPc eigensolver in Firedrake so I’ve done some of what you mentioned.
For instance: creating an eigensolver:
es = SLEPc.EPS().create(comm=COMM_WORLD)
setting type of eigenvalues:
es.setWhichEigenpairs(self.eigensolver.Which.SMALLEST_REAL)
setting number of eigenvalues to compute:
es.setDimensions(self.num_eigenvalues)
setting the matrix (petsc_a is a PETSc matrix)
es.setOperators(petsc_a)
getting eigenvectors (in my case the eigenvalues are all real, petsc_a is the matrix)
vr, vi = petsc_a.getVecs()
for i in range(num_eigenvalues):
lmbda = es.getEigenpair(i, vr, vi)
eigenvalue[i] = lmbda.real
eigenfunction[i].vector()[:] = vr
 
Here eigenfunction is an array of Firedrake functions in the right space that has already been allocated.
Hope that helps,
Dave
Dear
 Francis,
On
 30/11/16 03:23, Francis Poulin wrote:
Hello,
A while ago someone pointed me towards test_slepc.py for an example
on how to solve an eigenvalue problem.  I have been using this as a
model to compute Quasi-Geostrophic basin modes in a closed basin.
We have a version working in FENicS and I am very keen to get it
working in Firedrake as well.
I have a good start, which can be found on github,
https://github.com/francispoulin/firedrakeQG/blob/master/basin_modes_qg.py
Unfortunately, I don't know how to  do a lot of things.  For
example specify number of eigenvalues to determine, specify a
target eigenvalue, what method to use, or even just look at the
eigenvector.
Doing
 this is strictly in the realm of SLEPc, rather than Firedrake
itself.
  We're "just" providing the operators.  The SLEPc
documentation
 is therefore probably helpful:
http://slepc.upv.es/documentation/manual.htm
It
 describes everything in terms of the C API, but fortunately,
slepc4py
 has  somewhat useful documentation strings (type
help(some_slepc_function)
 in ipython).  The translation from C
functions
 to slepc4py functions is quite mechanical:
EPSFooBar(eps,
 ...)
translates
 to:
eps.fooBar(...)
and
 so forth.
Thanks,
Lawrence
©Melior
 Innovations Inc, 2016 all rights reserved. This e-mail is covered by the Electronic Communication Privacy Act, 18 U.S.C. Section 2510-2521 and may be legally privileged. If you have received this transmission in error, please notify the original sender immediately
 by return e-mail or by telephone at the telephone number provided above and delete/trash the original message from your system. Thank you for your assistance. The information contained herein is confidential information provided only for the use of the individual
 or entity for whom it was intended. If the reader of this message received it in error, is not the intended recipient, or if an attachment is made in error, the reader is hereby notified that any dissemination, distribution or copying of this communication
 is strictly prohibited.
_______________________________________________
firedrake
 mailing list
firedrake@imperial.ac.uk
https://mailman.ic.ac.uk/mailman/listinfo/firedrake 
 
 
©Melior Innovations Inc, 2016 all rights reserved. This e-mail is covered by the Electronic Communication Privacy Act, 18 U.S.C. Section 2510-2521 and may be legally privileged. If you have received this transmission in error, please notify the original sender
 immediately by return e-mail or by telephone at the telephone number provided above and delete/trash the original message from your system. Thank you for your assistance. The information contained herein is confidential information provided only for the use
 of the individual or entity for whom it was intended. If the reader of this message received it in error, is not the intended recipient, or if an attachment is made in error, the reader is hereby notified that any dissemination, distribution or copying of
 this communication is strictly prohibited.