Skip to content
Snippets Groups Projects
Commit bbf130f7 authored by oliver.sander_at_tu-dresden.de's avatar oliver.sander_at_tu-dresden.de
Browse files

New abstract base class LocalEnergy

This is supposed to the base class for everything that can
integrate energy densities over one grid element.

LocalFEStiffness inherits from the new LocalEnergy.
parent 51754ca8
Branches
No related tags found
No related merge requests found
install(FILES install(FILES
feassembler.hh feassembler.hh
localadolcstiffness.hh localadolcstiffness.hh
localenergy.hh
localfestiffness.hh localfestiffness.hh
neohookefunctionalassembler.hh neohookefunctionalassembler.hh
neohookeoperatorassembler.hh neohookeoperatorassembler.hh
......
#ifndef DUNE_ELASTICITY_ASSEMBLERS_LOCALENERGY_HH
#define DUNE_ELASTICITY_ASSEMBLERS_LOCALENERGY_HH
#include <vector>
namespace Dune {
namespace Elasticity {
/** \brief Base class for energies defined by integrating over one grid element */
template<class GridView, class LocalFiniteElement, class VectorType>
class LocalEnergy
{
typedef typename VectorType::value_type::field_type RT;
typedef typename GridView::template Codim<0>::Entity Element;
public:
/** \brief Compute the energy
*
* \param element A grid element
* \param LocalFiniteElement A finite element on the reference element
* \param localConfiguration The coefficients of a FE function on the current element
*/
/** \brief Compute the energy at the current configuration */
virtual RT energy (const Element& element,
const LocalFiniteElement& localFiniteElement,
const VectorType& localConfiguration) const = 0;
};
} // namespace Elasticity
} // namespace Dune
#endif // DUNE_ELASTICITY_ASSEMBLERS_LOCALENERGY_HH
...@@ -4,9 +4,11 @@ ...@@ -4,9 +4,11 @@
#include <dune/common/fmatrix.hh> #include <dune/common/fmatrix.hh>
#include <dune/istl/matrix.hh> #include <dune/istl/matrix.hh>
#include <dune/elasticity/assemblers/localenergy.hh>
template<class GridView, class LocalFiniteElement, class VectorType> template<class GridView, class LocalFiniteElement, class VectorType>
class LocalFEStiffness class LocalFEStiffness
: public Dune::Elasticity::LocalEnergy<GridView, LocalFiniteElement, VectorType>
{ {
// grid types // grid types
typedef typename GridView::Grid::ctype DT; typedef typename GridView::Grid::ctype DT;
...@@ -28,11 +30,6 @@ public: ...@@ -28,11 +30,6 @@ public:
const VectorType& localConfiguration, const VectorType& localConfiguration,
VectorType& localGradient); VectorType& localGradient);
/** \brief Compute the energy at the current configuration */
virtual RT energy (const Entity& e,
const LocalFiniteElement& localFiniteElement,
const VectorType& localConfiguration) const = 0;
// assembled data // assembled data
Dune::Matrix<Dune::FieldMatrix<RT,blocksize,blocksize> > A_; Dune::Matrix<Dune::FieldMatrix<RT,blocksize,blocksize> > A_;
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment