Skip to content
Snippets Groups Projects
Commit b014e93b authored by graeser's avatar graeser Committed by graeser
Browse files

Geometries are no longer returned by reference.

Thus one should store a copy. Especially the following
is in general forbidden:

  const JIT& j = e.geometry().jacobianInverseTransposed()

[[Imported from SVN: r8215]]
parent af52b941
No related branches found
No related tags found
No related merge requests found
...@@ -76,6 +76,7 @@ void assemblePQ1Stiffness(const GridView& gridView, Matrix& matrix) ...@@ -76,6 +76,7 @@ void assemblePQ1Stiffness(const GridView& gridView, Matrix& matrix)
typedef typename GridView::IndexSet IndexSet; typedef typename GridView::IndexSet IndexSet;
typedef typename GridView::template Codim<0>::Iterator ElementIterator; typedef typename GridView::template Codim<0>::Iterator ElementIterator;
typedef typename GridView::template Codim<0>::Entity Element; typedef typename GridView::template Codim<0>::Entity Element;
typedef typename GridView::template Codim<0>::Entity::Geometry Geometry;
typedef typename Dune::PQkLocalFiniteElementCache<double, double, dim, 1> FiniteElementCache; typedef typename Dune::PQkLocalFiniteElementCache<double, double, dim, 1> FiniteElementCache;
typedef typename FiniteElementCache::FiniteElementType FiniteElement; typedef typename FiniteElementCache::FiniteElementType FiniteElement;
...@@ -93,6 +94,7 @@ void assemblePQ1Stiffness(const GridView& gridView, Matrix& matrix) ...@@ -93,6 +94,7 @@ void assemblePQ1Stiffness(const GridView& gridView, Matrix& matrix)
for (; it != end; ++it) for (; it != end; ++it)
{ {
Element& e = *it; Element& e = *it;
const Geometry geometry = e.geometry();
const FiniteElement& fe = cache.get(e.type()); const FiniteElement& fe = cache.get(e.type());
int localSize = fe.localBasis().size(); int localSize = fe.localBasis().size();
...@@ -115,10 +117,10 @@ void assemblePQ1Stiffness(const GridView& gridView, Matrix& matrix) ...@@ -115,10 +117,10 @@ void assemblePQ1Stiffness(const GridView& gridView, Matrix& matrix)
const LocalCoordinate& quadPos = quad[pt].position(); const LocalCoordinate& quadPos = quad[pt].position();
// get transposed inverse of Jacobian of transformation // get transposed inverse of Jacobian of transformation
const JacobianInverseTransposed& invJacobian = e.geometry().jacobianInverseTransposed(quadPos); const JacobianInverseTransposed invJacobian = geometry.jacobianInverseTransposed(quadPos);
// get integration factor // get integration factor
const double integrationElement = e.geometry().integrationElement(quadPos); const double integrationElement = geometry.integrationElement(quadPos);
// evaluate gradients of shape functions // evaluate gradients of shape functions
fe.localBasis().evaluateJacobian(quadPos, referenceGradients); fe.localBasis().evaluateJacobian(quadPos, referenceGradients);
...@@ -156,6 +158,7 @@ void assemblePQ1Mass(const GridView& gridView, Matrix& matrix) ...@@ -156,6 +158,7 @@ void assemblePQ1Mass(const GridView& gridView, Matrix& matrix)
typedef typename GridView::IndexSet IndexSet; typedef typename GridView::IndexSet IndexSet;
typedef typename GridView::template Codim<0>::Iterator ElementIterator; typedef typename GridView::template Codim<0>::Iterator ElementIterator;
typedef typename GridView::template Codim<0>::Entity Element; typedef typename GridView::template Codim<0>::Entity Element;
typedef typename GridView::template Codim<0>::Entity::Geometry Geometry;
typedef typename Dune::PQkLocalFiniteElementCache<double, double, dim, 1> FiniteElementCache; typedef typename Dune::PQkLocalFiniteElementCache<double, double, dim, 1> FiniteElementCache;
typedef typename FiniteElementCache::FiniteElementType FiniteElement; typedef typename FiniteElementCache::FiniteElementType FiniteElement;
...@@ -171,6 +174,7 @@ void assemblePQ1Mass(const GridView& gridView, Matrix& matrix) ...@@ -171,6 +174,7 @@ void assemblePQ1Mass(const GridView& gridView, Matrix& matrix)
for (; it != end; ++it) for (; it != end; ++it)
{ {
Element& e = *it; Element& e = *it;
const Geometry geometry = e.geometry();
const FiniteElement& fe = cache.get(e.type()); const FiniteElement& fe = cache.get(e.type());
int localSize = fe.localBasis().size(); int localSize = fe.localBasis().size();
...@@ -192,7 +196,7 @@ void assemblePQ1Mass(const GridView& gridView, Matrix& matrix) ...@@ -192,7 +196,7 @@ void assemblePQ1Mass(const GridView& gridView, Matrix& matrix)
const LocalCoordinate& quadPos = quad[pt].position(); const LocalCoordinate& quadPos = quad[pt].position();
// get integration factor // get integration factor
const double integrationElement = e.geometry().integrationElement(quadPos); const double integrationElement = geometry.integrationElement(quadPos);
// evaluate basis functions // evaluate basis functions
fe.localBasis().evaluateFunction(quadPos, values); fe.localBasis().evaluateFunction(quadPos, values);
...@@ -227,6 +231,7 @@ void assemblePQ1RHS(const GridView& gridView, Vector& r, const Function& f) ...@@ -227,6 +231,7 @@ void assemblePQ1RHS(const GridView& gridView, Vector& r, const Function& f)
typedef typename GridView::IndexSet IndexSet; typedef typename GridView::IndexSet IndexSet;
typedef typename GridView::template Codim<0>::Iterator ElementIterator; typedef typename GridView::template Codim<0>::Iterator ElementIterator;
typedef typename GridView::template Codim<0>::Entity Element; typedef typename GridView::template Codim<0>::Entity Element;
typedef typename GridView::template Codim<0>::Entity::Geometry Geometry;
typedef typename Dune::PQkLocalFiniteElementCache<double, double, dim, 1> FiniteElementCache; typedef typename Dune::PQkLocalFiniteElementCache<double, double, dim, 1> FiniteElementCache;
typedef typename FiniteElementCache::FiniteElementType FiniteElement; typedef typename FiniteElementCache::FiniteElementType FiniteElement;
...@@ -243,6 +248,7 @@ void assemblePQ1RHS(const GridView& gridView, Vector& r, const Function& f) ...@@ -243,6 +248,7 @@ void assemblePQ1RHS(const GridView& gridView, Vector& r, const Function& f)
for (; it != end; ++it) for (; it != end; ++it)
{ {
Element& e = *it; Element& e = *it;
const Geometry geometry = e.geometry();
const FiniteElement& fe = cache.get(e.type()); const FiniteElement& fe = cache.get(e.type());
int localSize = fe.localBasis().size(); int localSize = fe.localBasis().size();
...@@ -264,14 +270,14 @@ void assemblePQ1RHS(const GridView& gridView, Vector& r, const Function& f) ...@@ -264,14 +270,14 @@ void assemblePQ1RHS(const GridView& gridView, Vector& r, const Function& f)
const LocalCoordinate& quadPos = quad[pt].position(); const LocalCoordinate& quadPos = quad[pt].position();
// get integration factor // get integration factor
const double integrationElement = e.geometry().integrationElement(quadPos); const double integrationElement = geometry.integrationElement(quadPos);
// evaluate basis functions // evaluate basis functions
fe.localBasis().evaluateFunction(quadPos, values); fe.localBasis().evaluateFunction(quadPos, values);
// evaluate function // evaluate function
FunctionRangeType fAtPos; FunctionRangeType fAtPos;
f.evaluate(e.geometry().global(quadPos), fAtPos); f.evaluate(geometry.global(quadPos), fAtPos);
// add vector entries // add vector entries
double z = quad[pt].weight() * integrationElement; double z = quad[pt].weight() * integrationElement;
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment