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 Hyperelastic | ||
10 | * \brief Local residual for the hyperelastic model | ||
11 | */ | ||
12 | #ifndef DUMUX_GEOMECHANICS_HYPERELASTIC_LOCAL_RESIDUAL_HH | ||
13 | #define DUMUX_GEOMECHANICS_HYPERELASTIC_LOCAL_RESIDUAL_HH | ||
14 | |||
15 | #include <dune/common/fvector.hh> | ||
16 | #include <dune/common/fmatrix.hh> | ||
17 | |||
18 | #include <dumux/common/math.hh> | ||
19 | #include <dumux/common/properties.hh> | ||
20 | #include <dumux/common/numeqvector.hh> | ||
21 | #include <dumux/discretization/extrusion.hh> | ||
22 | |||
23 | namespace Dumux { | ||
24 | |||
25 | /*! | ||
26 | * \ingroup Hyperelastic | ||
27 | * \brief Local residual for the hyperelastic model | ||
28 | */ | ||
29 | template<class TypeTag> | ||
30 | class HyperelasticLocalResidual | ||
31 | : public GetPropType<TypeTag, Properties::BaseLocalResidual> | ||
32 | { | ||
33 | using ParentType = GetPropType<TypeTag, Properties::BaseLocalResidual>; | ||
34 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
35 | using Problem = GetPropType<TypeTag, Properties::Problem>; | ||
36 | using NumEqVector = Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>; | ||
37 | using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>; | ||
38 | using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView; | ||
39 | using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView; | ||
40 | using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>; | ||
41 | using FVElementGeometry = typename GridGeometry::LocalView; | ||
42 | using SubControlVolume = typename GridGeometry::SubControlVolume; | ||
43 | using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace; | ||
44 | using GridView = typename GridGeometry::GridView; | ||
45 | using Element = typename GridView::template Codim<0>::Entity; | ||
46 | using Extrusion = Extrusion_t<GridGeometry>; | ||
47 | using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>; | ||
48 | using Indices = typename ModelTraits::Indices; | ||
49 | static constexpr int dimWorld = GridView::dimensionworld; | ||
50 | public: | ||
51 | 98304 | using ParentType::ParentType; | |
52 | using ElementResidualVector = typename ParentType::ElementResidualVector; | ||
53 | |||
54 | /*! | ||
55 | * \brief Evaluate the rate of change of all conserved quantities | ||
56 | */ | ||
57 | ✗ | NumEqVector computeStorage(const Problem& problem, | |
58 | const SubControlVolume& scv, | ||
59 | const VolumeVariables& volVars) const | ||
60 | { | ||
61 | // only the static model is implemented so far | ||
62 | ✗ | return NumEqVector(0.0); | |
63 | } | ||
64 | |||
65 | /*! | ||
66 | * \brief Evaluate the fluxes over a face of a sub control volume | ||
67 | * flux: -div(P) where P is the 1st Piola-Kirchoff stress tensor (stress in reference configuration) | ||
68 | */ | ||
69 | 14745600 | NumEqVector computeFlux(const Problem& problem, | |
70 | const Element& element, | ||
71 | const FVElementGeometry& fvGeometry, | ||
72 | const ElementVolumeVariables& elemVolVars, | ||
73 | const SubControlVolumeFace& scvf, | ||
74 | const ElementFluxVariablesCache& elemFluxVarsCache) const | ||
75 | { | ||
76 | // the deformation gradient F = grad(d) + I | ||
77 | 14745600 | const auto F = [&]() | |
78 | { | ||
79 | 14745600 | const auto& fluxVarCache = elemFluxVarsCache[scvf]; | |
80 | 14745600 | Dune::FieldMatrix<Scalar, dimWorld, dimWorld> F(0.0); | |
81 |
4/4✓ Branch 0 taken 117964800 times.
✓ Branch 1 taken 14745600 times.
✓ Branch 2 taken 117964800 times.
✓ Branch 3 taken 14745600 times.
|
147456000 | for (const auto& scv : scvs(fvGeometry)) |
82 | { | ||
83 | 117964800 | const auto& volVars = elemVolVars[scv]; | |
84 |
2/2✓ Branch 0 taken 353894400 times.
✓ Branch 1 taken 117964800 times.
|
471859200 | for (int dir = 0; dir < dimWorld; ++dir) |
85 | 1769472000 | F[dir].axpy(volVars.displacement(dir), fluxVarCache.gradN(scv.indexInElement())); | |
86 | } | ||
87 | |||
88 |
2/2✓ Branch 0 taken 44236800 times.
✓ Branch 1 taken 14745600 times.
|
58982400 | for (int dir = 0; dir < dimWorld; ++dir) |
89 | 132710400 | F[dir][dir] += 1; | |
90 | |||
91 | 14745600 | return F; | |
92 | 14745600 | }(); | |
93 | |||
94 | // numerical flux ==> -P*n*dA | ||
95 | 14745600 | NumEqVector flux(0.0); | |
96 | |||
97 | // problem implements the constitutive law P=P(F) | ||
98 | 14745600 | const auto P = problem.firstPiolaKirchhoffStressTensor(F); | |
99 | 29491200 | P.mv(scvf.unitOuterNormal(), flux); | |
100 | 14745600 | flux *= -scvf.area(); | |
101 | 14745600 | return flux; | |
102 | } | ||
103 | |||
104 | NumEqVector computeSource(const Problem& problem, | ||
105 | const Element& element, | ||
106 | const FVElementGeometry& fvGeometry, | ||
107 | const ElementVolumeVariables& elemVolVars, | ||
108 | const SubControlVolume& scv) const | ||
109 | { | ||
110 | // loading and body forces are implemented in the problem | ||
111 | 9830400 | return ParentType::computeSource(problem, element, fvGeometry, elemVolVars, scv); | |
112 | } | ||
113 | }; | ||
114 | |||
115 | } // end namespace Dumux | ||
116 | |||
117 | #endif | ||
118 |