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 LowReKEpsilonModel | ||
10 | * \copydoc Dumux::LowReKEpsilonResidualImpl | ||
11 | */ | ||
12 | #ifndef DUMUX_STAGGERED_LOWREKEPSILON_LOCAL_RESIDUAL_HH | ||
13 | #define DUMUX_STAGGERED_LOWREKEPSILON_LOCAL_RESIDUAL_HH | ||
14 | |||
15 | #include <dune/common/hybridutilities.hh> | ||
16 | #include <dumux/common/properties.hh> | ||
17 | #include <dumux/discretization/method.hh> | ||
18 | #include <dumux/freeflow/navierstokes/localresidual.hh> | ||
19 | |||
20 | namespace Dumux { | ||
21 | |||
22 | /*! | ||
23 | * \ingroup LowReKEpsilonModel | ||
24 | * \brief Element-wise calculation of the residual for low-Reynolds k-epsilon models using the staggered discretization | ||
25 | */ | ||
26 | |||
27 | // forward declaration | ||
28 | template<class TypeTag, class BaseLocalResidual, class DiscretizationMethod> | ||
29 | class LowReKEpsilonResidualImpl; | ||
30 | |||
31 | template<class TypeTag, class BaseLocalResidual> | ||
32 | class LowReKEpsilonResidualImpl<TypeTag, BaseLocalResidual, DiscretizationMethods::Staggered> | ||
33 | : public BaseLocalResidual | ||
34 | { | ||
35 | using ParentType = BaseLocalResidual; | ||
36 | |||
37 | using GridVariables = GetPropType<TypeTag, Properties::GridVariables>; | ||
38 | |||
39 | using GridVolumeVariables = typename GridVariables::GridVolumeVariables; | ||
40 | using ElementVolumeVariables = typename GridVolumeVariables::LocalView; | ||
41 | using VolumeVariables = typename GridVolumeVariables::VolumeVariables; | ||
42 | |||
43 | using GridFluxVariablesCache = typename GridVariables::GridFluxVariablesCache; | ||
44 | using ElementFluxVariablesCache = typename GridFluxVariablesCache::LocalView; | ||
45 | using FluxVariables = GetPropType<TypeTag, Properties::FluxVariables>; | ||
46 | |||
47 | using GridFaceVariables = typename GridVariables::GridFaceVariables; | ||
48 | using ElementFaceVariables = typename GridFaceVariables::LocalView; | ||
49 | |||
50 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
51 | using Problem = GetPropType<TypeTag, Properties::Problem>; | ||
52 | using GridView = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView; | ||
53 | using Element = typename GridView::template Codim<0>::Entity; | ||
54 | using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView; | ||
55 | using SubControlVolume = typename FVElementGeometry::SubControlVolume; | ||
56 | using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>; | ||
57 | using Indices = typename GetPropType<TypeTag, Properties::ModelTraits>::Indices; | ||
58 | using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>; | ||
59 | |||
60 | static constexpr int turbulentKineticEnergyEqIdx = Indices::turbulentKineticEnergyEqIdx - ModelTraits::dim(); | ||
61 | static constexpr int dissipationEqIdx = Indices::dissipationEqIdx - ModelTraits::dim(); | ||
62 | |||
63 | public: | ||
64 |
4/8✓ Branch 1 taken 115432 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 115432 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 115432 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 115432 times.
✗ Branch 11 not taken.
|
461728 | using ParentType::ParentType; |
65 | |||
66 | //! Evaluate fluxes entering or leaving the cell center control volume. | ||
67 | ✗ | CellCenterPrimaryVariables computeStorageForCellCenter(const Problem& problem, | |
68 | const SubControlVolume& scv, | ||
69 | const VolumeVariables& volVars) const | ||
70 | { | ||
71 | 7332800 | CellCenterPrimaryVariables storage = ParentType::computeStorageForCellCenter(problem, scv, volVars); | |
72 | |||
73 | 14665600 | storage[turbulentKineticEnergyEqIdx] = volVars.turbulentKineticEnergy() * volVars.density(); | |
74 | 21998400 | storage[dissipationEqIdx] = volVars.dissipationTilde() * volVars.density(); | |
75 | |||
76 | ✗ | return storage; | |
77 | } | ||
78 | |||
79 | 5066880 | CellCenterPrimaryVariables computeSourceForCellCenter(const Problem& problem, | |
80 | const Element &element, | ||
81 | const FVElementGeometry& fvGeometry, | ||
82 | const ElementVolumeVariables& elemVolVars, | ||
83 | const ElementFaceVariables& elemFaceVars, | ||
84 | const SubControlVolume &scv) const | ||
85 | { | ||
86 | 5066880 | CellCenterPrimaryVariables source = ParentType::computeSourceForCellCenter(problem, element, fvGeometry, | |
87 | elemVolVars, elemFaceVars, scv); | ||
88 | |||
89 | 5066880 | const auto& volVars = elemVolVars[scv]; | |
90 | |||
91 | // production | ||
92 | 15200640 | source[turbulentKineticEnergyEqIdx] += 2.0 * volVars.dynamicEddyViscosity() | |
93 | 5066880 | * volVars.stressTensorScalarProduct(); | |
94 | 15200640 | source[dissipationEqIdx] += volVars.cOneEpsilon() * volVars.fOne() | |
95 | 5066880 | * volVars.dissipationTilde() / volVars.turbulentKineticEnergy() | |
96 | 5066880 | * 2.0 * volVars.dynamicEddyViscosity() | |
97 | 5066880 | * volVars.stressTensorScalarProduct(); | |
98 | |||
99 | // destruction | ||
100 | 10133760 | source[turbulentKineticEnergyEqIdx] -= volVars.dissipationTilde() * volVars.density(); | |
101 | 5066880 | source[dissipationEqIdx] -= volVars.cTwoEpsilon() * volVars.fTwo() * volVars.density() | |
102 | 5066880 | * volVars.dissipationTilde() * volVars.dissipationTilde() | |
103 | 5066880 | / volVars.turbulentKineticEnergy(); | |
104 | |||
105 | // dampening functions | ||
106 | 15200640 | source[turbulentKineticEnergyEqIdx] -= volVars.dValue() * volVars.density(); | |
107 | 5066880 | source[dissipationEqIdx] += volVars.eValue() * volVars.density(); | |
108 | |||
109 | 5066880 | return source; | |
110 | } | ||
111 | }; | ||
112 | } // end namespace Dumux | ||
113 | |||
114 | #endif | ||
115 |