26#ifndef DUMUX_BLOOD_FLOW_PROBLEM_HH
27#define DUMUX_BLOOD_FLOW_PROBLEM_HH
29#include <dune/foamgrid/foamgrid.hh>
49template <
class TypeTag>
class BloodFlowProblem;
61template<
class TypeTag>
62struct Grid<TypeTag, TTag::BloodFlow> {
using type = Dune::FoamGrid<1, 3>; };
64template<
class TypeTag>
66template<
class TypeTag>
68template<
class TypeTag>
70template<
class TypeTag>
72template<
class TypeTag>
74template<
class TypeTag>
78template<
class TypeTag>
82template<
class TypeTag>
90template<
class TypeTag>
94template<
class TypeTag>
106template <
class TypeTag>
116 using GridView =
typename GridGeometry::GridView;
117 using FVElementGeometry =
typename GridGeometry::LocalView;
118 using SubControlVolume =
typename GridGeometry::SubControlVolume;
121 using Element =
typename GridView::template Codim<0>::Entity;
122 using GlobalPosition =
typename GridGeometry::GlobalCoordinate;
132 name_ = getParam<std::string>(
"Vtk.OutputName") +
"_" + getParamFromGroup<std::string>(this->
paramGroup(),
"Problem.Name");
133 p_in_ = getParam<Scalar>(
"BoundaryConditions1D.PressureInput");
134 delta_p_ = getParam<Scalar>(
"BoundaryConditions1D.DeltaPressure");
137 for (
const auto& element : elements(this->
gridGeometry().gridView()))
140 fvGeometry.bindElement(element);
142 for (
auto&& scv : scvs(fvGeometry))
143 exactPressure_[scv.dofIndex()] =
exactSolution(scv.dofPosition());
153 template<
class ElementSolution>
155 const SubControlVolume &scv,
156 const ElementSolution& elemSol)
const
158 const auto eIdx = this->
gridGeometry().elementMapper().index(element);
160 return M_PI*radius*radius;
173 const std::string&
name()
const
181 {
return 273.15 + 37.0; }
197 BoundaryTypes values;
198 values.setAllDirichlet();
211 PrimaryVariables values;
212 if(globalPos[2] > 0.5)
213 values[Indices::pressureIdx] = p_in_;
215 values[Indices::pressureIdx] = p_in_-delta_p_;
258 template<
class ElementVolumeVariables,
260 std::enable_if_t<!enable, int> = 0>
262 const Element &element,
263 const FVElementGeometry& fvGeometry,
264 const ElementVolumeVariables& elemVolVars,
265 const SubControlVolume &scv)
const
268 const Scalar pressure1D = this->
couplingManager().lowDimPriVars(source.id())[Indices::pressureIdx];
269 const Scalar pressure3D = this->
couplingManager().bulkPriVars(source.id())[Indices::pressureIdx];
273 const Scalar beta = 2*M_PI/(2*M_PI + std::log(radius));
274 const Scalar sourceValue = beta*(pressure3D - pressure1D);
280 template<
class ElementVolumeVariables,
282 std::enable_if_t<enable, int> = 0>
284 const Element &element,
285 const FVElementGeometry& fvGeometry,
286 const ElementVolumeVariables& elemVolVars,
287 const SubControlVolume &scv)
const
289 static const Scalar kernelWidth = getParam<Scalar>(
"MixedDimension.KernelWidth");
290 static const Scalar kernelFactor = std::log(kernelWidth)-9.0/10.0;
293 const Scalar pressure1D = this->
couplingManager().lowDimPriVars(source.id())[Indices::pressureIdx];
294 const Scalar pressure3D = this->
couplingManager().bulkPriVars(source.id())[Indices::pressureIdx];
298 const Scalar beta = 2*M_PI/(2*M_PI + std::log(radius));
299 const Scalar sourceValue = beta*(pressure3D / kernelFactor * std::log(radius) - pressure1D);
306 template<
class MatrixBlock,
class VolumeVariables>
308 const Element& element,
309 const FVElementGeometry& fvGeometry,
310 const VolumeVariables& curElemVolVars,
311 const SubControlVolume& scv)
const
313 const auto eIdx = this->
gridGeometry().elementMapper().index(element);
315 auto key = std::make_pair(eIdx, 0);
323 for (
const auto&
source : pointSources)
324 block[0][0] -= this->
couplingManager().pointSourceDerivative(source, Dune::index_constant<1>{}, Dune::index_constant<1>{});
335 {
return PrimaryVariables(0.0); }
341 {
return 1 + globalPos[2]; }
347 PrimaryVariables
source(0.0);
348 for (
const auto& element : elements(this->
gridGeometry().gridView()))
351 fvGeometry.bindElement(element);
353 auto elemVolVars =
localView(gridVars.curGridVolVars());
354 elemVolVars.bindElement(element, fvGeometry, sol);
356 for (
auto&& scv : scvs(fvGeometry))
358 auto pointSources = this->
scvPointSources(element, fvGeometry, elemVolVars, scv);
359 pointSources *= scv.volume()*elemVolVars[scv].extrusionFactor();
364 std::cout <<
"Global integrated source (1D): " <<
source <<
'\n';
372 template<
class VtkOutputModule>
375 vtk.
addField(exactPressure_,
"exact pressure");
380 {
return *couplingManager_; }
383 std::vector<Scalar> exactPressure_;
385 static constexpr Scalar eps_ = 1.5e-7;
391 std::shared_ptr<CouplingManager> couplingManager_;
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Defines a type tag and some properties for models using the box scheme.
Properties for all models using cell-centered finite volume scheme with TPFA.
Setting constant fluid properties via the input file.
A liquid phase consisting of a single component.
Coupling manager for low-dimensional domains embedded in the bulk domain. Point sources on each integ...
Definition of the spatial parameters for the blood flow problem.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
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
Base class for all finite-volume problems.
Definition: common/fvproblem.hh:50
NumEqVector scvPointSources(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv) const
Adds contribution of point sources for a specific sub control volume to the values....
Definition: common/fvproblem.hh:435
const std::string & paramGroup() const
The parameter group in which to retrieve runtime parameters.
Definition: common/fvproblem.hh:592
NumEqVector source(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv) const
Evaluate the source term for all phases within a given sub-control-volume.
Definition: common/fvproblem.hh:327
const PointSourceMap & pointSourceMap() const
Get the point source map. It stores the point sources per scv.
Definition: common/fvproblem.hh:505
const GridGeometry & gridGeometry() const
The finite volume grid geometry.
Definition: common/fvproblem.hh:588
The DUNE grid type.
Definition: common/properties.hh:57
Property to specify the type of a problem which has to be solved.
Definition: common/properties.hh:69
Definition: common/properties.hh:91
Definition: common/properties.hh:169
If disabled, the volume variables are not stored (reduces memory, but is slower)
Definition: common/properties.hh:178
specifies if data on flux vars should be saved (faster, but more memory consuming)
Definition: common/properties.hh:188
specifies if the parameters for the advective fluxes depend on the solution
Definition: common/properties.hh:210
specifies if the parameters for the diffusive fluxes depend on the solution
Definition: common/properties.hh:214
specifies if the parameters for the heat conduction fluxes depend on the solution
Definition: common/properties.hh:218
The type of the spatial parameters object.
Definition: common/properties.hh:221
The type of the fluid system to use.
Definition: common/properties.hh:223
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: io/vtkoutputmodule.hh:66
void addField(const Vector &v, const std::string &name, FieldType fieldType=FieldType::automatic)
Definition: io/vtkoutputmodule.hh:161
A liquid phase consisting of a single component.
Definition: 1pliquid.hh:46
Element-wise calculation of the residual and its derivatives for a single-phase, incompressible,...
Definition: 1p/incompressiblelocalresidual.hh:41
Base class for all fully implicit porous media problems.
Definition: dumux/porousmediumflow/problem.hh:39
SpatialParams & spatialParams()
Returns the spatial parameters object.
Definition: dumux/porousmediumflow/problem.hh:146
Exact solution 1D-3D.
Definition: problem_bloodflow.hh:108
Scalar extrusionFactor(const Element &element, const SubControlVolume &scv, const ElementSolution &elemSol) const
Returns how much the domain is extruded at a given sub-control volume.
Definition: problem_bloodflow.hh:154
const std::string & name() const
The problem name.
Definition: problem_bloodflow.hh:173
Scalar exactSolution(const GlobalPosition &globalPos) const
The exact pressure solution.
Definition: problem_bloodflow.hh:340
void addVtkOutputFields(VtkOutputModule &vtk) const
Adds additional VTK output data to the VTKWriter.
Definition: problem_bloodflow.hh:373
PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
Evaluates the initial value for a control volume.
Definition: problem_bloodflow.hh:334
void addPointSources(std::vector< PointSource > &pointSources) const
Applies a vector of point sources which are possibly solution dependent.
Definition: problem_bloodflow.hh:237
BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
Specifies which kind of boundary condition should be used for which equation on a given boundary segm...
Definition: problem_bloodflow.hh:195
const CouplingManager & couplingManager() const
Get the coupling manager.
Definition: problem_bloodflow.hh:379
BloodFlowProblem(std::shared_ptr< const GridGeometry > gridGeometry, std::shared_ptr< CouplingManager > couplingManager)
Definition: problem_bloodflow.hh:126
Scalar temperature() const
Returns the temperature within the domain in [K].
Definition: problem_bloodflow.hh:180
void computeSourceIntegral(const SolutionVector &sol, const GridVariables &gridVars)
Definition: problem_bloodflow.hh:345
void pointSource(PointSource &source, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv) const
Evaluates the point sources (added by addPointSources) for all phases within a given sub-control-volu...
Definition: problem_bloodflow.hh:261
void addSourceDerivatives(MatrixBlock &block, const Element &element, const FVElementGeometry &fvGeometry, const VolumeVariables &curElemVolVars, const SubControlVolume &scv) const
Definition: problem_bloodflow.hh:307
PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
Evaluates the boundary conditions for a Dirichlet control volume.
Definition: problem_bloodflow.hh:209
Definition: problem_bloodflow.hh:55
std::tuple< OneP > InheritsFrom
Definition: problem_bloodflow.hh:55
Definition: problem_bloodflow.hh:56
std::tuple< BloodFlow, CCTpfaModel > InheritsFrom
Definition: problem_bloodflow.hh:56
Definition: problem_bloodflow.hh:57
std::tuple< BloodFlow, BoxModel > InheritsFrom
Definition: problem_bloodflow.hh:57
Dune::FoamGrid< 1, 3 > type
Definition: problem_bloodflow.hh:62
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: problem_bloodflow.hh:85
Definition of the spatial parameters for the blood flow problem.
Definition: spatialparams_bloodflow.hh:39
Declares all properties used in Dumux.
A single-phase, isothermal flow model using the fully implicit scheme.
Base class for all porous media problems.
Element-wise calculation of the residual and its derivatives for a single-phase, incompressible,...