version 3.8
geomechanics/hyperelastic/localresidual.hh
Go to the documentation of this file.
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//
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>
22
23namespace Dumux {
24
29template<class TypeTag>
31: public GetPropType<TypeTag, Properties::BaseLocalResidual>
32{
38 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
39 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
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>;
48 using Indices = typename ModelTraits::Indices;
49 static constexpr int dimWorld = GridView::dimensionworld;
50public:
51 using ParentType::ParentType;
52 using ElementResidualVector = typename ParentType::ElementResidualVector;
53
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
69 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 const auto F = [&]()
78 {
79 const auto& fluxVarCache = elemFluxVarsCache[scvf];
80 Dune::FieldMatrix<Scalar, dimWorld, dimWorld> F(0.0);
81 for (const auto& scv : scvs(fvGeometry))
82 {
83 const auto& volVars = elemVolVars[scv];
84 for (int dir = 0; dir < dimWorld; ++dir)
85 F[dir].axpy(volVars.displacement(dir), fluxVarCache.gradN(scv.indexInElement()));
86 }
87
88 for (int dir = 0; dir < dimWorld; ++dir)
89 F[dir][dir] += 1;
90
91 return F;
92 }();
93
94 // numerical flux ==> -P*n*dA
95 NumEqVector flux(0.0);
96
97 // problem implements the constitutive law P=P(F)
98 const auto P = problem.firstPiolaKirchhoffStressTensor(F);
99 P.mv(scvf.unitOuterNormal(), flux);
100 flux *= -scvf.area();
101 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 return ParentType::computeSource(problem, element, fvGeometry, elemVolVars, scv);
112 }
113};
114
115} // end namespace Dumux
116
117#endif
Local residual for the hyperelastic model.
Definition: geomechanics/hyperelastic/localresidual.hh:32
typename ParentType::ElementResidualVector ElementResidualVector
Definition: geomechanics/hyperelastic/localresidual.hh:52
NumEqVector computeStorage(const Problem &problem, const SubControlVolume &scv, const VolumeVariables &volVars) const
Evaluate the rate of change of all conserved quantities.
Definition: geomechanics/hyperelastic/localresidual.hh:57
NumEqVector computeSource(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv) const
Definition: geomechanics/hyperelastic/localresidual.hh:104
NumEqVector computeFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache) const
Evaluate the fluxes over a face of a sub control volume flux: -div(P) where P is the 1st Piola-Kircho...
Definition: geomechanics/hyperelastic/localresidual.hh:69
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
typename NumEqVectorTraits< PrimaryVariables >::type NumEqVector
A vector with the same size as numbers of equations This is the default implementation and has to be ...
Definition: numeqvector.hh:34
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Define some often used mathematical functions.
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
A helper to deduce a vector with the same size as numbers of equations.