105 const Element& element,
106 const FVElementGeometry& fvGeometry,
107 const ElementVolumeVariables& elemVolVars,
108 const ElementFluxVarsCache& elemFluxVarsCache,
111 using Velocity =
typename VelocityVector::value_type;
113 const auto geometry = element.geometry();
114 const Dune::GeometryType geomType = geometry.type();
117 auto upwindTerm = [phaseIdx](
const auto& volVars) {
return volVars.mobility(phaseIdx); };
119 if(isBox && dim == 1)
121 Velocity tmpVelocity(0.0);
122 tmpVelocity = (geometry.corner(1) - geometry.corner(0));
123 tmpVelocity /= tmpVelocity.two_norm();
125 for (
auto&& scvf : scvfs(fvGeometry))
131 FluxVariables fluxVars;
132 fluxVars.init(problem_, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
136 Scalar localArea = scvfReferenceArea_(geomType, scvf.index());
137 Scalar flux = fluxVars.advectiveFlux(phaseIdx, upwindTerm) / localArea;
138 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
139 flux /= insideVolVars.extrusionFactor();
142 const int eIdxGlobal = gridGeometry_.elementMapper().index(element);
143 velocity[eIdxGlobal] = tmpVelocity;
149 const auto referenceElement = ReferenceElements::general(geomType);
150 const auto& localPos = referenceElement.position(0, 0);
151 const auto jacobianT2 = geometry.jacobianTransposed(localPos);
155 using ScvVelocities = Dune::BlockVector<Velocity>;
156 ScvVelocities scvVelocities(fvGeometry.numScv());
159 for (
auto&& scvf : scvfs(fvGeometry))
165 const auto localPosIP = geometry.local(scvf.ipGlobal());
168 const auto jacobianT1 = geometry.jacobianTransposed(localPosIP);
169 const auto globalNormal = scvf.unitOuterNormal();
170 GlobalPosition localNormal(0);
171 jacobianT1.mv(globalNormal, localNormal);
172 localNormal /= localNormal.two_norm();
175 FluxVariables fluxVars;
176 fluxVars.init(problem_, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
180 Scalar localArea = scvfReferenceArea_(geomType, scvf.index());
181 Scalar flux = fluxVars.advectiveFlux(phaseIdx, upwindTerm) / localArea;
182 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
183 flux /= insideVolVars.extrusionFactor();
186 Velocity tmpVelocity = localNormal;
189 scvVelocities[scvf.insideScvIdx()] += tmpVelocity;
190 scvVelocities[scvf.outsideScvIdx()] += tmpVelocity;
194 for (
auto&& scv : scvs(fvGeometry))
196 int vIdxGlobal = scv.dofIndex();
199 Velocity scvVelocity(0);
201 jacobianT2.mtv(scvVelocities[scv.indexInElement()], scvVelocity);
202 scvVelocity /= geometry.integrationElement(localPos)*cellNum_[vIdxGlobal];
204 velocity[vIdxGlobal] += scvVelocity;
213 const int numScvfsPerFace = isMpfa ? element.template subEntity<1>(0).geometry().corners() : 1;
215 if (fvGeometry.numScvf() != element.subEntities(1)*numScvfsPerFace)
216 DUNE_THROW(Dune::NotImplemented,
"Velocity output for non-conforming grids");
218 if (!geomType.isCube() && !geomType.isSimplex())
219 DUNE_THROW(Dune::NotImplemented,
"Velocity output for other geometry types than cube and simplex");
225 std::vector<bool> handledScvf;
227 handledScvf.resize(element.subEntities(1),
false);
230 std::vector<unsigned int> scvfIndexInInside(fvGeometry.numScvf());
231 int localScvfIdx = 0;
232 for (
const auto& intersection : intersections(gridGeometry_.gridView(), element))
235 if (handledScvf[intersection.indexInInside()])
238 if (intersection.neighbor() || intersection.boundary())
240 for (
int i = 0; i < numScvfsPerFace; ++i)
241 scvfIndexInInside[localScvfIdx++] = intersection.indexInInside();
245 handledScvf[intersection.indexInInside()] =
true;
249 std::vector<Scalar> scvfFluxes(element.subEntities(1), 0.0);
251 for (
auto&& scvf : scvfs(fvGeometry))
253 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
254 if (!scvf.boundary())
256 FluxVariables fluxVars;
257 fluxVars.init(problem_, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
258 scvfFluxes[scvfIndexInInside[localScvfIdx]] += fluxVars.advectiveFlux(phaseIdx, upwindTerm)/insideVolVars.extrusionFactor();
262 auto bcTypes = problemBoundaryTypes_(element, scvf);
263 if (bcTypes.hasOnlyDirichlet())
265 FluxVariables fluxVars;
266 fluxVars.init(problem_, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
267 scvfFluxes[scvfIndexInInside[localScvfIdx]] += fluxVars.advectiveFlux(phaseIdx, upwindTerm)/insideVolVars.extrusionFactor();
282 for (
auto&& scvf : scvfs(fvGeometry))
286 auto bcTypes = problemBoundaryTypes_(element, scvf);
287 if (bcTypes.hasNeumann())
293 if (stationaryVelocityField)
295 const auto flux = AdvectionType::flux(problem_, element, fvGeometry, elemVolVars,
296 scvf, phaseIdx, elemFluxVarsCache);
298 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
299 scvfFluxes[scvfIndexInInside[localScvfIdx]] += flux / insideVolVars.extrusionFactor();
304 const auto neumannFlux = problem_.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
305 using NumEqVector = std::decay_t<
decltype(neumannFlux)>;
306 if (Dune::FloatCmp::eq<NumEqVector, Dune::FloatCmp::CmpStyle::absolute>(neumannFlux, NumEqVector(0.0), 1e-30))
307 scvfFluxes[scvfIndexInInside[localScvfIdx]] = 0;
311 else if (dim == 1 || geomType.isCube())
313 const auto fIdx = scvfIndexInInside[localScvfIdx];
314 const auto fIdxOpposite = fIdx%2 ? fIdx-1 : fIdx+1;
315 scvfFluxes[fIdx] = -scvfFluxes[fIdxOpposite];
319 else if (geomType.isSimplex())
320 scvfFluxes[scvfIndexInInside[localScvfIdx]] = 0;
323 DUNE_THROW(Dune::NotImplemented,
"Velocity computation at Neumann boundaries for cell-centered and prism/pyramid");
333 Velocity refVelocity;
336 if (dim == 1 || geomType.isCube())
338 for (
int i = 0; i < dim; i++)
339 refVelocity[i] = 0.5 * (scvfFluxes[2*i + 1] - scvfFluxes[2*i]);
342 else if (geomType.isSimplex())
344 for (
int dimIdx = 0; dimIdx < dim; dimIdx++)
346 refVelocity[dimIdx] = -scvfFluxes[dim - 1 - dimIdx];
347 for (
int fIdx = 0; fIdx < dim + 1; fIdx++)
348 refVelocity[dimIdx] += scvfFluxes[fIdx]/(dim + 1);
353 DUNE_THROW(Dune::NotImplemented,
"Velocity computation for cell-centered and prism/pyramid");
355 Velocity scvVelocity(0);
356 jacobianT2.mtv(refVelocity, scvVelocity);
358 scvVelocity /= geometry.integrationElement(localPos);
360 int eIdxGlobal = gridGeometry_.elementMapper().index(element);
362 velocity[eIdxGlobal] = scvVelocity;