2 #ifndef DUNE_PDELAB_LINEARELASTICITY_HH 3 #define DUNE_PDELAB_LINEARELASTICITY_HH 7 #include <dune/common/exceptions.hh> 8 #include <dune/common/fvector.hh> 10 #include <dune/geometry/type.hh> 11 #include <dune/geometry/referenceelements.hh> 12 #include <dune/geometry/quadraturerules.hh> 60 template<
typename EG,
typename LFSU,
typename X,
typename LFSV,
typename M>
61 void jacobian_volume (
const EG& eg,
const LFSU& lfsu,
const X& x,
const LFSV& lfsv, M & mat)
const 64 using namespace TypeTree::Indices;
65 using LFSU_SUB = TypeTree::Child<LFSU,0>;
66 using RF =
typename M::value_type;
67 using JacobianType =
typename LFSU_SUB::Traits::FiniteElementType::
68 Traits::LocalBasisType::Traits::JacobianType;
69 using size_type =
typename LFSU_SUB::Traits::SizeType;
72 const int dim = EG::Entity::dimension;
73 const int dimw = EG::Geometry::coorddimension;
74 static_assert(dim == dimw,
"doesn't work on manifolds");
77 const auto& cell = eg.entity();
80 auto geo = eg.geometry();
83 typename EG::Geometry::JacobianInverseTransposed jac;
86 std::vector<JacobianType> js(lfsu.child(0).size());
87 std::vector<FieldVector<RF,dim> > gradphi(lfsu.child(0).size());
93 lfsu.child(0).finiteElement().localBasis().evaluateJacobian(qp.position(),js);
96 jac = geo.jacobianInverseTransposed(qp.position());
97 for (size_type i=0; i<lfsu.child(0).size(); i++)
100 jac.umv(js[i][0],gradphi[i]);
104 auto mu =
param_.mu(cell,qp.position());
105 auto lambda =
param_.lambda(cell,qp.position());
108 auto factor = qp.weight() * geo.integrationElement(qp.position());
110 for(
int d=0; d<
dim; ++d)
112 for (size_type i=0; i<lfsu.child(0).size(); i++)
114 for (
int k=0; k<
dim; k++)
116 for (size_type j=0; j<lfsv.child(k).size(); j++)
119 mat.accumulate(lfsv.child(k),j,lfsu.child(k),i,
121 mu * gradphi[i][d] * gradphi[j][d]
123 mat.accumulate(lfsv.child(k),j,lfsu.child(d),i,
125 mu * gradphi[i][k] * gradphi[j][d]
128 mat.accumulate(lfsv.child(k),j,lfsu.child(d),i,
129 lambda * gradphi[i][d]*gradphi[j][k]
139 template<
typename EG,
typename LFSU_HAT,
typename X,
typename LFSV,
typename R>
140 void alpha_volume (
const EG& eg,
const LFSU_HAT& lfsu_hat,
const X& x,
const LFSV& lfsv, R& r)
const 143 using namespace TypeTree::Indices;
144 using LFSU = TypeTree::Child<LFSU_HAT,0>;
145 using RF =
typename R::value_type;
146 using JacobianType =
typename LFSU::Traits::FiniteElementType::
147 Traits::LocalBasisType::Traits::JacobianType;
148 using size_type =
typename LFSU::Traits::SizeType;
151 const int dim = EG::Entity::dimension;
152 const int dimw = EG::Geometry::coorddimension;
153 static_assert(dim == dimw,
"doesn't work on manifolds");
156 const auto& cell = eg.entity();
159 auto geo = eg.geometry();
162 typename EG::Geometry::JacobianInverseTransposed jac;
165 std::vector<JacobianType> js(lfsu_hat.child(0).size());
166 std::vector<FieldVector<RF,dim> > gradphi(lfsu_hat.child(0).size());
167 Dune::FieldVector<RF,dim> gradu(0.0);
173 lfsu_hat.child(0).finiteElement().localBasis().evaluateJacobian(qp.position(),js);
176 jac = geo.jacobianInverseTransposed(qp.position());
177 for (size_type i=0; i<lfsu_hat.child(0).size(); i++)
180 jac.umv(js[i][0],gradphi[i]);
184 auto mu =
param_.mu(cell,qp.position());
185 auto lambda =
param_.lambda(cell,qp.position());
188 auto factor = qp.weight() * geo.integrationElement(qp.position());
190 for(
int d=0; d<
dim; ++d)
192 const LFSU & lfsu = lfsu_hat.child(d);
196 for (
size_t i=0; i<lfsu.size(); i++)
198 gradu.axpy(x(lfsu,i),gradphi[i]);
201 for (size_type i=0; i<lfsv.child(d).size(); i++)
203 for (
int k=0; k<
dim; k++)
206 r.accumulate(lfsv.child(d),i,
208 mu * gradu[k] * gradphi[i][k]
210 r.accumulate(lfsv.child(k),i,
212 mu * gradu[k] * gradphi[i][d]
215 r.accumulate(lfsv.child(k),i,
216 lambda * gradu[d]*gradphi[i][k]
225 template<
typename EG,
typename LFSV_HAT,
typename R>
229 using namespace TypeTree::Indices;
230 using LFSV = TypeTree::Child<LFSV_HAT,0>;
231 using RF =
typename R::value_type;
232 using RangeType =
typename LFSV::Traits::FiniteElementType::
233 Traits::LocalBasisType::Traits::RangeType;
234 using size_type =
typename LFSV::Traits::SizeType;
237 const int dim = EG::Entity::dimension;
240 const auto& cell = eg.entity();
243 auto geo = eg.geometry();
246 std::vector<RangeType> phi(lfsv_hat.child(0).size());
247 FieldVector<RF,dim> y(0.0);
253 lfsv_hat.child(0).finiteElement().localBasis().evaluateFunction(qp.position(),phi);
257 param_.f(cell,qp.position(),y);
260 auto factor = qp.weight() * geo.integrationElement(qp.position());
262 for(
int d=0; d<
dim; ++d)
264 const auto& lfsv = lfsv_hat.child(d);
267 for (size_type i=0; i<lfsv.size(); i++)
268 r.accumulate(lfsv,i, -y[d]*phi[i] * factor);
274 template<
typename IG,
typename LFSV_HAT,
typename R>
278 using namespace TypeTree::Indices;
279 using LFSV = TypeTree::Child<LFSV_HAT,0>;
280 using RF =
typename R::value_type;
281 using RangeType =
typename LFSV::Traits::FiniteElementType::
282 Traits::LocalBasisType::Traits::RangeType;
283 using size_type =
typename LFSV::Traits::SizeType;
286 const int dim = IG::Entity::dimension;
289 auto geo = ig.geometry();
292 auto geo_in_inside = ig.geometryInInside();
295 std::vector<RangeType> phi(lfsv_hat.child(0).size());
296 FieldVector<RF,dim> y(0.0);
302 auto local = geo_in_inside.global(qp.position());
306 if(
param_.isDirichlet( ig.intersection(), qp.position() ) )
310 lfsv_hat.child(0).finiteElement().localBasis().evaluateFunction(local,phi);
318 auto factor = qp.weight() * geo.integrationElement(qp.position());
320 for(
int d=0; d<
dim; ++d)
322 const auto& lfsv = lfsv_hat.child(d);
325 for (size_type i=0; i<lfsv.size(); i++)
326 r.accumulate(lfsv,i, y[d]*phi[i] * factor);
void jacobian_volume(const EG &eg, const LFSU &lfsu, const X &x, const LFSV &lfsv, M &mat) const
Definition: linearelasticity.hh:61
int intorder_
Definition: linearelasticity.hh:332
sparsity pattern generator
Definition: pattern.hh:13
Definition: linearelasticity.hh:40
const IG & ig
Definition: constraints.hh:148
static const int dim
Definition: adaptivity.hh:83
Definition: adaptivity.hh:27
Default flags for all local operators.
Definition: flags.hh:18
QuadratureRuleWrapper< QuadratureRule< typename Geometry::ctype, Geometry::mydimension > > quadratureRule(const Geometry &geo, std::size_t order, QuadratureType::Enum quadrature_type=QuadratureType::GaussLegendre)
Returns a quadrature rule for the given geometry.
Definition: quadraturerules.hh:117
Default class for additional methods in instationary local operators.
Definition: idefault.hh:89
T ParameterType
Definition: linearelasticity.hh:46
void lambda_volume(const EG &eg, const LFSV_HAT &lfsv_hat, R &r) const
Definition: linearelasticity.hh:226
const P & p
Definition: constraints.hh:147
LinearElasticity(const ParameterType &p, int intorder=4)
Definition: linearelasticity.hh:56
void lambda_boundary(const IG &ig, const LFSV_HAT &lfsv_hat, R &r) const
Definition: linearelasticity.hh:275
Definition: linearelasticity.hh:49
Definition: linearelasticity.hh:54
Definition: linearelasticity.hh:53
Definition: linearelasticity.hh:52
void alpha_volume(const EG &eg, const LFSU_HAT &lfsu_hat, const X &x, const LFSV &lfsv, R &r) const
Definition: linearelasticity.hh:140
const ParameterType & param_
Definition: linearelasticity.hh:333