25#ifndef DUMUX_FVVELOCITY1P_HH
26#define DUMUX_FVVELOCITY1P_HH
46template<
class TypeTag>
59 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
62 using Element =
typename GridView::Traits::template Codim<0>::Entity;
63 using Intersection =
typename GridView::Intersection;
67 dim = GridView::dimension, dimWorld = GridView::dimensionworld
72 pressEqIdx = Indices::pressureEqIdx
75 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
76 using VelocityVector = Dune::FieldVector<Scalar, dimWorld>;
77 using GravityVector = Dune::FieldVector<Scalar, dimWorld>;
78 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
87 : problem_(problem), gravity_(problem.gravity())
89 const auto element = *problem_.gridView().template begin<0>();
91 Scalar referencePress = problem_.referencePressure(element);
111 template<
class MultiWriter>
114 Dune::BlockVector<Dune::FieldVector<Scalar, dim> > &velocity = *(writer.template allocateManagedBuffer<Scalar, dim> (
115 problem_.gridView().size(0)));
118 for (
const auto& element : elements(problem_.gridView()))
121 int eIdxGlobal = problem_.variables().index(element);
123 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
125 const typename Element::Geometry& geometry = element.geometry();
128 const auto refElement = referenceElement(geometry);
130 const int numberOfFaces = refElement.size(1);
131 std::vector<Scalar> flux(numberOfFaces,0);
134 for (
const auto& intersection : intersections(problem_.gridView(), element))
136 int isIndex = intersection.indexInInside();
138 flux[isIndex] = intersection.geometry().volume()
139 * (intersection.centerUnitOuterNormal() * cellData.fluxData().velocity(isIndex));
144 Dune::FieldVector<Scalar, dim> refVelocity;
146 if (refElement.type().isSimplex()) {
147 for (
int dimIdx = 0; dimIdx < dim; dimIdx++)
149 refVelocity[dimIdx] = -flux[dim - 1 - dimIdx];
150 for (
int fIdx = 0; fIdx < dim + 1; fIdx++)
152 refVelocity[dimIdx] += flux[fIdx]/(dim + 1);
157 else if (refElement.type().isCube()){
158 for (
int i = 0; i < dim; i++)
159 refVelocity[i] = 0.5 * (flux[2*i + 1] - flux[2*i]);
163 DUNE_THROW(Dune::NotImplemented,
"velocity output for prism/pyramid not implemented");
166 const Dune::FieldVector<Scalar, dim>& localPos
167 = refElement.position(0, 0);
170 const typename Element::Geometry::JacobianTransposed& jacobianT =
171 geometry.jacobianTransposed(localPos);
174 Dune::FieldVector<Scalar, dim> elementVelocity(0);
175 jacobianT.umtv(refVelocity, elementVelocity);
176 elementVelocity /= geometry.integrationElement(localPos);
178 velocity[eIdxGlobal] = elementVelocity;
181 writer.attachCellData(velocity,
"velocity", dim);
187 const GravityVector& gravity_;
201template<
class TypeTag>
204 auto elementI = intersection.inside();
205 auto elementJ = intersection.outside();
207 int eIdxGlobalJ = problem_.variables().index(elementJ);
209 CellData& cellDataJ = problem_.variables().cellData(eIdxGlobalJ);
212 const GlobalPosition& globalPosI = elementI.geometry().center();
213 const GlobalPosition& globalPosJ = elementJ.geometry().center();
216 int isIndexI = intersection.indexInInside();
217 int isIndexJ = intersection.indexInOutside();
220 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
223 GlobalPosition distVec = globalPosJ - globalPosI;
226 Scalar dist = distVec.two_norm();
229 DimMatrix meanPermeability(0);
231 problem_.spatialParams().meanK(meanPermeability, problem_.spatialParams().intrinsicPermeability(elementI),
232 problem_.spatialParams().intrinsicPermeability(elementJ));
240 Scalar potential = (cellData.pressure() - cellDataJ.pressure()) / dist;
242 potential += density_ * (unitOuterNormal * gravity_);
245 cellData.fluxData().setPotential(isIndexI, potential);
246 cellDataJ.fluxData().setPotential(isIndexJ, -potential);
250 velocity *= (cellData.pressure() - cellDataJ.pressure()) / dist;
252 GravityVector gravityTerm(unitOuterNormal);
255 velocity += gravityTerm;
258 cellData.fluxData().setVelocity(isIndexI, velocity);
259 cellData.fluxData().setVelocityMarker(isIndexI);
261 cellDataJ.fluxData().setVelocity(isIndexJ, velocity);
262 cellDataJ.fluxData().setVelocityMarker(isIndexJ);
275template<
class TypeTag>
278 auto element = intersection.inside();
281 int isIndex = intersection.indexInInside();
284 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
286 BoundaryTypes bcType;
288 problem_.boundaryTypes(bcType, intersection);
289 PrimaryVariables boundValues(0.0);
291 if (bcType.isDirichlet(pressEqIdx))
293 problem_.dirichlet(boundValues, intersection);
296 const GlobalPosition& globalPosI = element.geometry().center();
299 const GlobalPosition& globalPosJ = intersection.geometry().center();
302 GlobalPosition distVec = globalPosJ - globalPosI;
305 Scalar dist = distVec.two_norm();
309 DimMatrix meanPermeability(0);
311 problem_.spatialParams().meanK(meanPermeability, problem_.spatialParams().intrinsicPermeability(element));
318 Scalar pressBound = boundValues;
321 Scalar potential = (cellData.pressure() - pressBound) / dist;
323 potential += density_ * (unitOuterNormal * gravity_);
326 cellData.fluxData().setPotential(isIndex, potential);
330 velocity *= (cellData.pressure() - pressBound) / dist;
332 GravityVector gravityTerm(unitOuterNormal);
335 velocity += gravityTerm;
338 cellData.fluxData().setVelocity(isIndex, velocity);
339 cellData.fluxData().setVelocityMarker(isIndex);
345 problem_.neumann(boundValues, intersection);
346 VelocityVector velocity(unitOuterNormal);
348 velocity *= boundValues[pressEqIdx] / density_;
351 cellData.fluxData().setPotential(isIndex, boundValues[pressEqIdx]);
354 cellData.fluxData().setVelocity(isIndex, velocity);
355 cellData.fluxData().setVelocityMarker(isIndex);
void calculateVelocityOnBoundary(const Intersection &, CellData &)
Calculates the velocity at a boundary.
Definition: 1p/sequential/diffusion/cellcentered/velocity.hh:276
void calculateVelocity(const Intersection &, CellData &)
Calculates the velocity at a cell-cell interface.
Definition: 1p/sequential/diffusion/cellcentered/velocity.hh:202
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type GetProp
get the type of a property (equivalent to old macro GET_PROP(...))
Definition: propertysystem.hh:140
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
std::string temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:51
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:74
std::string permeability() noexcept
I/O name of permeability.
Definition: name.hh:143
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
Single phase finite volume velocity reconstruction.
Definition: 1p/sequential/diffusion/cellcentered/velocity.hh:48
FVVelocity1P(Problem &problem)
Constructs a FVVelocity1P object.
Definition: 1p/sequential/diffusion/cellcentered/velocity.hh:86
void addOutputVtkFields(MultiWriter &writer)
Adds velocity output to the output file.
Definition: 1p/sequential/diffusion/cellcentered/velocity.hh:112
Defines the properties required for the single phase sequential model.