libMesh
rb_eim_assembly.h
Go to the documentation of this file.
1 // rbOOmit: An implementation of the Certified Reduced Basis method.
2 // Copyright (C) 2009, 2010 David J. Knezevic
3 
4 // This file is part of rbOOmit.
5 
6 // rbOOmit is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License, or (at your option) any later version.
10 
11 // rbOOmit is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 // Lesser General Public License for more details.
15 
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
19 
20 #ifndef LIBMESH_RB_EIM_ASSEMBLY_H
21 #define LIBMESH_RB_EIM_ASSEMBLY_H
22 
23 // rbOOmit includes
24 #include "libmesh/elem_assembly.h"
25 
26 // libMesh includes
27 #include "libmesh/auto_ptr.h"
28 #include "libmesh/numeric_vector.h"
29 #include "libmesh/point.h"
30 #include "libmesh/fe.h"
31 
32 namespace libMesh
33 {
34 
35 class Elem;
36 class RBParameters;
37 class RBEIMConstruction;
38 
48 {
49 public:
50 
54  RBEIMAssembly(RBEIMConstruction & rb_eim_con_in,
55  unsigned int basis_function_index_in);
56 
60  virtual ~RBEIMAssembly();
61 
66  virtual void evaluate_basis_function(unsigned int var,
67  const Elem & element,
68  const QBase & element_qrule,
69  std::vector<Number> & values);
70 
75 
80 
84  FEBase & get_fe();
85 
86 private:
87 
91  void initialize_fe();
92 
97 
101  unsigned int _basis_function_index;
102 
109 
115 };
116 
117 }
118 
119 #endif // LIBMESH_RB_EIM_ASSEMBLY_H
This class is part of the rbOOmit framework.
This is the base class from which all geometric element types are derived.
Definition: elem.h:89
virtual ~RBEIMAssembly()
Destructor.
UniquePtr< QBase > _qrule
The libMesh namespace provides an interface to certain functionality in the library.
RBEIMConstruction & get_rb_eim_construction()
Get a reference to the RBEIMConstruction object.
std::unique_ptr< T > UniquePtr
Definition: auto_ptr.h:46
This class provides functionality required to define an assembly object that arises from an "Empirica...
void initialize_fe()
Initialize the FE object.
RBEIMAssembly(RBEIMConstruction &rb_eim_con_in, unsigned int basis_function_index_in)
Constructor.
NumericVector< Number > & get_ghosted_basis_function()
Get a reference to the ghosted_basis_function.
unsigned int _basis_function_index
The EIM basis function index (from rb_eim_eval) for this assembly object.
virtual void evaluate_basis_function(unsigned int var, const Elem &element, const QBase &element_qrule, std::vector< Number > &values)
Evaluate variable var_number of this object&#39;s EIM basis function at the points qpoints.
FEBase & get_fe()
Retrieve the FE object.
ElemAssembly provides a per-element (interior and boundary) assembly functionality.
Definition: elem_assembly.h:38
RBEIMConstruction & _rb_eim_con
The RBEIMConstruction object that this RBEIMAssembly is based on.
UniquePtr< FEBase > _fe
We store an FE object and an associated quadrature rule.
The QBase class provides the basic functionality from which various quadrature rules can be derived...
Definition: quadrature.h:53
This class forms the foundation from which generic finite elements may be derived.
UniquePtr< NumericVector< Number > > _ghosted_basis_function
The basis function that we sample to evaluate the empirical interpolation approximation.