25 #include <queso/Defines.h>
27 #ifdef QUESO_HAVE_LIBMESH_SLEPC
32 #include <queso/FunctionOperatorBuilder.h>
33 #include <queso/LibMeshFunction.h>
34 #include <queso/LibMeshOperatorBase.h>
36 #include <libmesh/libmesh.h>
37 #include <libmesh/mesh.h>
38 #include <libmesh/mesh_generation.h>
39 #include <libmesh/equation_systems.h>
40 #include <libmesh/explicit_system.h>
41 #include <libmesh/eigen_system.h>
42 #include <libmesh/condensed_eigen_system.h>
43 #include <libmesh/numeric_vector.h>
44 #include <libmesh/exodusII_io.h>
45 #include <libmesh/parallel.h>
46 #include <libmesh/parallel_implementation.h>
53 equation_systems(new libMesh::EquationSystems(m)),
56 #ifndef LIBMESH_HAVE_SLEPC
57 if (m.processor_id() == 0)
58 std::cerr <<
"ERROR: This example requires libMesh to be\n"
59 <<
"compiled with SLEPc eigen solvers support!"
63 #ifdef LIBMESH_DEFAULT_SINGLE_PRECISION
65 libmesh_example_assert(
false,
"--disable-singleprecision");
68 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
71 libmesh_example_assert(
false,
"--disable-complex");
76 this->
equation_systems->add_system<libMesh::CondensedEigenSystem>(
"Eigensystem");
78 #endif // LIBMESH_HAVE_SLEPC
88 std::ofstream evals_file(filename.c_str());
90 std::pair<libMesh::Real, libMesh::Real> eval;
91 for (i = 0; i < this->
nconv; i++) {
93 ->get_system<libMesh::EigenSystem>(
"Eigensystem").get_eigenpair(i);
95 evals_file << eval.first <<
" " << eval.second << std::endl;
104 unsigned int i)
const
106 if (i < this->
nconv) {
108 es->get_system<libMesh::EigenSystem>(
"Eigensystem").get_eigenpair(i);
109 libMesh::ExodusII_IO(es->get_mesh()).write_equation_systems(filename, *es);
112 std::cerr <<
"Warning: eigenpair " << i
113 <<
" did not converge. Not saving."
124 if (i < this->
nconv) {
125 std::pair<libMesh::Real, libMesh::Real> eval;
127 eval = es->get_system<libMesh::EigenSystem>(
"Eigensystem").get_eigenpair(i);
160 std::pair<libMesh::Real, libMesh::Real> eval;
162 eval = es->get_system<libMesh::EigenSystem>(
"Eigensystem").get_eigenpair(i);
163 kl_eq_sys->get_system<libMesh::ExplicitSystem>(
"Function").solution->add(
164 xi[i] / pow(eval.first, alpha / 2.0),
165 *es->get_system<libMesh::EigenSystem>(
"Eigensystem").solution);
174 #endif // QUESO_HAVE_LIBMESH_SLEPC
std::shared_ptr< T > Type
virtual void save_converged_evec(const std::string &filename, unsigned int i) const
Save converged eigenfunction i to filename.
virtual void save_converged_evals(const std::string &filename) const
Save the eigenvalues to file filename.
Abstract base class for operator objects. Operators are assumed to be symmetric and positive-definite...
const FunctionOperatorBuilder & builder
~LibMeshOperatorBase()
Destructor.
virtual SharedPtr< libMesh::EquationSystems >::Type get_equation_systems() const
Return the internal libmesh equation systems object.
LibMeshOperatorBase(const FunctionOperatorBuilder &builder, libMesh::MeshBase &m)
Constuct an operator on the mesh m using a builder builder.
SharedPtr< libMesh::EquationSystems >::Type equation_systems
virtual double get_inverted_eigenvalue(unsigned int i) const
Return the reciprocal of eigenvalue i.
virtual unsigned int get_num_converged() const
Return the number of converged eigenpairs.
virtual double get_eigenvalue(unsigned int i) const
Return eigenvalue i.
unsigned int nconv
The number of converged eigenpairs.
virtual SharedPtr< FunctionBase >::Type inverse_kl_transform(std::vector< double > &xi, double alpha) const
Given coefficients xi, computes the Karhunen-Loeve transform.
Function objects using libMesh for the backend.
virtual libMesh::EquationSystems & get_equation_systems() const
Return the internal libmesh equation systems object.