3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
1p/sequential/diffusion/cellcentered/velocity.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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_FVVELOCITY1P_HH
26#define DUMUX_FVVELOCITY1P_HH
27
29
30namespace Dumux {
46template<class TypeTag>
48{
52
54
56
59 using PrimaryVariables = typename SolutionTypes::PrimaryVariables;
61
62 using Element = typename GridView::Traits::template Codim<0>::Entity;
63 using Intersection = typename GridView::Intersection;
64
65 enum
66 {
67 dim = GridView::dimension, dimWorld = GridView::dimensionworld
68 };
69
70 enum
71 {
72 pressEqIdx = Indices::pressureEqIdx // only one equation!
73 };
74
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>;
79
80public:
86 FVVelocity1P(Problem& problem)
87 : problem_(problem), gravity_(problem.gravity())
88 {
89 const auto element = *problem_.gridView().template begin<0>();
90 Scalar temperature = problem_.temperature(element);
91 Scalar referencePress = problem_.referencePressure(element);
92
93 density_ = FluidSystem::density(temperature, referencePress);
94 viscosity_ = FluidSystem::viscosity(temperature, referencePress);
95 }
96
97 // Calculates the velocity at a cell-cell interface.
98 void calculateVelocity(const Intersection&, CellData&);
99
100 // Calculates the velocity at a boundary.
101 void calculateVelocityOnBoundary(const Intersection&, CellData&);
102
111 template<class MultiWriter>
112 void addOutputVtkFields(MultiWriter &writer)
113 {
114 Dune::BlockVector<Dune::FieldVector<Scalar, dim> > &velocity = *(writer.template allocateManagedBuffer<Scalar, dim> (
115 problem_.gridView().size(0)));
116
117 // compute update vector
118 for (const auto& element : elements(problem_.gridView()))
119 {
120 // cell index
121 int eIdxGlobal = problem_.variables().index(element);
122
123 CellData& cellData = problem_.variables().cellData(eIdxGlobal);
124
125 const typename Element::Geometry& geometry = element.geometry();
126
127 // get corresponding reference element
128 const auto refElement = referenceElement(geometry);
129
130 const int numberOfFaces = refElement.size(1);
131 std::vector<Scalar> flux(numberOfFaces,0);
132
133 // run through all intersections with neighbors and boundary
134 for (const auto& intersection : intersections(problem_.gridView(), element))
135 {
136 int isIndex = intersection.indexInInside();
137
138 flux[isIndex] = intersection.geometry().volume()
139 * (intersection.centerUnitOuterNormal() * cellData.fluxData().velocity(isIndex));
140 }
141
142 // calculate velocity on reference element as the Raviart-Thomas-0
143 // interpolant of the fluxes
144 Dune::FieldVector<Scalar, dim> refVelocity;
145 // simplices
146 if (refElement.type().isSimplex()) {
147 for (int dimIdx = 0; dimIdx < dim; dimIdx++)
148 {
149 refVelocity[dimIdx] = -flux[dim - 1 - dimIdx];
150 for (int fIdx = 0; fIdx < dim + 1; fIdx++)
151 {
152 refVelocity[dimIdx] += flux[fIdx]/(dim + 1);
153 }
154 }
155 }
156 // cubes
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]);
160 }
161 // 3D prism and pyramids
162 else {
163 DUNE_THROW(Dune::NotImplemented, "velocity output for prism/pyramid not implemented");
164 }
165
166 const Dune::FieldVector<Scalar, dim>& localPos
167 = refElement.position(0, 0);
168
169 // get the transposed Jacobian of the element mapping
170 const typename Element::Geometry::JacobianTransposed& jacobianT =
171 geometry.jacobianTransposed(localPos);
172
173 // calculate the element velocity by the Piola transformation
174 Dune::FieldVector<Scalar, dim> elementVelocity(0);
175 jacobianT.umtv(refVelocity, elementVelocity);
176 elementVelocity /= geometry.integrationElement(localPos);
177
178 velocity[eIdxGlobal] = elementVelocity;
179 }
180
181 writer.attachCellData(velocity, "velocity", dim);
182
183 return;
184 }
185private:
186 Problem &problem_;
187 const GravityVector& gravity_;
188 Scalar density_;
189 Scalar viscosity_;
190};
191
201template<class TypeTag>
202void FVVelocity1P<TypeTag>::calculateVelocity(const Intersection& intersection, CellData& cellData)
203{
204 auto elementI = intersection.inside();
205 auto elementJ = intersection.outside();
206
207 int eIdxGlobalJ = problem_.variables().index(elementJ);
208
209 CellData& cellDataJ = problem_.variables().cellData(eIdxGlobalJ);
210
211 // get global coordinates of cell centers
212 const GlobalPosition& globalPosI = elementI.geometry().center();
213 const GlobalPosition& globalPosJ = elementJ.geometry().center();
214
215 //get face index
216 int isIndexI = intersection.indexInInside();
217 int isIndexJ = intersection.indexInOutside();
218
219 //get face normal
220 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
221
222 // distance vector between barycenters
223 GlobalPosition distVec = globalPosJ - globalPosI;
224
225 // compute distance between cell centers
226 Scalar dist = distVec.two_norm();
227
228 // compute vectorized permeabilities
229 DimMatrix meanPermeability(0);
230
231 problem_.spatialParams().meanK(meanPermeability, problem_.spatialParams().intrinsicPermeability(elementI),
232 problem_.spatialParams().intrinsicPermeability(elementJ));
233
234 Dune::FieldVector < Scalar, dim > permeability(0);
235 meanPermeability.mv(unitOuterNormal, permeability);
236
237 permeability /= viscosity_;
238
239 // calculate potential gradients
240 Scalar potential = (cellData.pressure() - cellDataJ.pressure()) / dist;
241
242 potential += density_ * (unitOuterNormal * gravity_);
243
244 // store potentials for further calculations (velocity, saturation, ...)
245 cellData.fluxData().setPotential(isIndexI, potential);
246 cellDataJ.fluxData().setPotential(isIndexJ, -potential);
247
248 // calculate the gravity term
249 VelocityVector velocity(permeability);
250 velocity *= (cellData.pressure() - cellDataJ.pressure()) / dist;
251
252 GravityVector gravityTerm(unitOuterNormal);
253 gravityTerm *= (gravity_ * permeability) * density_;
254
255 velocity += gravityTerm;
256
257 // store velocities
258 cellData.fluxData().setVelocity(isIndexI, velocity);
259 cellData.fluxData().setVelocityMarker(isIndexI);
260
261 cellDataJ.fluxData().setVelocity(isIndexJ, velocity);
262 cellDataJ.fluxData().setVelocityMarker(isIndexJ);
263 return;
264}
265
275template<class TypeTag>
276void FVVelocity1P<TypeTag>::calculateVelocityOnBoundary(const Intersection& intersection, CellData& cellData)
277{
278 auto element = intersection.inside();
279
280 // get face index
281 int isIndex = intersection.indexInInside();
282
283 // get face normal
284 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
285
286 BoundaryTypes bcType;
287 // get boundary type
288 problem_.boundaryTypes(bcType, intersection);
289 PrimaryVariables boundValues(0.0);
290
291 if (bcType.isDirichlet(pressEqIdx))
292 {
293 problem_.dirichlet(boundValues, intersection);
294
295 // get global coordinates of cell centers
296 const GlobalPosition& globalPosI = element.geometry().center();
297
298 // center of face in global coordinates
299 const GlobalPosition& globalPosJ = intersection.geometry().center();
300
301 // distance vector between barycenters
302 GlobalPosition distVec = globalPosJ - globalPosI;
303
304 // compute distance between cell centers
305 Scalar dist = distVec.two_norm();
306
307 // permeability vector at boundary
308 // compute vectorized permeabilities
309 DimMatrix meanPermeability(0);
310
311 problem_.spatialParams().meanK(meanPermeability, problem_.spatialParams().intrinsicPermeability(element));
312
313 // multiply with normal vector at the boundary
314 Dune::FieldVector < Scalar, dim > permeability(0);
315 meanPermeability.mv(unitOuterNormal, permeability);
316 permeability /= viscosity_;
317
318 Scalar pressBound = boundValues;
319
320 // calculate potential gradients
321 Scalar potential = (cellData.pressure() - pressBound) / dist;
322
323 potential += density_ * (unitOuterNormal * gravity_);
324
325 // store potentials for further calculations (velocity, saturation, ...)
326 cellData.fluxData().setPotential(isIndex, potential);
327
328 // calculate the gravity term
329 VelocityVector velocity(permeability);
330 velocity *= (cellData.pressure() - pressBound) / dist;
331
332 GravityVector gravityTerm(unitOuterNormal);
333 gravityTerm *= (gravity_ * permeability) * density_;
334
335 velocity += gravityTerm;
336
337 // store velocities
338 cellData.fluxData().setVelocity(isIndex, velocity);
339 cellData.fluxData().setVelocityMarker(isIndex);
340
341 } // end Dirichlet boundary
342
343 else
344 {
345 problem_.neumann(boundValues, intersection);
346 VelocityVector velocity(unitOuterNormal);
347
348 velocity *= boundValues[pressEqIdx] / density_;
349
350 // store potential gradients for further calculations
351 cellData.fluxData().setPotential(isIndex, boundValues[pressEqIdx]);
352
353 //store velocity
354 cellData.fluxData().setVelocity(isIndex, velocity);
355 cellData.fluxData().setVelocityMarker(isIndex);
356 } // end Neumann boundary
357 return;
358}
359} // end namespace Dumux
360#endif
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
Definition: adapt.hh:29
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.