4 #ifndef DUNE_PDELAB_LOCALOPERATOR_LAPLACE_HH
5 #define DUNE_PDELAB_LOCALOPERATOR_LAPLACE_HH
9 #include <dune/common/fvector.hh>
11 #include <dune/geometry/type.hh>
12 #include <dune/geometry/quadraturerules.hh>
14 #include <dune/localfunctions/common/interfaceswitch.hh>
65 template<
typename EG,
typename LFSU,
typename X,
typename LFSV,
typename R>
66 void alpha_volume (
const EG&
eg,
const LFSU& lfsu,
const X& x,
const LFSV& lfsv, R& r)
const
69 typedef FiniteElementInterfaceSwitch<
70 typename LFSU::Traits::FiniteElementType
72 typedef BasisInterfaceSwitch<
73 typename FESwitch::Basis
75 typedef typename BasisSwitch::DomainField DF;
76 typedef typename BasisSwitch::RangeField RF;
79 static const int dimLocal = EG::Geometry::mydimension;
80 static const int dimGlobal = EG::Geometry::coorddimension;
83 Dune::GeometryType gt = eg.geometry().type();
84 const Dune::QuadratureRule<DF,dimLocal>& rule =
85 Dune::QuadratureRules<DF,dimLocal>::rule(gt,
quadOrder_);
88 for(
typename Dune::QuadratureRule<DF,dimLocal>::const_iterator it =
89 rule.begin(); it!=rule.end(); ++it)
93 std::vector<Dune::FieldMatrix<RF,1,dimGlobal> >
94 gradphiu(lfsu.size());
95 BasisSwitch::gradient(FESwitch::basis(lfsu.finiteElement()),
96 eg.geometry(), it->position(), gradphiu);
97 std::vector<Dune::FieldMatrix<RF,1,dimGlobal> >
98 gradphiv(lfsv.size());
99 BasisSwitch::gradient(FESwitch::basis(lfsv.finiteElement()),
100 eg.geometry(), it->position(), gradphiv);
103 Dune::FieldVector<RF,dimGlobal> gradu(0.0);
104 for (
size_t i=0; i<lfsu.size(); i++)
105 gradu.axpy(x(lfsu,i),gradphiu[i][0]);
108 RF factor = r.weight() * it->weight() * eg.geometry().integrationElement(it->position());
109 for (
size_t i=0; i<lfsv.size(); i++)
110 r.rawAccumulate(lfsv,i,(gradu*gradphiv[i][0])*factor);
124 template<
typename EG,
typename LFSU,
typename X,
typename LFSV,
typename M>
125 void jacobian_volume (
const EG&
eg,
const LFSU& lfsu,
const X& x,
const LFSV& lfsv, M & matrix)
const
128 typedef FiniteElementInterfaceSwitch<
129 typename LFSU::Traits::FiniteElementType
131 typedef BasisInterfaceSwitch<
132 typename FESwitch::Basis
136 typedef typename BasisSwitch::DomainField DF;
137 typedef typename BasisSwitch::RangeField RF;
138 typedef typename LFSU::Traits::SizeType size_type;
141 const int dim = EG::Entity::dimension;
144 Dune::GeometryType gt = eg.geometry().type();
145 const Dune::QuadratureRule<DF,dim>& rule = Dune::QuadratureRules<DF,dim>::rule(gt,
quadOrder_);
148 for (
typename QuadratureRule<DF,dim>::const_iterator it=rule.begin(); it!=rule.end(); ++it)
150 std::vector<Dune::FieldMatrix<RF,1,dim> > gradphi(lfsu.size());
151 BasisSwitch::gradient(FESwitch::basis(lfsu.finiteElement()),
152 eg.geometry(), it->position(), gradphi);
155 RF factor = it->weight() * eg.geometry().integrationElement(it->position());
157 for (size_type i=0; i<lfsu.size(); i++)
159 for (size_type j=0; j<lfsv.size(); j++)
162 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:125
Laplace(unsigned int quadOrder)
Constructor.
Definition: laplace.hh:51
Definition: laplace.hh:42
Default flags for all local operators.
Definition: flags.hh:18
sparsity pattern generator
Definition: pattern.hh:13
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:66
Definition: adaptivity.hh:27
Definition: laplace.hh:36
static const int dim
Definition: adaptivity.hh:83
unsigned int quadOrder_
Definition: laplace.hh:170
Definition: laplace.hh:45
const EG & eg
Definition: constraints.hh:280