4 #ifndef DUNE_PDELAB_LOCALOPERATOR_LAPLACE_HH
5 #define DUNE_PDELAB_LOCALOPERATOR_LAPLACE_HH
9 #include <dune/common/fvector.hh>
10 #include <dune/common/static_assert.hh>
12 #include <dune/geometry/type.hh>
13 #include <dune/geometry/quadraturerules.hh>
15 #include <dune/localfunctions/common/interfaceswitch.hh>
66 template<
typename EG,
typename LFSU,
typename X,
typename LFSV,
typename R>
67 void alpha_volume (
const EG&
eg,
const LFSU& lfsu,
const X& x,
const LFSV& lfsv, R& r)
const
70 typedef FiniteElementInterfaceSwitch<
71 typename LFSU::Traits::FiniteElementType
73 typedef BasisInterfaceSwitch<
74 typename FESwitch::Basis
76 typedef typename BasisSwitch::DomainField DF;
77 typedef typename BasisSwitch::RangeField RF;
80 static const int dimLocal = EG::Geometry::mydimension;
81 static const int dimGlobal = EG::Geometry::coorddimension;
84 Dune::GeometryType gt = eg.geometry().type();
85 const Dune::QuadratureRule<DF,dimLocal>& rule =
86 Dune::QuadratureRules<DF,dimLocal>::rule(gt,
quadOrder_);
89 for(
typename Dune::QuadratureRule<DF,dimLocal>::const_iterator it =
90 rule.begin(); it!=rule.end(); ++it)
94 std::vector<Dune::FieldMatrix<RF,1,dimGlobal> >
95 gradphiu(lfsu.size());
96 BasisSwitch::gradient(FESwitch::basis(lfsu.finiteElement()),
97 eg.geometry(), it->position(), gradphiu);
98 std::vector<Dune::FieldMatrix<RF,1,dimGlobal> >
99 gradphiv(lfsv.size());
100 BasisSwitch::gradient(FESwitch::basis(lfsv.finiteElement()),
101 eg.geometry(), it->position(), gradphiv);
104 Dune::FieldVector<RF,dimGlobal> gradu(0.0);
105 for (
size_t i=0; i<lfsu.size(); i++)
106 gradu.axpy(x(lfsu,i),gradphiu[i][0]);
109 RF factor = r.weight() * it->weight() * eg.geometry().integrationElement(it->position());
110 for (
size_t i=0; i<lfsv.size(); i++)
111 r.rawAccumulate(lfsv,i,(gradu*gradphiv[i][0])*factor);
125 template<
typename EG,
typename LFSU,
typename X,
typename LFSV,
typename M>
126 void jacobian_volume (
const EG&
eg,
const LFSU& lfsu,
const X& x,
const LFSV& lfsv, M & matrix)
const
129 typedef FiniteElementInterfaceSwitch<
130 typename LFSU::Traits::FiniteElementType
132 typedef BasisInterfaceSwitch<
133 typename FESwitch::Basis
137 typedef typename BasisSwitch::DomainField DF;
138 typedef typename BasisSwitch::RangeField RF;
139 typedef typename LFSU::Traits::SizeType size_type;
142 const int dim = EG::Geometry::dimension;
145 Dune::GeometryType gt = eg.geometry().type();
146 const Dune::QuadratureRule<DF,dim>& rule = Dune::QuadratureRules<DF,dim>::rule(gt,
quadOrder_);
149 for (
typename QuadratureRule<DF,dim>::const_iterator it=rule.begin(); it!=rule.end(); ++it)
151 std::vector<Dune::FieldMatrix<RF,1,dim> > gradphi(lfsu.size());
152 BasisSwitch::gradient(FESwitch::basis(lfsu.finiteElement()),
153 eg.geometry(), it->position(), gradphi);
156 RF factor = it->weight() * eg.geometry().integrationElement(it->position());
158 for (size_type i=0; i<lfsu.size(); i++)
160 for (size_type j=0; j<lfsv.size(); j++)
163 matrix.accumulate(lfsv,j,lfsu,i, gradphi[i][0] * gradphi[j][0] * factor);
void jacobian_volume(const EG &eg, const LFSU &lfsu, const X &x, const LFSV &lfsv, M &matrix) const
Compute the Laplace stiffness matrix for the element given in 'eg'.
Definition: laplace.hh:126
unsigned int quadOrder_
Definition: laplace.hh:171
Definition: laplace.hh:43
Default flags for all local operators.
Definition: flags.hh:18
static const int dim
Definition: adaptivity.hh:82
void alpha_volume(const EG &eg, const LFSU &lfsu, const X &x, const LFSV &lfsv, R &r) const
Compute Laplace matrix times a given vector for one element.
Definition: laplace.hh:67
Laplace(unsigned int quadOrder)
Constructor.
Definition: laplace.hh:52
sparsity pattern generator
Definition: pattern.hh:14
const EG & eg
Definition: common/constraints.hh:277
Definition: laplace.hh:46
Definition: laplace.hh:37