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 Hyperelastic | ||
10 | * \brief Local residual for the hyperelastic model | ||
11 | */ | ||
12 | #ifndef DUMUX_SOLIDMECHANICS_HYPERELASTIC_LOCAL_RESIDUAL_HH | ||
13 | #define DUMUX_SOLIDMECHANICS_HYPERELASTIC_LOCAL_RESIDUAL_HH | ||
14 | |||
15 | #include <dune/common/fvector.hh> | ||
16 | #include <dune/common/fmatrix.hh> | ||
17 | #include <dune/common/std/type_traits.hh> | ||
18 | #include <dune/common/exceptions.hh> | ||
19 | |||
20 | #include <dumux/common/math.hh> | ||
21 | #include <dumux/common/properties.hh> | ||
22 | #include <dumux/common/numeqvector.hh> | ||
23 | #include <dumux/discretization/extrusion.hh> | ||
24 | #include <dumux/discretization/defaultlocaloperator.hh> | ||
25 | |||
26 | namespace Dumux { | ||
27 | |||
28 | #ifndef DOXYGEN | ||
29 | namespace Detail::Hyperelastic { | ||
30 | // helper struct detecting if the user-defined problem class has a solidDensity method | ||
31 | template <typename T, typename ...Ts> | ||
32 | using SolidDensityDetector = decltype(std::declval<T>().solidDensity(std::declval<Ts>()...)); | ||
33 | |||
34 | template<class T, typename ...Args> | ||
35 | static constexpr bool hasSolidDensity() | ||
36 | { return Dune::Std::is_detected<SolidDensityDetector, T, Args...>::value; } | ||
37 | |||
38 | // helper struct detecting if the user-defined problem class has a acceleration method | ||
39 | template <typename T, typename ...Ts> | ||
40 | using AccelerationDetector = decltype(std::declval<T>().acceleration(std::declval<Ts>()...)); | ||
41 | |||
42 | template<class T, typename ...Args> | ||
43 | static constexpr bool hasAcceleration() | ||
44 | { return Dune::Std::is_detected<AccelerationDetector, T, Args...>::value; } | ||
45 | } // end namespace Detail::Hyperelastic | ||
46 | #endif | ||
47 | |||
48 | /*! | ||
49 | * \ingroup Hyperelastic | ||
50 | * \brief Local residual for the hyperelastic model | ||
51 | */ | ||
52 | template<class TypeTag> | ||
53 | class HyperelasticLocalResidual | ||
54 | : public DiscretizationDefaultLocalOperator<TypeTag> | ||
55 | { | ||
56 | using ParentType = DiscretizationDefaultLocalOperator<TypeTag>; | ||
57 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
58 | using Problem = GetPropType<TypeTag, Properties::Problem>; | ||
59 | using NumEqVector = Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>; | ||
60 | using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>; | ||
61 | using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView; | ||
62 | using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView; | ||
63 | using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>; | ||
64 | using FVElementGeometry = typename GridGeometry::LocalView; | ||
65 | using SubControlVolume = typename GridGeometry::SubControlVolume; | ||
66 | using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace; | ||
67 | using GridView = typename GridGeometry::GridView; | ||
68 | using Element = typename GridView::template Codim<0>::Entity; | ||
69 | using Extrusion = Extrusion_t<GridGeometry>; | ||
70 | using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>; | ||
71 | using Indices = typename ModelTraits::Indices; | ||
72 | static constexpr int dimWorld = GridView::dimensionworld; | ||
73 | public: | ||
74 | 595284 | using ParentType::ParentType; | |
75 | using ElementResidualVector = typename ParentType::ElementResidualVector; | ||
76 | |||
77 | using ParentType::evalStorage; | ||
78 | |||
79 | /*! | ||
80 | * \brief Evaluate the storage contribution to the residual: rho*d^2u/dt^2 | ||
81 | */ | ||
82 | 11468772 | void evalStorage(ElementResidualVector& residual, | |
83 | const Problem& problem, | ||
84 | const Element& element, | ||
85 | const FVElementGeometry& fvGeometry, | ||
86 | const ElementVolumeVariables& prevElemVolVars, | ||
87 | const ElementVolumeVariables& curElemVolVars, | ||
88 | const SubControlVolume& scv) const | ||
89 | { | ||
90 | 11468772 | const auto& curVolVars = curElemVolVars[scv]; | |
91 | |||
92 | if constexpr (Detail::Hyperelastic::hasSolidDensity<Problem, const Element&, const SubControlVolume&>() | ||
93 | && Detail::Hyperelastic::hasAcceleration<Problem, const Element&, const SubControlVolume&, Scalar, decltype(curVolVars.displacement())>()) | ||
94 | { | ||
95 | 11468772 | const auto& curVolVars = curElemVolVars[scv]; | |
96 | 11468772 | NumEqVector storage = problem.solidDensity(element, scv)*problem.acceleration( | |
97 | 11468772 | element, scv, this->timeLoop().timeStepSize(), curVolVars.displacement() | |
98 | ); | ||
99 | 11468772 | storage *= curVolVars.extrusionFactor()*Extrusion::volume(fvGeometry, scv); | |
100 | 11468772 | residual[scv.localDofIndex()] += storage; | |
101 | } | ||
102 | else | ||
103 | ✗ | DUNE_THROW(Dune::InvalidStateException, | |
104 | "Problem class must implement solidDensity and acceleration" | ||
105 | " methods to evaluate storage term" | ||
106 | ); | ||
107 | 11468772 | } | |
108 | |||
109 | /*! | ||
110 | * \brief Evaluate the fluxes over a face of a sub control volume | ||
111 | * flux: -div(P) where P is the 1st Piola-Kirchoff stress tensor (stress in reference configuration) | ||
112 | */ | ||
113 | 26214372 | NumEqVector computeFlux(const Problem& problem, | |
114 | const Element& element, | ||
115 | const FVElementGeometry& fvGeometry, | ||
116 | const ElementVolumeVariables& elemVolVars, | ||
117 | const SubControlVolumeFace& scvf, | ||
118 | const ElementFluxVariablesCache& elemFluxVarsCache) const | ||
119 | { | ||
120 | // the deformation gradient F = grad(d) + I | ||
121 | 52428744 | const auto F = [&]() | |
122 | { | ||
123 | 26214372 | const auto& fluxVarCache = elemFluxVarsCache[scvf]; | |
124 | 26214372 | Dune::FieldMatrix<Scalar, dimWorld, dimWorld> F(0.0); | |
125 |
2/2✓ Branch 0 taken 152371116 times.
✓ Branch 1 taken 26214372 times.
|
178585488 | for (const auto& scv : scvs(fvGeometry)) |
126 | { | ||
127 | 152371116 | const auto& volVars = elemVolVars[scv]; | |
128 |
2/2✓ Branch 0 taken 422707032 times.
✓ Branch 1 taken 152371116 times.
|
575078148 | for (int dir = 0; dir < dimWorld; ++dir) |
129 | 845414064 | F[dir].axpy(volVars.displacement(dir), fluxVarCache.gradN(scv.indexInElement())); | |
130 | } | ||
131 | |||
132 |
2/2✓ Branch 0 taken 67174344 times.
✓ Branch 1 taken 26214372 times.
|
93388716 | for (int dir = 0; dir < dimWorld; ++dir) |
133 | 67174344 | F[dir][dir] += 1; | |
134 | |||
135 | 26214372 | return F; | |
136 | 26214372 | }(); | |
137 | |||
138 | // numerical flux ==> -P*n*dA | ||
139 | 26214372 | NumEqVector flux(0.0); | |
140 | |||
141 | // problem implements the constitutive law P=P(F) | ||
142 | 26214372 | const auto P = problem.firstPiolaKirchhoffStressTensor(F); | |
143 | 52428744 | P.mv(scvf.unitOuterNormal(), flux); | |
144 | 26214372 | flux *= -scvf.area(); | |
145 | 26214372 | return flux; | |
146 | } | ||
147 | |||
148 | 21299172 | NumEqVector computeSource(const Problem& problem, | |
149 | const Element& element, | ||
150 | const FVElementGeometry& fvGeometry, | ||
151 | const ElementVolumeVariables& elemVolVars, | ||
152 | const SubControlVolume& scv) const | ||
153 | { | ||
154 | // loading and body forces are implemented in the problem | ||
155 | 21299172 | return ParentType::computeSource(problem, element, fvGeometry, elemVolVars, scv); | |
156 | } | ||
157 | }; | ||
158 | |||
159 | } // end namespace Dumux | ||
160 | |||
161 | #endif | ||
162 |