27#ifndef EWOMS_QUADRATURE_GEOMETRIES_HH
28#define EWOMS_QUADRATURE_GEOMETRIES_HH
30#include <dune/common/fmatrix.hh>
31#include <dune/common/fvector.hh>
32#include <dune/geometry/type.hh>
42template <
class Scalar,
unsigned dim>
46 enum { numCorners = (1 << dim) };
48 using LocalPosition = Dune::FieldVector<Scalar, dim>;
49 using GlobalPosition = Dune::FieldVector<Scalar, dim>;
51 Dune::GeometryType type()
const
52 {
return Dune::GeometryType((1 << dim) - 1, dim); }
54 template <
class CornerContainer>
55 void setCorners(
const CornerContainer& corners,
unsigned nCorners)
57 for (
unsigned cornerIdx = 0; cornerIdx < nCorners; ++cornerIdx) {
58 for (
unsigned j = 0; j < dim; ++j) {
59 corners_[cornerIdx][j] = corners[cornerIdx][j];
64 for (
unsigned cornerIdx = 0; cornerIdx < nCorners; ++cornerIdx) {
65 center_ += corners_[cornerIdx];
79 GlobalPosition
global(
const LocalPosition& localPos)
const
81 GlobalPosition globalPos(0.0);
83 for (
unsigned cornerIdx = 0; cornerIdx < numCorners; ++cornerIdx) {
95 void jacobian(Dune::FieldMatrix<Scalar, dim, dim>& jac,
96 const LocalPosition& localPos)
const
99 for (
unsigned cornerIdx = 0; cornerIdx < numCorners; ++cornerIdx) {
100 for (
unsigned k = 0; k < dim; ++k) {
101 Scalar dWeight_dk = (cornerIdx & (1 << k)) ? 1 : -1;
102 for (
unsigned j = 0; j < dim; ++j) {
104 if (cornerIdx & (1 << j)) {
105 dWeight_dk *= localPos[j];
108 dWeight_dk *= 1 - localPos[j];
113 jac[k].axpy(dWeight_dk, corners_[cornerIdx]);
125 Dune::FieldMatrix<Scalar, dim, dim> jac;
127 return jac.determinant();
133 const GlobalPosition&
corner(
unsigned cornerIdx)
const
134 {
return corners_[cornerIdx]; }
140 Scalar
cornerWeight(
const LocalPosition& localPos,
unsigned cornerIdx)
const
145 for (
unsigned j = 0; j < dim; ++j) {
146 weight *= (cornerIdx & (1 << j)) ? localPos[j] : (1 - localPos[j]);
153 std::array<GlobalPosition, numCorners> corners_{};
154 GlobalPosition center_;
Quadrature geometry for quadrilaterals.
Definition quadraturegeometries.hh:44
Scalar integrationElement(const LocalPosition &localPos) const
Return the determinant of the Jacobian of the mapping from local to global coordinates at a given loc...
Definition quadraturegeometries.hh:123
Scalar cornerWeight(const LocalPosition &localPos, unsigned cornerIdx) const
Return the weight of an individual corner for the local to global mapping.
Definition quadraturegeometries.hh:140
const GlobalPosition & corner(unsigned cornerIdx) const
Return the position of the corner with a given index.
Definition quadraturegeometries.hh:133
GlobalPosition global(const LocalPosition &localPos) const
Convert a local coordinate into a global one.
Definition quadraturegeometries.hh:79
const GlobalPosition & center() const
Returns the center of weight of the polyhedron.
Definition quadraturegeometries.hh:73
void jacobian(Dune::FieldMatrix< Scalar, dim, dim > &jac, const LocalPosition &localPos) const
Returns the Jacobian matrix of the local to global mapping at a given local position.
Definition quadraturegeometries.hh:95
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:45