version 3.11-dev
Loading...
Searching...
No Matches
felocalresidual.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-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_NAVIERSTOKES_MOMENTUM_CVFE_FE_LOCAL_RESIDUAL_HELPER_HH
13#define DUMUX_NAVIERSTOKES_MOMENTUM_CVFE_FE_LOCAL_RESIDUAL_HELPER_HH
14
15#include <dune/geometry/quadraturerules.hh>
16#include <dumux/common/typetraits/localdofs_.hh>
20
21#include "flux.hh"
22
23namespace Dumux {
24
29template<class Scalar, class NumEqVector, class LocalBasis, class Extrusion>
31{
32 using RangeType = typename LocalBasis::Traits::RangeType;
33
34public:
45 template<class ResidualVector, class Problem, class FVElementGeometry, class ElementVariables>
46 static void addStorageTerms(ResidualVector& residual,
47 const Problem& problem,
48 const FVElementGeometry& fvGeometry,
49 const ElementVariables& prevElemVars,
50 const ElementVariables& curElemVars,
51 const Scalar timeStepSize)
52 {
53 if constexpr (Detail::LocalDofs::hasNonCVLocalDofsInterface<FVElementGeometry>())
54 {
55 // Make sure we don't iterate over quadrature points if there are no hybrid dofs
56 if (nonCVLocalDofs(fvGeometry).empty())
57 return;
58
59 const auto& localBasis = fvGeometry.feLocalBasis();
60 std::vector<RangeType> integralShapeFunctions(localBasis.size(), RangeType(0.0));
61
62 // We apply mass lumping such that we only need to calculate the integral of basis functions
63 const auto& geometry = fvGeometry.elementGeometry();
64 const auto& element = fvGeometry.element();
65 using GlobalPosition = typename FVElementGeometry::GridGeometry::GlobalCoordinate;
67
68 for (const auto& qpData : CVFE::quadratureRule(fvGeometry, element))
69 {
70 const auto& ipData = qpData.ipData();
71 // Obtain and store shape function values and gradients at the current quad point
72 FeIpData feIpData(geometry, ipData.local(), ipData.global(), localBasis);
73
74 // get density from the problem
75 for (const auto& localDof : nonCVLocalDofs(fvGeometry))
76 integralShapeFunctions[localDof.index()] += qpData.weight() * feIpData.shapeValue(localDof.index());
77 }
78
79 for (const auto& localDof : nonCVLocalDofs(fvGeometry))
80 {
81 const auto localDofIdx = localDof.index();
82 const auto& data = ipData(fvGeometry, localDof);
83 const auto curDensity = problem.density(element, fvGeometry, data, false);
84 const auto prevDensity = problem.density(element, fvGeometry, data, true);
85 const auto curVelocity = curElemVars[localDofIdx].velocity();
86 const auto prevVelocity = prevElemVars[localDofIdx].velocity();
87 auto timeDeriv = (curDensity*curVelocity - prevDensity*prevVelocity);
88 timeDeriv /= timeStepSize;
89
90 // add storage to residual
91 for (int eqIdx = 0; eqIdx < NumEqVector::dimension; ++eqIdx)
92 residual[localDofIdx][eqIdx] += integralShapeFunctions[localDofIdx]*timeDeriv[eqIdx];
93 }
94 }
95 }
96
105 template<class ResidualVector, class Problem, class FVElementGeometry, class ElementVariables>
106 static void addFluxAndSourceTerms(ResidualVector& residual,
107 const Problem& problem,
108 const FVElementGeometry& fvGeometry,
109 const ElementVariables& elemVars)
110 {
111 if constexpr (Detail::LocalDofs::hasNonCVLocalDofsInterface<FVElementGeometry>())
112 {
113 // Make sure we don't iterate over quadrature points if there are no hybrid dofs
114 if (nonCVLocalDofs(fvGeometry).empty())
115 return;
116
117 if constexpr (requires { problem.pointSources(); })
118 {
119 if (!problem.pointSources().empty())
120 DUNE_THROW(Dune::NotImplemented, "Point sources are not implemented for hybrid momentum schemes.");
121 }
122
123 static const bool enableUnsymmetrizedVelocityGradient
124 = getParamFromGroup<bool>(problem.paramGroup(), "FreeFlow.EnableUnsymmetrizedVelocityGradient", false);
125
126 const auto& element = fvGeometry.element();
127 using Cache = typename ElementVariables::InterpolationPointData;
129 for (const auto& qpData : CVFE::quadratureRule(fvGeometry, element))
130 {
131 const auto& ipData = qpData.ipData();
132 // Obtain and store shape function values and gradients at the current quad point
133 const auto& ipCache = cache(elemVars, ipData);
134 FluxFunctionContext context(problem, fvGeometry, elemVars, ipCache);
135 const auto& v = context.velocity();
136 const auto& gradV = context.gradVelocity();
137
138 // get viscosity from the problem
139 const Scalar mu = problem.effectiveViscosity(element, fvGeometry, ipData);
140 // get density from the problem
141 const Scalar density = problem.density(element, fvGeometry, ipData);
142
143 for (const auto& localDof : nonCVLocalDofs(fvGeometry))
144 {
145 const auto localDofIdx = localDof.index();
146 NumEqVector fluxAndSourceTerm(0.0);
147 // add advection term
148 if (problem.enableInertiaTerms())
149 fluxAndSourceTerm -= density*(v*ipCache.gradN(localDofIdx))*v;
150
151 // add diffusion term
152 fluxAndSourceTerm += enableUnsymmetrizedVelocityGradient ?
153 mu*mv(gradV, ipCache.gradN(localDofIdx))
154 : mu*mv(gradV + getTransposed(gradV), ipCache.gradN(localDofIdx));
155
156 // add pressure term
157 fluxAndSourceTerm -= problem.pressure(element, fvGeometry, ipData) * ipCache.gradN(localDofIdx);
158
159 // finally add source and flux boundary term and add everything to residual
160 auto sourceAtIp = problem.source(fvGeometry, elemVars, ipData);
161 // add gravity term rho*g (note that gravity might be zero in case it's disabled in the problem)
162 sourceAtIp += density * problem.gravity();
163
164 const auto& shapeValues = ipCache.shapeValues();
165 for (int eqIdx = 0; eqIdx < NumEqVector::dimension; ++eqIdx)
166 {
167 fluxAndSourceTerm[eqIdx] -= shapeValues[localDofIdx] * sourceAtIp[eqIdx];
168 residual[localDofIdx][eqIdx] += qpData.weight()*fluxAndSourceTerm[eqIdx];
169 }
170 }
171 }
172 }
173 }
174
175};
176
177} // end namespace Dumux
178
179#endif
Boundary flag to store e.g. in sub control volume faces.
Definition fem/interpolationpointdata.hh:21
Helper class for evaluating FE-based local residuals.
Definition felocalresidual.hh:31
static void addFluxAndSourceTerms(ResidualVector &residual, const Problem &problem, const FVElementGeometry &fvGeometry, const ElementVariables &elemVars)
Add flux and source residual contribution for non-CV local dofs.
Definition felocalresidual.hh:106
static void addStorageTerms(ResidualVector &residual, const Problem &problem, const FVElementGeometry &fvGeometry, const ElementVariables &prevElemVars, const ElementVariables &curElemVars, const Scalar timeStepSize)
Add storage residual contribution for non-CV local dofs.
Definition felocalresidual.hh:46
Context for interpolating data on interpolation points.
Definition flux.hh:99
Helper classes to compute the integration elements.
Shape functions and gradients at an interpolation point.
The flux variables class for the Navier-Stokes model using control-volume finite element schemes.
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
Dune::DenseVector< V >::derived_type mv(const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V > &v)
Returns the result of the projection of a vector v with a Matrix M.
Definition math.hh:829
Dune::FieldMatrix< Scalar, n, m > getTransposed(const Dune::FieldMatrix< Scalar, m, n > &M)
Transpose a FieldMatrix.
Definition math.hh:712
T getParamFromGroup(Args &&... args)
A free function to get a parameter from the parameter tree singleton with a model group.
Definition parameters.hh:149
auto quadratureRule(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolume &scv, QuadratureRules::MidpointQuadrature)
Midpoint quadrature for scv.
Definition quadraturerules.hh:159
Definition adapt.hh:17