Line | Branch | Exec | Source |
---|---|---|---|
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 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup CCMpfaDiscretization | ||
10 | * \brief Helper class providing functionality to compute the geometry | ||
11 | * of the interaction-volume local sub-control volumes involved | ||
12 | * in the mpfa-o scheme. | ||
13 | */ | ||
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 | |||
23 | namespace Dumux { | ||
24 | |||
25 | /*! | ||
26 | * \ingroup CCMpfaDiscretization | ||
27 | * \brief Helper class providing functionality to compute the geometry | ||
28 | * of the interaction-volume local sub-control volumes of mpfa-o type. | ||
29 | */ | ||
30 | template<class LocalScvType> | ||
31 | class CCMpfaOScvGeometryHelper | ||
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 > | ||
43 | struct CornerStorage | ||
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 > | ||
48 | struct hasSingleGeometryType | ||
49 | { | ||
50 | static const bool v = true; | ||
51 | static const unsigned int topologyId = Dune::GeometryTypes::cube(d).id(); | ||
52 | }; | ||
53 | }; | ||
54 | |||
55 | public: | ||
56 | //! export the geometry type of the local scvs | ||
57 | using ScvGeometry = Dune::MultiLinearGeometry<ctype, dim, dimWorld, MLGTraits>; | ||
58 | |||
59 | //! returns the geometry of the i-th local scv | ||
60 | template<class InteractionVolume, class FVElementGeometry> | ||
61 | 400 | static ScvGeometry computeScvGeometry(LocalIndexType ivLocalScvIdx, | |
62 | const InteractionVolume& iv, | ||
63 | const FVElementGeometry& fvGeometry) | ||
64 | { | ||
65 | 400 | const auto& scv = iv.localScv(ivLocalScvIdx); | |
66 | |||
67 | if (dim == 2) | ||
68 | { | ||
69 | 400 | const auto& firstGridScvf = fvGeometry.scvf(iv.localScvf(scv.localScvfIndex(0)).gridScvfIndex()); | |
70 | 400 | const auto& secondGridScvf = fvGeometry.scvf(iv.localScvf(scv.localScvfIndex(1)).gridScvfIndex()); | |
71 | |||
72 | 400 | typename MLGTraits::template CornerStorage<dim, dimWorld>::Type corners; | |
73 |
2/2✓ Branch 0 taken 121 times.
✓ Branch 1 taken 279 times.
|
800 | corners[0] = fvGeometry.scv( scv.gridScvIndex() ).center(); |
74 | 400 | corners[1] = fvGeometry.facetCorner(firstGridScvf); | |
75 | 400 | corners[2] = fvGeometry.facetCorner(secondGridScvf); | |
76 | 400 | corners[3] = fvGeometry.vertexCorner(secondGridScvf); | |
77 | |||
78 | using std::swap; | ||
79 | 400 | typename LocalScvType::LocalBasis basis; | |
80 | 1600 | basis[0] = corners[1] - corners[0]; | |
81 |
2/2✓ Branch 0 taken 200 times.
✓ Branch 1 taken 200 times.
|
1600 | basis[1] = corners[2] - corners[0]; |
82 |
4/4✓ Branch 0 taken 200 times.
✓ Branch 1 taken 200 times.
✓ Branch 2 taken 200 times.
✓ Branch 3 taken 200 times.
|
800 | if ( !fvGeometry.gridGeometry().mpfaHelper().isRightHandSystem(basis) ) |
83 | 600 | swap(corners[1], corners[2]); | |
84 | |||
85 | 400 | 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 | ||
97 |