3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
scvgeometryhelper.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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program 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 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program 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 this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
26#ifndef DUMUX_DISCRETIZATION_CC_MPFA_O_SCV_GEOMETRY_HELPER_HH
27#define DUMUX_DISCRETIZATION_CC_MPFA_O_SCV_GEOMETRY_HELPER_HH
28
29#include <array>
30
31#include <dune/common/exceptions.hh>
32#include <dune/geometry/type.hh>
33#include <dune/geometry/multilineargeometry.hh>
34
35namespace Dumux {
36
42template<class LocalScvType>
44{
45 using ctype = typename LocalScvType::ctype;
46 using LocalIndexType = typename LocalScvType::LocalIndexType;
47
48 static constexpr int dim = LocalScvType::myDimension;
49 static constexpr int dimWorld = LocalScvType::worldDimension;
50
51 struct MLGTraits : public Dune::MultiLinearGeometryTraits<ctype>
52 {
53 // we know the number of corners is always (2^(dim) corners (1<<dim))
54 template< int mydim, int cdim >
56 { using Type = std::array< typename LocalScvType::GlobalCoordinate, (1<<dim) >; };
57
58 // we know all scvs will have the same geometry type
59 template< int d >
61 {
62 static const bool v = true;
63 static const unsigned int topologyId = Dune::GeometryTypes::cube(d).id();
64 };
65 };
66
67public:
69 using ScvGeometry = Dune::MultiLinearGeometry<ctype, dim, dimWorld, MLGTraits>;
70
72 template<class InteractionVolume, class FVElementGeometry>
73 static ScvGeometry computeScvGeometry(LocalIndexType ivLocalScvIdx,
74 const InteractionVolume& iv,
75 const FVElementGeometry& fvGeometry)
76 {
77 const auto& scv = iv.localScv(ivLocalScvIdx);
78
79 if (dim == 2)
80 {
81 const auto& firstGridScvf = fvGeometry.scvf(iv.localScvf(scv.localScvfIndex(0)).gridScvfIndex());
82 const auto& secondGridScvf = fvGeometry.scvf(iv.localScvf(scv.localScvfIndex(1)).gridScvfIndex());
83
84 typename MLGTraits::template CornerStorage<dim, dimWorld>::Type corners;
85 corners[0] = fvGeometry.scv( scv.gridScvIndex() ).center();
86 corners[1] = firstGridScvf.facetCorner();
87 corners[2] = secondGridScvf.facetCorner();
88 corners[3] = secondGridScvf.vertexCorner();
89
90 using std::swap;
91 typename LocalScvType::LocalBasis basis;
92 basis[0] = corners[1] - corners[0];
93 basis[1] = corners[2] - corners[0];
94 if ( !fvGeometry.gridGeometry().mpfaHelper().isRightHandSystem(basis) )
95 swap(corners[1], corners[2]);
96
97 return ScvGeometry(Dune::GeometryTypes::cube(ScvGeometry::mydimension), corners);
98 }
99 else if (dim == 3)
100 DUNE_THROW(Dune::NotImplemented, "Mpfa-o local scv geometry computation in 3d");
101 else
102 DUNE_THROW(Dune::InvalidStateException, "Mpfa only works in 2d or 3d");
103 }
104};
105
106} // end namespace Dumux
107
108#endif
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
Helper class providing functionality to compute the geometry of the interaction-volume local sub-cont...
Definition: scvgeometryhelper.hh:44
static ScvGeometry computeScvGeometry(LocalIndexType ivLocalScvIdx, const InteractionVolume &iv, const FVElementGeometry &fvGeometry)
returns the geometry of the i-th local scv
Definition: scvgeometryhelper.hh:73
Dune::MultiLinearGeometry< ctype, dim, dimWorld, MLGTraits > ScvGeometry
export the geometry type of the local scvs
Definition: scvgeometryhelper.hh:69
std::array< typename LocalScvType::GlobalCoordinate,(1<< dim) > Type
Definition: scvgeometryhelper.hh:56
static const bool v
Definition: scvgeometryhelper.hh:62
static const unsigned int topologyId
Definition: scvgeometryhelper.hh:63