Skip to content
Snippets Groups Projects

Cosserat-Continuum-Nonplanar

Compare and Show latest version
2 files
+ 27
12
Compare changes
  • Side-by-side
  • Inline
Files
2
@@ -23,7 +23,14 @@
#include <dune/localfunctions/lagrange/lfecache.hh>
#endif
template<class Basis, int dim, class field_type=double, class GridFunction =
/** \brief Assembles the cosserat energy for a single element.
*
* \tparam Basis Type of the Basis used for assembling
* \tparam dim Dimension of the Targetspace, 3
* \tparam field_type The coordinate type of the TargetSpace
* \tparam StressFreeStateGridFunction Type of the GridFunction representing the Cosserat shell in a stress free state
*/
template<class Basis, int dim, class field_type=double, class StressFreeStateGridFunction =
Dune::Functions::DiscreteGlobalBasisFunction<Basis,std::vector<Dune::FieldVector<double, Basis::GridView::dimensionworld>> > >
class NonplanarCosseratShellEnergy
: public Dune::GFE::LocalEnergy<Basis,RigidBodyMotion<field_type,dim> >
@@ -43,14 +50,15 @@ class NonplanarCosseratShellEnergy
public:
/** \brief Constructor with a set of material parameters
* \param parameters The material parameters
* \param parameters The material parameters
* \param stressFreeStateGridFunction Pointer to a parametrization representing the Cosserat shell in a stress-free state
*/
NonplanarCosseratShellEnergy(const Dune::ParameterTree& parameters,
const GridFunction* gridFunctionPtr,
const StressFreeStateGridFunction* stressFreeStateGridFunction,
const BoundaryPatch<GridView>* neumannBoundary,
const std::function<Dune::FieldVector<double,3>(Dune::FieldVector<double,dimworld>)> neumannFunction,
const std::function<Dune::FieldVector<double,3>(Dune::FieldVector<double,dimworld>)> volumeLoad)
: gridFunctionPtr_(gridFunctionPtr),
: stressFreeStateGridFunction_(stressFreeStateGridFunction),
neumannBoundary_(neumannBoundary),
neumannFunction_(neumannFunction),
volumeLoad_(volumeLoad)
@@ -117,7 +125,7 @@ public:
double b1_, b2_, b3_;
/** \brief The geometry used for assembling */
const GridFunction* gridFunctionPtr_;
const StressFreeStateGridFunction* stressFreeStateGridFunction_;
/** \brief The Neumann boundary */
const BoundaryPatch<GridView>* neumannBoundary_;
@@ -129,9 +137,9 @@ public:
const std::function<Dune::FieldVector<double,3>(Dune::FieldVector<double,dimworld>)> volumeLoad_;
};
template <class Basis, int dim, class field_type, class GridFunction>
typename NonplanarCosseratShellEnergy<Basis, dim, field_type, GridFunction>::RT
NonplanarCosseratShellEnergy<Basis,dim,field_type, GridFunction>::
template <class Basis, int dim, class field_type, class StressFreeStateGridFunction>
typename NonplanarCosseratShellEnergy<Basis, dim, field_type, StressFreeStateGridFunction>::RT
NonplanarCosseratShellEnergy<Basis,dim,field_type, StressFreeStateGridFunction>::
energy(const typename Basis::LocalView& localView,
const std::vector<RigidBodyMotion<field_type,dim> >& localSolution) const
{
@@ -139,15 +147,18 @@ energy(const typename Basis::LocalView& localView,
auto element = localView.element();
#if HAVE_DUNE_CURVEDGEOMETRY
// Construct a curved geometry to use a better approximation of the geometry of the element
// Construct a curved geometry of this element of the Cosserat shell in stress-free state
// When using element.geometry(), then the curvatures on the element are zero, when using a curved geometry, they are not
// If a parametrization representing the Cosserat shell in a stress-free state is given,
// this is used for the curved geometry approximation.
// The variable local holds the local coordinates in the reference element
// and localGeometry.global maps them to the world coordinates
Dune::CurvedGeometry<DT, gridDim, dimworld, Dune::CurvedGeometryTraits<DT, Dune::LagrangeLFECache<DT,DT,gridDim>>> geometry(referenceElement(element),
[this,element](const auto& local) {
if (not gridFunctionPtr_) {
if (not stressFreeStateGridFunction_) {
return element.geometry().global(local);
}
auto localGridFunction = localFunction(*gridFunctionPtr_);
auto localGridFunction = localFunction(*stressFreeStateGridFunction_);
localGridFunction.bind(element);
return localGridFunction(local);
}, 2); /*order*/
Loading