37#ifndef OPM_GEOMETRY_HEADER
38#define OPM_GEOMETRY_HEADER
43#include <opm/grid/utility/platform_dependent/disable_warnings.h>
45#include <dune/common/version.hh>
46#include <dune/geometry/referenceelements.hh>
47#include <dune/grid/common/geometry.hh>
49#include <dune/geometry/type.hh>
51#include <opm/grid/cpgrid/EntityRep.hpp>
52#include <opm/grid/cpgrid/DefaultGeometryPolicy.hpp>
53#include <opm/grid/cpgrid/OrientedEntityTable.hpp>
55#include <opm/grid/common/Volumes.hpp>
56#include <opm/grid/utility/platform_dependent/reenable_warnings.h>
57#include <opm/grid/utility/SparseTable.hpp>
59#include <opm/grid/utility/ErrorMacros.hpp>
74 template <
int mydim,
int cdim>
86 static_assert(cdim == 3,
"");
89 enum { dimension = 3 };
91 enum { mydimension = 0};
93 enum { coorddimension = cdim };
95 enum { dimensionworld = 3 };
108 typedef FieldMatrix< ctype, coorddimension, mydimension >
Jacobian;
152 return Dune::GeometryTypes::cube(mydimension);
164 static_cast<void>(cor);
191 JacobianInverseTransposed
219 GlobalCoordinate pos_;
230 static_assert(cdim == 3,
"");
233 enum { dimension = 3 };
235 enum { mydimension = 2 };
237 enum { coorddimension = cdim };
239 enum { dimensionworld = 3 };
252 typedef FieldMatrix< ctype, coorddimension, mydimension >
Jacobian;
265 : pos_(pos), vol_(vol)
271 : pos_(0.0), vol_(0.0)
278 OPM_THROW(std::runtime_error,
"Geometry::global() meaningless on singular geometry.");
284 OPM_THROW(std::runtime_error,
"Geometry::local() meaningless on singular geometry.");
297 return Dune::GeometryTypes::none(mydimension);
329 const FieldMatrix<ctype, mydimension, coorddimension>&
332 OPM_THROW(std::runtime_error,
"Meaningless to call jacobianTransposed() on singular geometries.");
336 const FieldMatrix<ctype, coorddimension, mydimension>&
339 OPM_THROW(std::runtime_error,
"Meaningless to call jacobianInverseTransposed() on singular geometries.");
363 GlobalCoordinate pos_;
375 static_assert(cdim == 3,
"");
378 enum { dimension = 3 };
380 enum { mydimension = 3 };
382 enum { coorddimension = cdim };
384 enum { dimensionworld = 3 };
397 typedef FieldMatrix< ctype, coorddimension, mydimension >
Jacobian;
405 typedef Dune::Impl::FieldMatrixHelper< double > MatrixHelperType;
418 const int* corner_indices)
419 : pos_(pos), vol_(vol),
420 allcorners_(allcorners_ptr), cor_idx_(corner_indices)
422 assert(allcorners_ && corner_indices);
427 : pos_(0.0), vol_(0.0), allcorners_(0), cor_idx_(0)
438 static_assert(mydimension == 3,
"");
439 static_assert(coorddimension == 3,
"");
442 uvw[0] -= local_coord;
444 const int pat[8][3] = { { 0, 0, 0 },
453 for (
int i = 0; i < 8; ++i) {
456 for (
int j = 0; j < 3; ++j) {
457 factor *= uvw[pat[i][j]][j];
459 corner_contrib *= factor;
460 xyz += corner_contrib;
469 static_assert(mydimension == 3,
"");
470 static_assert(coorddimension == 3,
"");
473 const ctype epsilon = 1e-12;
474 auto refElement = Dune::ReferenceElements<ctype, 3>::cube();
482 MatrixHelperType::template xTRightInvA<3, 3>(JT, z, dx );
484 }
while (dx.two_norm2() > epsilon*epsilon);
495 return MatrixHelperType::template sqrtDetAAT<3, 3>(Jt);
502 return Dune::GeometryTypes::cube(mydimension);
515 assert(allcorners_ && cor_idx_);
516 return (allcorners_->data())[cor_idx_[cor]].center();
525 void set_volume(ctype
volume) {
541 const JacobianTransposed
544 static_assert(mydimension == 3,
"");
545 static_assert(coorddimension == 3,
"");
549 uvw[0] -= local_coord;
551 const int pat[8][3] = { { 0, 0, 0 },
560 for (
int i = 0; i < 8; ++i) {
561 for (
int deriv = 0; deriv < 3; ++deriv) {
564 for (
int j = 0; j < 3; ++j) {
565 factor *= (j != deriv) ? uvw[pat[i][j]][j]
566 : (pat[i][j] == 0 ? -1.0 : 1.0);
569 corner_contrib *= factor;
570 Jt[deriv] += corner_contrib;
577 const JacobianInverseTransposed
626 std::vector<std::array<int,8>>& refined_cell_to_point,
632 const std::array<int,3>& patch_dim,
633 const std::vector<double>& widthsX,
634 const std::vector<double>& lengthsY,
635 const std::vector<double>& heightsZ)
const
638 *(all_geom.
geomVector(std::integral_constant<int,3>()));
640 *(all_geom.
geomVector(std::integral_constant<int,1>()));
642 *(all_geom.
geomVector(std::integral_constant<int,0>()));
652 const std::array<int,3>& refined_dim = { cells_per_dim[0]*patch_dim[0],
653 cells_per_dim[1]*patch_dim[1],
654 cells_per_dim[2]*patch_dim[2]};
655 refined_corners.resize((refined_dim[0] + 1)*(refined_dim[1] + 1)*(refined_dim[2] + 1));
662 assert(
static_cast<int>(widthsX.size()) == patch_dim[0]);
663 assert(
static_cast<int>(lengthsY.size()) == patch_dim[1]);
664 assert(
static_cast<int>(heightsZ.size()) == patch_dim[2]);
665 const auto localCoordNumerator = [](
const std::vector<double>& vec,
int sumLimit,
double multiplier) {
667 assert(!vec.empty());
668 assert(sumLimit <
static_cast<int>(vec.size()));
669 lcn += multiplier*vec[sumLimit];
670 for (
int m = 0; m < sumLimit; ++m) {
677 const double sumWidths = std::accumulate(widthsX.begin(), widthsX.end(),
double(0));
679 const double sumLengths = std::accumulate(lengthsY.begin(), lengthsY.end(),
double(0));
681 const double sumHeights = std::accumulate(heightsZ.begin(), heightsZ.end(),
double(0));
684 for (
int j = 0; j < refined_dim[1] +1; ++j) {
686 for (
int i = 0; i < refined_dim[0] +1; ++i) {
688 for (
int k = 0; k < refined_dim[2] +1; ++k) {
692 int refined_corner_idx =
693 (j*(refined_dim[0]+1)*(refined_dim[2]+1)) + (i*(refined_dim[2]+1)) + k;
696 if ( i == refined_dim[0]) {
699 local_x = localCoordNumerator(widthsX, i/cells_per_dim[0],
double((i % cells_per_dim[0])) / cells_per_dim[0]);
701 if ( j == refined_dim[1]) {
702 local_y = sumLengths;
704 local_y = localCoordNumerator(lengthsY, j/cells_per_dim[1],
double((j % cells_per_dim[1])) / cells_per_dim[1]);
706 if ( k == refined_dim[2]) {
707 local_z = sumHeights;
709 local_z = localCoordNumerator(heightsZ, k/cells_per_dim[2],
double((k % cells_per_dim[2])) /cells_per_dim[2]);
712 const LocalCoordinate& local_refined_corner = { local_x/sumWidths, local_y/sumLengths, local_z/sumHeights };
713 assert(local_x/sumWidths <= 1.);
714 assert(local_y/sumLengths <= 1.);
715 assert(local_z/sumHeights <= 1.);
726 const int refined_faces_size =
727 (refined_dim[0]*refined_dim[1]*(refined_dim[2]+1))
728 + ((refined_dim[0]+1)*refined_dim[1]*refined_dim[2])
729 + (refined_dim[0]*(refined_dim[1]+1)*refined_dim[2]);
730 refined_faces.resize(refined_faces_size);
731 refined_face_tags.resize(refined_faces_size);
732 refined_face_normals.resize(refined_faces_size);
769 for (
int constant_direction = 0; constant_direction < 3; ++constant_direction){
774 std::array<int,3> refined_dim_mixed = {
775 refined_dim[(2+constant_direction)%3],
776 refined_dim[(1+constant_direction)%3],
777 refined_dim[constant_direction % 3] };
778 for (
int l = 0; l < refined_dim_mixed[0] + 1; ++l) {
779 for (
int m = 0; m < refined_dim_mixed[1]; ++m) {
780 for (
int n = 0; n < refined_dim_mixed[2]; ++n) {
782 auto [
face_tag, idx, face_to_point, face_to_cell, wrong_local_centroid] =
783 getIndicesFace(l, m, n, constant_direction, refined_dim);
787 refined_face_to_point.
appendRow(face_to_point.begin(), face_to_point.end());
789 refined_face_to_cell.
appendRow(face_to_cell.begin(), face_to_cell.end());
792 for (
int corn = 0; corn < 4; ++corn){
793 face_center += refined_corners[face_to_point[corn]].center();
800 GlobalCoordinate face_vector0 = refined_corners[face_to_point[0]].center() - face_center;
801 GlobalCoordinate face_vector1 = refined_corners[face_to_point[1]].center() - face_center;
802 mutable_face_normals[idx] = {
803 (face_vector0[1]*face_vector1[2]) - (face_vector0[2]*face_vector1[1]),
804 (face_vector0[2]*face_vector1[0]) - (face_vector0[0]*face_vector1[2]),
805 (face_vector0[0]*face_vector1[1]) - (face_vector0[1]*face_vector1[0])};
806 mutable_face_normals[idx] /= mutable_face_normals[idx].two_norm();
808 mutable_face_normals[idx] *= -1;
812 std::vector<std::array<int,2>> refined_face_to_edges = {
813 { face_to_point[0], face_to_point[1] },
814 { face_to_point[0], face_to_point[2] },
815 { face_to_point[1], face_to_point[3] },
816 { face_to_point[2], face_to_point[3] }
820 double refined_face_area = 0.0;
821 for (
int edge = 0; edge < 4; ++edge) {
825 refined_corners[refined_face_to_edges[edge][0]].center(),
826 refined_corners[refined_face_to_edges[edge][1]].center(),
828 refined_face_area += std::fabs(
area(trian_corners));
849 refined_cells.resize(refined_dim[0] * refined_dim[1] * refined_dim[2]);
852 refined_cell_to_point.resize(refined_dim[0] * refined_dim[1] * refined_dim[2]);
881 double sum_all_refined_cell_volumes = 0.0;
894 for (
int k = 0; k < refined_dim[2]; ++k) {
895 for (
int j = 0; j < refined_dim[1]; ++j) {
896 for (
int i = 0; i < refined_dim[0]; ++i) {
898 int refined_cell_idx = (k*refined_dim[0]*refined_dim[1]) + (j*refined_dim[0]) +i;
901 double refined_cell_volume = 0.0;
904 std::array<int,8> cell_to_point = {
905 (j*(refined_dim[0]+1)*(refined_dim[2]+1)) + (i*(refined_dim[2]+1)) +k,
906 (j*(refined_dim[0]+1)*(refined_dim[2]+1)) + ((i+1)*(refined_dim[2]+1)) +k,
907 ((j+1)*(refined_dim[0]+1)*(refined_dim[2]+1)) + (i*(refined_dim[2]+1)) +k,
908 ((j+1)*(refined_dim[0]+1)*(refined_dim[2]+1)) + ((i+1)*(refined_dim[2]+1)) +k,
909 (j*(refined_dim[0]+1)*(refined_dim[2]+1)) + (i*(refined_dim[2]+1)) +k+1,
910 (j*(refined_dim[0]+1)*(refined_dim[2]+1)) + ((i+1)*(refined_dim[2]+1)) +k+1,
911 ((j+1)*(refined_dim[0]+1)*(refined_dim[2]+1)) + (i*(refined_dim[2]+1)) +k+1,
912 ((j+1)*(refined_dim[0]+1)*(refined_dim[2]+1)) + ((i+1)*(refined_dim[2]+1)) +k+1
915 refined_cell_to_point[refined_cell_idx] = cell_to_point;
919 for (
int corn = 0; corn < 8; ++corn) {
920 refined_cell_center += refined_corners[cell_to_point[corn]].center();
922 refined_cell_center /= 8.;
926 std::vector<int> hexa_to_face = {
928 (k*refined_dim[0]*refined_dim[1]) + (j*refined_dim[0]) + i,
930 (refined_dim[0]*refined_dim[1]*(refined_dim[2]+1))
931 + ((refined_dim[0]+1)*refined_dim[1]*refined_dim[2])
932 + (j*refined_dim[0]*refined_dim[2]) + (i*refined_dim[2]) + k,
934 (refined_dim[0]*refined_dim[1]*(refined_dim[2]+1))
935 + (i*refined_dim[1]*refined_dim[2]) + (k*refined_dim[1]) + j,
937 (refined_dim[0]*refined_dim[1]*(refined_dim[2]+1))
938 + ((i+1)*refined_dim[1]*refined_dim[2]) + (k*refined_dim[1]) + j,
940 (refined_dim[0]*refined_dim[1]*(refined_dim[2]+1)) +
941 ((refined_dim[0]+1)*refined_dim[1]*refined_dim[2])
942 + ((j+1)*refined_dim[0]*refined_dim[2]) + (i*refined_dim[2]) + k,
944 ((k+1)*refined_dim[0]*refined_dim[1]) + (j*refined_dim[0]) + i};
954 std::vector<cpgrid::EntityRep<1>> cell_to_face = {
955 { hexa_to_face[0],
false}, {hexa_to_face[1],
false},
956 { hexa_to_face[2],
false}, {hexa_to_face[3],
true},
957 { hexa_to_face[4],
true}, {hexa_to_face[5],
true} };
958 refined_cell_to_face.
appendRow(cell_to_face.begin(), cell_to_face.end());
962 std::vector<Geometry<0,3>::GlobalCoordinate> hexa_face_centroids;
963 for (
auto& idx : hexa_to_face) {
964 hexa_face_centroids.push_back(refined_faces[idx].
center());
978 std::vector<std::array<int,4>> cell_face4corners;
979 cell_face4corners.reserve(6);
980 for (
int face = 0; face < 6; ++face) {
981 cell_face4corners.push_back({
982 refined_face_to_point[hexa_to_face[face]][0],
983 refined_face_to_point[hexa_to_face[face]][1],
984 refined_face_to_point[hexa_to_face[face]][2],
985 refined_face_to_point[hexa_to_face[face]][3] });
989 std::vector<std::vector<std::array<int,2>>> tetra_edge_indices;
990 tetra_edge_indices.reserve(6);
991 for (
auto& face_indices : cell_face4corners)
993 std::vector<std::array<int,2>> face4edges_indices = {
994 { face_indices[0], face_indices[1]},
995 { face_indices[0], face_indices[2]},
996 { face_indices[1], face_indices[3]},
997 { face_indices[2], face_indices[3]} };
998 tetra_edge_indices.push_back(face4edges_indices);
1004 for (
int face = 0; face < 6; ++face) {
1005 for (
int edge = 0; edge < 4; ++edge) {
1009 refined_corners[tetra_edge_indices[face][edge][0]].center(),
1010 refined_corners[tetra_edge_indices[face][edge][1]].center(),
1011 hexa_face_centroids[face],
1012 refined_cell_center };
1018 sum_all_refined_cell_volumes += refined_cell_volume;
1021 int* indices_storage_ptr = refined_cell_to_point[refined_cell_idx].data();
1023 refined_cells[refined_cell_idx] =
1025 refined_cell_volume,
1026 all_geom.
geomVector(std::integral_constant<int,3>()),
1027 indices_storage_ptr);
1034 if (std::fabs(sum_all_refined_cell_volumes - this->
volume())
1037 for(
auto& cell: refined_cells){
1038 cell.vol_ *= correction;
1045 GlobalCoordinate pos_;
1047 std::shared_ptr<const EntityVariable<Geometry<0, 3>,3>> allcorners_;
1048 const int* cor_idx_;
1063 const std::tuple<
enum face_tag, int,
1064 std::array<int, 4>, std::vector<cpgrid::EntityRep<0>>,
1066 getIndicesFace(
int l,
int m,
int n,
int constant_direction,
const std::array<int, 3>& cells_per_dim)
const
1069 std::vector<cpgrid::EntityRep<0>> neighboring_cells_of_one_face;
1070 switch(constant_direction) {
1075 neighboring_cells_of_one_face.push_back({((l-1)*cells_per_dim[0]*cells_per_dim[1])
1076 + (m*cells_per_dim[0]) + n,
true});
1078 if (l != cells_per_dim[2]) {
1079 neighboring_cells_of_one_face.push_back({ (l*cells_per_dim[0]*cells_per_dim[1])
1080 + (m*cells_per_dim[0]) + n,
false});
1082 return {
face_tag::K_FACE, (l*cells_per_dim[0]*cells_per_dim[1]) + (m*cells_per_dim[0]) + n,
1083 {(m*(cells_per_dim[0]+1)*(cells_per_dim[2]+1)) + (n*(cells_per_dim[2]+1)) +l,
1084 (m*(cells_per_dim[0]+1)*(cells_per_dim[2]+1)) + ((n+1)*(cells_per_dim[2]+1)) +l,
1085 ((m+1)*(cells_per_dim[0]+1)*(cells_per_dim[2]+1)) + (n*(cells_per_dim[2]+1)) +l,
1086 ((m+1)*(cells_per_dim[0]+1)*(cells_per_dim[2]+1)) + ((n+1)*(cells_per_dim[2]+1)) +l},
1087 neighboring_cells_of_one_face,
1088 {(.5 + n)/cells_per_dim[0], (.5 + m)/cells_per_dim[1], double(l)/cells_per_dim[2]}};
1093 neighboring_cells_of_one_face.push_back({(m*cells_per_dim[0]*cells_per_dim[1])
1094 + (n*cells_per_dim[0]) +l-1,
true});
1096 if (l != cells_per_dim[0]) {
1097 neighboring_cells_of_one_face.push_back({ (m*cells_per_dim[0]*cells_per_dim[1])
1098 + (n*cells_per_dim[0]) + l,
false});
1100 return {
face_tag::I_FACE, (cells_per_dim[0]*cells_per_dim[1]*(cells_per_dim[2]+1))
1101 + (l*cells_per_dim[1]*cells_per_dim[2]) + (m*cells_per_dim[1]) + n,
1102 {(n*(cells_per_dim[0]+1)*(cells_per_dim[2]+1)) + (l*(cells_per_dim[2]+1)) +m,
1103 ((n+1)*(cells_per_dim[0]+1)*(cells_per_dim[2]+1)) + (l*(cells_per_dim[2]+1)) +m,
1104 (n*(cells_per_dim[0]+1)*(cells_per_dim[2]+1)) + (l*(cells_per_dim[2]+1)) +m+1,
1105 ((n+1)*(cells_per_dim[0]+1)*(cells_per_dim[2]+1)) + (l*(cells_per_dim[2]+1)) +m+1},
1106 neighboring_cells_of_one_face,
1107 { double(l)/cells_per_dim[0], (.5 + n)/cells_per_dim[1], (.5 + m)/cells_per_dim[2]}};
1112 neighboring_cells_of_one_face.push_back({(n*cells_per_dim[0]*cells_per_dim[1])
1113 + ((l-1)*cells_per_dim[0]) +m,
true});
1115 if (l != cells_per_dim[1]) {
1116 neighboring_cells_of_one_face.push_back({(n*cells_per_dim[0]*cells_per_dim[1])
1117 + (l*cells_per_dim[0]) + m,
false});
1119 return {
face_tag::J_FACE, (cells_per_dim[0]*cells_per_dim[1]*(cells_per_dim[2] +1))
1120 + ((cells_per_dim[0]+1)*cells_per_dim[1]*cells_per_dim[2])
1121 + (l*cells_per_dim[0]*cells_per_dim[2]) + (m*cells_per_dim[2]) + n,
1122 {(l*(cells_per_dim[0]+1)*(cells_per_dim[2]+1)) + (m*(cells_per_dim[2]+1)) +n,
1123 (l*(cells_per_dim[0]+1)*(cells_per_dim[2]+1)) + ((m+1)*(cells_per_dim[2]+1)) +n,
1124 (l*(cells_per_dim[0]+1)*(cells_per_dim[2]+1)) + (m*(cells_per_dim[2]+1)) +n+1,
1125 (l*(cells_per_dim[0]+1)*(cells_per_dim[2]+1)) + ((m+1)*(cells_per_dim[2]+1)) +n+1},
1126 neighboring_cells_of_one_face,
1127 {(.5 + m)/cells_per_dim[0],
double(l)/cells_per_dim[1], (.5 + n)/cells_per_dim[2]}};
1130 OPM_THROW(std::logic_error,
"Unhandled dimension. This should never happen!");
1136 template<
int mydim,
int cdim >
1137 auto referenceElement(
const cpgrid::Geometry<mydim,cdim>& geo) ->
decltype(referenceElement<double,mydim>(geo.type()))
1139 return referenceElement<double,mydim>(geo.type());
Definition DefaultGeometryPolicy.hpp:53
const EntityVariable< cpgrid::Geometry< 3 - codim, 3 >, codim > & geomVector() const
Definition DefaultGeometryPolicy.hpp:86
Represents an entity of a given codim, with positive or negative orientation.
Definition EntityRep.hpp:98
Base class for EntityVariable and SignedEntityVariable.
Definition EntityRep.hpp:218
A class design to hold a variable with a value for each entity of the given codimension,...
Definition EntityRep.hpp:271
FieldVector< ctype, coorddimension > GlobalCoordinate
Range type of.
Definition Geometry.hpp:105
bool affine() const
The mapping implemented by this geometry is constant, therefore affine.
Definition Geometry.hpp:213
FieldMatrix< ctype, coorddimension, mydimension > Jacobian
Type of Jacobian matrix.
Definition Geometry.hpp:108
JacobianInverse jacobianInverse(const LocalCoordinate &) const
This method is meaningless for singular geometries.
Definition Geometry.hpp:207
Volume volume() const
Volume of vertex is arbitrarily set to 1.
Definition Geometry.hpp:170
FieldVector< ctype, mydimension > LocalCoordinate
Domain type of.
Definition Geometry.hpp:103
Jacobian jacobian(const LocalCoordinate &) const
This method is meaningless for singular geometries.
Definition Geometry.hpp:200
double ctype
Coordinate element type.
Definition Geometry.hpp:98
FieldMatrix< ctype, coorddimension, mydimension > JacobianInverseTransposed
Type of the inverse of the transposed Jacobian matrix.
Definition Geometry.hpp:114
Geometry(const GlobalCoordinate &pos)
Construct from vertex position.
Definition Geometry.hpp:119
FieldMatrix< ctype, coorddimension, mydimension > JacobianInverse
Type of inverse of Jacobian matrix.
Definition Geometry.hpp:110
ctype Volume
Number type used for the geometry volume.
Definition Geometry.hpp:100
JacobianTransposed jacobianTransposed(const LocalCoordinate &) const
This method is meaningless for singular geometries.
Definition Geometry.hpp:183
GeometryType type() const
Using the cube type for vertices.
Definition Geometry.hpp:150
FieldMatrix< ctype, mydimension, coorddimension > JacobianTransposed
Type of transposed Jacobian matrix.
Definition Geometry.hpp:112
GlobalCoordinate corner(int cor) const
Returns the single corner: the vertex itself.
Definition Geometry.hpp:162
int corners() const
A vertex is defined by a single corner.
Definition Geometry.hpp:156
const GlobalCoordinate & center() const
Returns the centroid of the geometry.
Definition Geometry.hpp:176
Geometry()
Default constructor, giving a non-valid geometry.
Definition Geometry.hpp:125
const GlobalCoordinate & global(const LocalCoordinate &) const
Returns the position of the vertex.
Definition Geometry.hpp:131
Volume integrationElement(const LocalCoordinate &) const
Returns 1 for the vertex geometry.
Definition Geometry.hpp:144
JacobianInverseTransposed jacobianInverseTransposed(const LocalCoordinate &) const
This method is meaningless for singular geometries.
Definition Geometry.hpp:192
LocalCoordinate local(const GlobalCoordinate &) const
Meaningless for the vertex geometry.
Definition Geometry.hpp:137
JacobianInverse jacobianInverse(const LocalCoordinate &) const
The inverse of the jacobian.
Definition Geometry.hpp:351
Volume volume() const
Volume (area, actually) of intersection.
Definition Geometry.hpp:317
const FieldMatrix< ctype, mydimension, coorddimension > & jacobianTransposed(const LocalCoordinate &) const
This method is meaningless for singular geometries.
Definition Geometry.hpp:330
const GlobalCoordinate & center() const
Returns the centroid of the geometry.
Definition Geometry.hpp:323
Volume integrationElement(const LocalCoordinate &) const
For the singular geometry, we return a constant integration element equal to the volume.
Definition Geometry.hpp:289
int corners() const
The number of corners of this convex polytope.
Definition Geometry.hpp:302
LocalCoordinate local(const GlobalCoordinate &) const
This method is meaningless for singular geometries.
Definition Geometry.hpp:282
double ctype
Coordinate element type.
Definition Geometry.hpp:242
FieldMatrix< ctype, coorddimension, mydimension > Jacobian
Type of Jacobian matrix.
Definition Geometry.hpp:252
bool affine() const
Since integrationElement() is constant, returns true.
Definition Geometry.hpp:357
ctype Volume
Number type used for the geometry volume.
Definition Geometry.hpp:244
FieldVector< ctype, mydimension > LocalCoordinate
Domain type of.
Definition Geometry.hpp:247
FieldMatrix< ctype, coorddimension, mydimension > JacobianInverseTransposed
Type of the inverse of the transposed Jacobian matrix.
Definition Geometry.hpp:258
const FieldMatrix< ctype, coorddimension, mydimension > & jacobianInverseTransposed(const LocalCoordinate &) const
This method is meaningless for singular geometries.
Definition Geometry.hpp:337
GeometryType type() const
We use the singular type (None) for intersections.
Definition Geometry.hpp:295
FieldMatrix< ctype, mydimension, coorddimension > JacobianTransposed
Type of transposed Jacobian matrix.
Definition Geometry.hpp:256
Geometry()
Default constructor, giving a non-valid geometry.
Definition Geometry.hpp:270
const GlobalCoordinate & global(const LocalCoordinate &) const
This method is meaningless for singular geometries.
Definition Geometry.hpp:276
FieldMatrix< ctype, coorddimension, mydimension > JacobianInverse
Type of inverse of Jacobian matrix.
Definition Geometry.hpp:254
Geometry(const GlobalCoordinate &pos, ctype vol)
Construct from centroid and volume (1- and 0-moments).
Definition Geometry.hpp:263
Jacobian jacobian(const LocalCoordinate &) const
The jacobian.
Definition Geometry.hpp:344
FieldVector< ctype, coorddimension > GlobalCoordinate
Range type of.
Definition Geometry.hpp:249
GlobalCoordinate corner(int) const
This method is meaningless for singular geometries.
Definition Geometry.hpp:308
const JacobianInverseTransposed jacobianInverseTransposed(const LocalCoordinate &local_coord) const
Inverse of Jacobian transposed.
Definition Geometry.hpp:578
FieldMatrix< ctype, coorddimension, mydimension > JacobianInverseTransposed
Type of the inverse of the transposed Jacobian matrix.
Definition Geometry.hpp:403
bool affine() const
The mapping implemented by this geometry is not generally affine.
Definition Geometry.hpp:600
Volume volume() const
Cell volume.
Definition Geometry.hpp:520
FieldVector< ctype, coorddimension > GlobalCoordinate
Range type of.
Definition Geometry.hpp:394
GlobalCoordinate corner(int cor) const
Get the cor-th of 8 corners of the hexahedral base cell.
Definition Geometry.hpp:513
FieldMatrix< ctype, mydimension, coorddimension > JacobianTransposed
Type of transposed Jacobian matrix.
Definition Geometry.hpp:401
Geometry(const GlobalCoordinate &pos, ctype vol, std::shared_ptr< const EntityVariable< cpgrid::Geometry< 0, 3 >, 3 > > allcorners_ptr, const int *corner_indices)
Construct from center, volume (1- and 0-moments) and corners.
Definition Geometry.hpp:415
GeometryType type() const
Using the cube type for all entities now (cells and vertices), but we use the singular type for inter...
Definition Geometry.hpp:500
void refineCellifiedPatch(const std::array< int, 3 > &cells_per_dim, DefaultGeometryPolicy &all_geom, std::vector< std::array< int, 8 > > &refined_cell_to_point, cpgrid::OrientedEntityTable< 0, 1 > &refined_cell_to_face, Opm::SparseTable< int > &refined_face_to_point, cpgrid::OrientedEntityTable< 1, 0 > &refined_face_to_cell, cpgrid::EntityVariable< enum face_tag, 1 > &refined_face_tags, cpgrid::SignedEntityVariable< PointType, 1 > &refined_face_normals, const std::array< int, 3 > &patch_dim, const std::vector< double > &widthsX, const std::vector< double > &lengthsY, const std::vector< double > &heightsZ) const
Definition Geometry.hpp:624
const JacobianTransposed jacobianTransposed(const LocalCoordinate &local_coord) const
Jacobian transposed.
Definition Geometry.hpp:542
GlobalCoordinate global(const LocalCoordinate &local_coord) const
Provide a trilinear mapping.
Definition Geometry.hpp:436
Geometry()
Default constructor, giving a non-valid geometry.
Definition Geometry.hpp:426
FieldMatrix< ctype, coorddimension, mydimension > JacobianInverse
Type of inverse of Jacobian matrix.
Definition Geometry.hpp:399
ctype Volume
Number type used for the geometry volume.
Definition Geometry.hpp:389
LocalCoordinate local(const GlobalCoordinate &y) const
Mapping from the cell to the reference domain.
Definition Geometry.hpp:467
int corners() const
The number of corners of this convex polytope.
Definition Geometry.hpp:507
const GlobalCoordinate & center() const
Returns the centroid of the geometry.
Definition Geometry.hpp:530
Jacobian jacobian(const LocalCoordinate &local_coord) const
The jacobian.
Definition Geometry.hpp:587
Dune::FieldVector< double, 3 > PointType
Refine a single cell considering different widths, lengths, and heights.
Definition Geometry.hpp:623
FieldVector< ctype, mydimension > LocalCoordinate
Domain type of.
Definition Geometry.hpp:392
Volume integrationElement(const LocalCoordinate &local_coord) const
Equal to \sqrt{\det{J^T J}} where J is the Jacobian.
Definition Geometry.hpp:492
JacobianInverse jacobianInverse(const LocalCoordinate &local_coord) const
The inverse of the jacobian.
Definition Geometry.hpp:594
FieldMatrix< ctype, coorddimension, mydimension > Jacobian
Type of Jacobian matrix.
Definition Geometry.hpp:397
double ctype
Coordinate element type.
Definition Geometry.hpp:387
This class encapsulates geometry for vertices, intersections, and cells.
Definition Geometry.hpp:76
Represents the topological relationships between sets of entities, for example cells and faces.
Definition OrientedEntityTable.hpp:139
void appendRow(DataIter row_beg, DataIter row_end)
Appends a row to the table.
Definition SparseTable.hpp:182
A class design to hold a variable with a value for each entity of the given codimension,...
Definition EntityRep.hpp:306
A SparseTable stores a table with rows of varying size as efficiently as possible.
Definition SparseTable.hpp:117
void appendRow(DataIter row_beg, DataIter row_end)
Appends a row to the table.
Definition SparseTable.hpp:182
Copyright 2019 Equinor AS.
Definition GridPartitioning.cpp:673
The namespace Dune is the main namespace for all Dune code.
Definition CartesianIndexMapper.hpp:10
T area(const Point< T, 2 > *c)
Computes the area of a 2-dimensional triangle.
Definition Volumes.hpp:119
T volume(const Point< T, 3 > *c)
Computes the volume of a 3D simplex (embedded i 3D space).
Definition Volumes.hpp:138
T simplex_volume(const Point< T, Dim > *a)
Computes the volume of a simplex consisting of (Dim+1) vertices embedded in Euclidean space of dimens...
Definition Volumes.hpp:104
Low-level corner-point processing routines and supporting data structures.
face_tag
Connection taxonomy.
Definition preprocess.h:66
@ K_FACE
Connection topologically normal to I-J plane.
Definition preprocess.h:69
@ J_FACE
Connection topologically normal to I-K plane.
Definition preprocess.h:68
@ I_FACE
Connection topologically normal to J-K plane.
Definition preprocess.h:67