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-FileCopyrightText: 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 Elastic | ||
10 | * \brief Element-wise calculation of the local residual for problems | ||
11 | * using the elastic model considering linear elasticity. | ||
12 | */ | ||
13 | #ifndef DUMUX_SOLIDMECHANICS_ELASTIC_LOCAL_RESIDUAL_HH | ||
14 | #define DUMUX_SOLIDMECHANICS_ELASTIC_LOCAL_RESIDUAL_HH | ||
15 | |||
16 | #include <dumux/common/parameters.hh> | ||
17 | #include <dumux/common/properties.hh> | ||
18 | #include <dumux/common/parameters.hh> | ||
19 | #include <dumux/common/numeqvector.hh> | ||
20 | #include <dumux/discretization/defaultlocaloperator.hh> | ||
21 | |||
22 | namespace Dumux { | ||
23 | |||
24 | /*! | ||
25 | * \ingroup Elastic | ||
26 | * \brief Element-wise calculation of the local residual for problems | ||
27 | * using the elastic model considering linear elasticity. | ||
28 | */ | ||
29 | template<class TypeTag> | ||
30 | class ElasticLocalResidual | ||
31 | : public DiscretizationDefaultLocalOperator<TypeTag> | ||
32 | { | ||
33 | using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>; | ||
34 | using ParentType = DiscretizationDefaultLocalOperator<TypeTag>; | ||
35 | |||
36 | using GridView = typename GridGeometry::GridView; | ||
37 | using Element = typename GridView::template Codim<0>::Entity; | ||
38 | |||
39 | using Problem = GetPropType<TypeTag, Properties::Problem>; | ||
40 | using Indices = typename GetPropType<TypeTag, Properties::ModelTraits>::Indices; | ||
41 | using FVElementGeometry = typename GridGeometry::LocalView; | ||
42 | using SubControlVolume = typename FVElementGeometry::SubControlVolume; | ||
43 | using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace; | ||
44 | using NumEqVector = Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>; | ||
45 | using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView; | ||
46 | using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView; | ||
47 | using VolumeVariables = typename ElementVolumeVariables::VolumeVariables; | ||
48 | |||
49 | // class assembling the stress tensor | ||
50 | using StressType = GetPropType<TypeTag, Properties::StressType>; | ||
51 | |||
52 | public: | ||
53 |
1/2✓ Branch 1 taken 5436 times.
✗ Branch 2 not taken.
|
5436 | using ParentType::ParentType; |
54 | |||
55 | /*! | ||
56 | * \brief For the elastic model the storage term is zero since | ||
57 | * we neglect inertia forces. | ||
58 | * | ||
59 | * \param problem The problem | ||
60 | * \param scv The sub control volume | ||
61 | * \param volVars The current or previous volVars | ||
62 | */ | ||
63 | 701440 | NumEqVector computeStorage(const Problem& problem, | |
64 | const SubControlVolume& scv, | ||
65 | const VolumeVariables& volVars) const | ||
66 | { | ||
67 | 1402880 | return NumEqVector(0.0); | |
68 | } | ||
69 | |||
70 | /*! | ||
71 | * \brief Evaluates the force in all grid directions acting | ||
72 | * on a face of a sub-control volume. | ||
73 | * | ||
74 | * \param problem The problem | ||
75 | * \param element The current element. | ||
76 | * \param fvGeometry The finite-volume geometry | ||
77 | * \param elemVolVars The volume variables of the current element | ||
78 | * \param scvf The sub control volume face to compute the flux on | ||
79 | * \param elemFluxVarsCache The cache related to flux computation | ||
80 | */ | ||
81 | 1204176 | NumEqVector computeFlux(const Problem& problem, | |
82 | const Element& element, | ||
83 | const FVElementGeometry& fvGeometry, | ||
84 | const ElementVolumeVariables& elemVolVars, | ||
85 | const SubControlVolumeFace& scvf, | ||
86 | const ElementFluxVariablesCache& elemFluxVarsCache) const | ||
87 | { | ||
88 | // obtain force on the face from stress type | ||
89 | 1204176 | return StressType::force(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache); | |
90 | } | ||
91 | |||
92 | /*! | ||
93 | * \brief Calculate the source term of the equation | ||
94 | * | ||
95 | * \param problem The problem to solve | ||
96 | * \param element The DUNE Codim<0> entity for which the residual | ||
97 | * ought to be calculated | ||
98 | * \param fvGeometry The finite-volume geometry of the element | ||
99 | * \param elemVolVars The volume variables associated with the element stencil | ||
100 | * \param scv The sub-control volume over which we integrate the source term | ||
101 | * \note This is the default implementation for geomechanical models adding to | ||
102 | * the user defined sources the source stemming from the gravitational acceleration. | ||
103 | * | ||
104 | */ | ||
105 | 7200 | NumEqVector computeSource(const Problem& problem, | |
106 | const Element& element, | ||
107 | const FVElementGeometry& fvGeometry, | ||
108 | const ElementVolumeVariables& elemVolVars, | ||
109 | const SubControlVolume &scv) const | ||
110 | { | ||
111 | 7200 | NumEqVector source(0.0); | |
112 | |||
113 | // add contributions from volume flux sources | ||
114 | 7200 | source += problem.source(element, fvGeometry, elemVolVars, scv); | |
115 | |||
116 | // add contribution from possible point sources | ||
117 | 7200 | source += problem.scvPointSources(element, fvGeometry, elemVolVars, scv); | |
118 | |||
119 | // maybe add gravitational acceleration | ||
120 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 7199 times.
✓ Branch 3 taken 1 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
|
7200 | static const bool gravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity"); |
121 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 7200 times.
|
7200 | if (gravity) |
122 | { | ||
123 | ✗ | const auto& g = problem.spatialParams().gravity(scv.center()); | |
124 | ✗ | for (int dir = 0; dir < GridView::dimensionworld; ++dir) | |
125 | ✗ | source[Indices::momentum(dir)] += elemVolVars[scv].solidDensity()*g[dir]; | |
126 | } | ||
127 | |||
128 | 7200 | return source; | |
129 | } | ||
130 | }; | ||
131 | |||
132 | } // end namespace Dumux | ||
133 | |||
134 | #endif | ||
135 |