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 This file contains free functions to evaluate the transmissibilities | ||
11 | * associated with flux evaluations across sub-control volume faces | ||
12 | * in the context of cell-centered Mpfa schemes. | ||
13 | */ | ||
14 | #ifndef DUMUX_DISCRETIZATION_CC_MPFA_COMPUTE_TRANSMISSIBILITY_HH | ||
15 | #define DUMUX_DISCRETIZATION_CC_MPFA_COMPUTE_TRANSMISSIBILITY_HH | ||
16 | |||
17 | #include <dune/common/typetraits.hh> | ||
18 | #include <dune/common/fvector.hh> | ||
19 | |||
20 | #include <dumux/common/math.hh> | ||
21 | #include <dumux/discretization/extrusion.hh> | ||
22 | |||
23 | namespace Dumux { | ||
24 | |||
25 | /*! | ||
26 | * \ingroup CCMpfaDiscretization | ||
27 | * \brief Free function to evaluate the Mpfa transmissibility associated | ||
28 | * with the flux (in the form of flux = t*gradU) across a | ||
29 | * sub-control volume face stemming from a given sub-control | ||
30 | * volume with corresponding tensor t. | ||
31 | * | ||
32 | * \param fvGeometry The element-centered control volume geometry | ||
33 | * \param scv The iv-local sub-control volume | ||
34 | * \param scvf The grid sub-control volume face | ||
35 | * \param t The tensor living in the scv | ||
36 | * \param extrusionFactor The extrusion factor of the scv | ||
37 | */ | ||
38 | template<class EG, class IVSubControlVolume, class Tensor> | ||
39 | Dune::FieldVector<typename Tensor::field_type, IVSubControlVolume::myDimension> | ||
40 | 3054080 | computeMpfaTransmissibility(const EG& fvGeometry, | |
41 | const IVSubControlVolume& scv, | ||
42 | const typename EG::SubControlVolumeFace& scvf, | ||
43 | const Tensor& t, | ||
44 | typename IVSubControlVolume::ctype extrusionFactor) | ||
45 | { | ||
46 | 3054080 | Dune::FieldVector<typename Tensor::field_type, IVSubControlVolume::myDimension> wijk; | |
47 |
2/2✓ Branch 0 taken 6108160 times.
✓ Branch 1 taken 3054080 times.
|
9162240 | for (unsigned int dir = 0; dir < IVSubControlVolume::myDimension; ++dir) |
48 | 30540800 | wijk[dir] = vtmv(scvf.unitOuterNormal(), t, scv.nu(dir)); | |
49 | |||
50 | using Extrusion = Extrusion_t<typename EG::GridGeometry>; | ||
51 | 6108160 | wijk *= Extrusion::area(fvGeometry, scvf)*extrusionFactor; | |
52 | 3054080 | wijk /= scv.detX(); | |
53 | |||
54 | 3054080 | return wijk; | |
55 | } | ||
56 | |||
57 | /*! | ||
58 | * \ingroup CCMpfaDiscretization | ||
59 | * \brief Free function to evaluate the Mpfa transmissibility associated | ||
60 | * with the flux (in the form of flux = t*gradU) across a | ||
61 | * sub-control volume face stemming from a given sub-control | ||
62 | * volume with corresponding tensor t, where t is a scalar. | ||
63 | * | ||
64 | * \param fvGeometry The element-centered control volume geometry | ||
65 | * \param scv The iv-local sub-control volume | ||
66 | * \param scvf The grid sub-control volume face | ||
67 | * \param t the scalar quantity living in the scv | ||
68 | * \param extrusionFactor The extrusion factor of the scv | ||
69 | */ | ||
70 | template< class EG, class IVSubControlVolume, class Tensor, std::enable_if_t< Dune::IsNumber<Tensor>::value, int > = 1 > | ||
71 | Dune::FieldVector<Tensor, IVSubControlVolume::myDimension> | ||
72 | 404860322 | computeMpfaTransmissibility(const EG& fvGeometry, | |
73 | const IVSubControlVolume& scv, | ||
74 | const typename EG::SubControlVolumeFace& scvf, | ||
75 | const Tensor& t, | ||
76 | typename IVSubControlVolume::ctype extrusionFactor) | ||
77 | { | ||
78 | 404860322 | Dune::FieldVector<Tensor, IVSubControlVolume::myDimension> wijk; | |
79 |
2/2✓ Branch 0 taken 809720644 times.
✓ Branch 1 taken 404860322 times.
|
1214580966 | for (unsigned int dir = 0; dir < IVSubControlVolume::myDimension; ++dir) |
80 | 4048603220 | wijk[dir] = vtmv(scvf.unitOuterNormal(), t, scv.nu(dir)); | |
81 | |||
82 | using Extrusion = Extrusion_t<typename EG::GridGeometry>; | ||
83 | 809720644 | wijk *= Extrusion::area(fvGeometry, scvf)*extrusionFactor; | |
84 | 404860322 | wijk /= scv.detX(); | |
85 | |||
86 | 404860322 | return wijk; | |
87 | } | ||
88 | |||
89 | } // end namespace Dumux | ||
90 | |||
91 | #endif | ||
92 |