queso-0.57.0
LibMeshOperatorBase.h
Go to the documentation of this file.
1 //-----------------------------------------------------------------------bl-
2 //--------------------------------------------------------------------------
3 //
4 // QUESO - a library to support the Quantification of Uncertainty
5 // for Estimation, Simulation and Optimization
6 //
7 // Copyright (C) 2008-2017 The PECOS Development Team
8 //
9 // This library is free software; you can redistribute it and/or
10 // modify it under the terms of the Version 2.1 GNU Lesser General
11 // Public License as published by the Free Software Foundation.
12 //
13 // This library is distributed in the hope that it will be useful,
14 // but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // Lesser General Public License for more details.
17 //
18 // You should have received a copy of the GNU Lesser General Public
19 // License along with this library; if not, write to the Free Software
20 // Foundation, Inc. 51 Franklin Street, Fifth Floor,
21 // Boston, MA 02110-1301 USA
22 //
23 //-----------------------------------------------------------------------el-
24 
25 #include <queso/Defines.h>
26 
27 #ifdef QUESO_HAVE_LIBMESH_SLEPC
28 
29 #ifndef QUESO_LIBMESHOPERATOR_BASE_H
30 #define QUESO_LIBMESHOPERATOR_BASE_H
31 
32 #include <string>
33 #include <set>
34 #include <vector>
35 #include <queso/SharedPtr.h>
36 #include <queso/FunctionBase.h>
37 #include <queso/OperatorBase.h>
38 #include <libmesh/system.h>
39 
40 namespace libMesh {
41  class Mesh;
42  class EquationSystems;
43 }
44 
45 namespace QUESO {
46 
47 class FunctionOperatorBuilder;
48 
58  public libMesh::System::Assembly {
59 public:
61 
62 
69  libMesh::MeshBase & m);
70 
74 
76 
79  virtual void assemble() = 0;
80 
82  virtual void print_info() const = 0;
83 
85  virtual void save_converged_evals(const std::string &filename) const;
86 
88  virtual void
89  save_converged_evec(const std::string &filename, unsigned int i) const;
90 
92  virtual unsigned int get_num_converged() const;
93 
95 
99  virtual double get_eigenvalue(unsigned int i) const;
100 
102  virtual double get_inverted_eigenvalue(unsigned int i) const;
103 
105  virtual libMesh::EquationSystems & get_equation_systems() const;
106 
108 
116  inverse_kl_transform(std::vector<double> & xi, double alpha) const;
117 
118 protected:
120 
122 
124  unsigned int num_req_pairs;
125 
127  unsigned int nconv;
128 };
129 
130 } // End namespace QUESO
131 
132 #endif // QUESO_LIBMESHOPERATOR_BASE_H
133 
134 #endif // QUESO_HAVE_LIBMESH_SLEPC
virtual void assemble()=0
Must implement this for the solve to work.
virtual unsigned int get_num_converged() const
Return the number of converged eigenpairs.
SharedPtr< libMesh::EquationSystems >::Type equation_systems
virtual SharedPtr< FunctionBase >::Type inverse_kl_transform(std::vector< double > &xi, double alpha) const
Given coefficients xi, computes the Karhunen-Loeve transform.
virtual libMesh::EquationSystems & get_equation_systems() const
Return the internal libmesh equation systems object.
unsigned int num_req_pairs
The number of requested eigenpairs.
Abstract base class for operator objects. Operators are assumed to be symmetric and positive-definite...
Definition: OperatorBase.h:45
LibMeshOperatorBase(const FunctionOperatorBuilder &builder, libMesh::MeshBase &m)
Constuct an operator on the mesh m using a builder builder.
const FunctionOperatorBuilder & builder
virtual void save_converged_evals(const std::string &filename) const
Save the eigenvalues to file filename.
unsigned int nconv
The number of converged eigenpairs.
virtual void save_converged_evec(const std::string &filename, unsigned int i) const
Save converged eigenfunction i to filename.
std::shared_ptr< T > Type
Definition: SharedPtr.h:69
virtual double get_inverted_eigenvalue(unsigned int i) const
Return the reciprocal of eigenvalue i.
virtual double get_eigenvalue(unsigned int i) const
Return eigenvalue i.
Abstract base class for operator objects using libmesh in the backend.
virtual void print_info() const =0
Print libmesh related information.

Generated on Sat Apr 22 2017 14:04:35 for queso-0.57.0 by  doxygen 1.8.5