quadraturegeometries.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
27#ifndef EWOMS_QUADRATURE_GEOMETRIES_HH
28#define EWOMS_QUADRATURE_GEOMETRIES_HH
29
30#include <dune/common/fmatrix.hh>
31#include <dune/common/fvector.hh>
32#include <dune/geometry/type.hh>
33
34#include <array>
35#include <cassert>
36
37namespace Opm {
38
42template <class Scalar, unsigned dim>
44{
45public:
46 enum { numCorners = (1 << dim) };
47
48 using LocalPosition = Dune::FieldVector<Scalar, dim>;
49 using GlobalPosition = Dune::FieldVector<Scalar, dim>;
50
51 Dune::GeometryType type() const
52 { return Dune::GeometryType(/*topologyId=*/(1 << dim) - 1, dim); }
53
54 template <class CornerContainer>
55 void setCorners(const CornerContainer& corners, unsigned nCorners)
56 {
57 for (unsigned cornerIdx = 0; cornerIdx < nCorners; ++cornerIdx) {
58 for (unsigned j = 0; j < dim; ++j) {
59 corners_[cornerIdx][j] = corners[cornerIdx][j];
60 }
61 }
62
63 center_ = 0;
64 for (unsigned cornerIdx = 0; cornerIdx < nCorners; ++cornerIdx) {
65 center_ += corners_[cornerIdx];
66 }
67 center_ /= nCorners;
68 }
69
73 const GlobalPosition& center() const
74 { return center_; }
75
79 GlobalPosition global(const LocalPosition& localPos) const
80 {
81 GlobalPosition globalPos(0.0);
82
83 for (unsigned cornerIdx = 0; cornerIdx < numCorners; ++cornerIdx) {
84 globalPos.axpy(cornerWeight(localPos, cornerIdx),
85 corners_[cornerIdx]);
86 }
87
88 return globalPos;
89 }
90
95 void jacobian(Dune::FieldMatrix<Scalar, dim, dim>& jac,
96 const LocalPosition& localPos) const
97 {
98 jac = 0.0;
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) {
103 if (k != j) {
104 if (cornerIdx & (1 << j)) {
105 dWeight_dk *= localPos[j];
106 }
107 else {
108 dWeight_dk *= 1 - localPos[j];
109 }
110 }
111 }
112
113 jac[k].axpy(dWeight_dk, corners_[cornerIdx]);
114 }
115 }
116 }
117
123 Scalar integrationElement(const LocalPosition& localPos) const
124 {
125 Dune::FieldMatrix<Scalar, dim, dim> jac;
126 jacobian(jac, localPos);
127 return jac.determinant();
128 }
129
133 const GlobalPosition& corner(unsigned cornerIdx) const
134 { return corners_[cornerIdx]; }
135
140 Scalar cornerWeight(const LocalPosition& localPos, unsigned cornerIdx) const
141 {
142 GlobalPosition globalPos(0.0);
143
144 // this code is based on the Q1 finite element code from
145 // dune-localfunctions
146 Scalar weight = 1.0;
147 for (unsigned j = 0; j < dim; ++j) {
148 weight *= (cornerIdx & (1 << j)) ? localPos[j] : (1 - localPos[j]);
149 }
150
151 return weight;
152 }
153
154private:
155 std::array<GlobalPosition, numCorners> corners_{};
156 GlobalPosition center_;
157};
158
159} // namespace Opm
160
161#endif // EWOMS_QUADRATURE_GEOMETRY_HH
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
Dune::GeometryType type() const
Definition: quadraturegeometries.hh:51
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
Dune::FieldVector< Scalar, dim > LocalPosition
Definition: quadraturegeometries.hh:48
const GlobalPosition & center() const
Returns the center of weight of the polyhedron.
Definition: quadraturegeometries.hh:73
void setCorners(const CornerContainer &corners, unsigned nCorners)
Definition: quadraturegeometries.hh:55
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
Dune::FieldVector< Scalar, dim > GlobalPosition
Definition: quadraturegeometries.hh:49
@ numCorners
Definition: quadraturegeometries.hh:46
static constexpr int dim
Definition: structuredgridvanguard.hh:68
Definition: blackoilboundaryratevector.hh:39