version 3.10-dev
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
14#ifndef DUMUX_DISCRETIZATION_CC_MPFA_O_SCV_GEOMETRY_HELPER_HH
15#define DUMUX_DISCRETIZATION_CC_MPFA_O_SCV_GEOMETRY_HELPER_HH
16
17#include <array>
18
19#include <dune/common/exceptions.hh>
20#include <dune/geometry/type.hh>
21#include <dune/geometry/multilineargeometry.hh>
22
23namespace Dumux {
24
30template<class LocalScvType>
32{
33 using ctype = typename LocalScvType::ctype;
34 using LocalIndexType = typename LocalScvType::LocalIndexType;
35
36 static constexpr int dim = LocalScvType::myDimension;
37 static constexpr int dimWorld = LocalScvType::worldDimension;
38
39 struct MLGTraits : public Dune::MultiLinearGeometryTraits<ctype>
40 {
41 // we know the number of corners is always (2^(dim) corners (1<<dim))
42 template< int mydim, int cdim >
44 { using Type = std::array< typename LocalScvType::GlobalCoordinate, (1<<dim) >; };
45
46 // we know all scvs will have the same geometry type
47 template< int d >
49 {
50 static const bool v = true;
51 static const unsigned int topologyId = Dune::GeometryTypes::cube(d).id();
52 };
53 };
54
55public:
57 using ScvGeometry = Dune::MultiLinearGeometry<ctype, dim, dimWorld, MLGTraits>;
58
60 template<class InteractionVolume, class FVElementGeometry>
61 static ScvGeometry computeScvGeometry(LocalIndexType ivLocalScvIdx,
62 const InteractionVolume& iv,
63 const FVElementGeometry& fvGeometry)
64 {
65 const auto& scv = iv.localScv(ivLocalScvIdx);
66
67 if (dim == 2)
68 {
69 const auto& firstGridScvf = fvGeometry.scvf(iv.localScvf(scv.localScvfIndex(0)).gridScvfIndex());
70 const auto& secondGridScvf = fvGeometry.scvf(iv.localScvf(scv.localScvfIndex(1)).gridScvfIndex());
71
72 typename MLGTraits::template CornerStorage<dim, dimWorld>::Type corners;
73 corners[0] = fvGeometry.scv( scv.gridScvIndex() ).center();
74 corners[1] = fvGeometry.facetCorner(firstGridScvf);
75 corners[2] = fvGeometry.facetCorner(secondGridScvf);
76 corners[3] = fvGeometry.vertexCorner(secondGridScvf);
77
78 using std::swap;
79 typename LocalScvType::LocalBasis basis;
80 basis[0] = corners[1] - corners[0];
81 basis[1] = corners[2] - corners[0];
82 if ( !fvGeometry.gridGeometry().mpfaHelper().isRightHandSystem(basis) )
83 swap(corners[1], corners[2]);
84
85 return ScvGeometry(Dune::GeometryTypes::cube(ScvGeometry::mydimension), corners);
86 }
87 else if (dim == 3)
88 DUNE_THROW(Dune::NotImplemented, "Mpfa-o local scv geometry computation in 3d");
89 else
90 DUNE_THROW(Dune::InvalidStateException, "Mpfa only works in 2d or 3d");
91 }
92};
93
94} // end namespace Dumux
95
96#endif
Helper class providing functionality to compute the geometry of the interaction-volume local sub-cont...
Definition: scvgeometryhelper.hh:32
static ScvGeometry computeScvGeometry(LocalIndexType ivLocalScvIdx, const InteractionVolume &iv, const FVElementGeometry &fvGeometry)
returns the geometry of the i-th local scv
Definition: scvgeometryhelper.hh:61
Dune::MultiLinearGeometry< ctype, dim, dimWorld, MLGTraits > ScvGeometry
export the geometry type of the local scvs
Definition: scvgeometryhelper.hh:57
Definition: adapt.hh:17
std::array< typename LocalScvType::GlobalCoordinate,(1<< dim) > Type
Definition: scvgeometryhelper.hh:44
static const bool v
Definition: scvgeometryhelper.hh:50
static const unsigned int topologyId
Definition: scvgeometryhelper.hh:51