25#ifndef DUMUX_MATRIX_ROBLEM_HH
26#define DUMUX_MATRIX_ROBLEM_HH
28#include <dune/geometry/quadraturerules.hh>
29#include <dune/geometry/referenceelements.hh>
30#include <dune/grid/yaspgrid.hh>
31#include <dune/localfunctions/lagrange/pqkfactory.hh>
49template <
class TypeTag>
60template<
class TypeTag>
61struct Grid<TypeTag, TTag::Matrix> {
using type = Dune::YaspGrid<3, Dune::EquidistantOffsetCoordinates<GetPropType<TypeTag, Properties::Scalar>, 3> >; };
63template<
class TypeTag>
65template<
class TypeTag>
67template<
class TypeTag>
69template<
class TypeTag>
71template<
class TypeTag>
73template<
class TypeTag>
77template<
class TypeTag>
81template<
class TypeTag>
85template<
class TypeTag>
93template<
class TypeTag>
106template <
class TypeTag>
112 using GridView =
typename GridGeometry::GridView;
113 using FVElementGeometry =
typename GridGeometry::LocalView;
114 using SubControlVolume =
typename GridGeometry::SubControlVolume;
125 dim = GridView::dimension,
126 dimWorld = GridView::dimensionworld
129 using Element =
typename GridView::template Codim<0>::Entity;
130 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
136 std::shared_ptr<typename ParentType::SpatialParams>
spatialParams,
143 name_ = getParam<std::string>(
"Vtk.OutputName") +
"_" + getParamFromGroup<std::string>(this->
paramGroup(),
"Problem.Name");
156 const std::string&
name()
const
165 {
return 273.15 + 37; }
182 BoundaryTypes values;
183 values.setAllNeumann();
225 template<
class ElementVolumeVariables>
227 const Element &element,
228 const FVElementGeometry& fvGeometry,
229 const ElementVolumeVariables& elemVolVars,
230 const SubControlVolume &scv)
const
233 const Scalar pressure3D = this->
couplingManager().bulkPriVars(source.id())[Indices::pressureIdx];
234 const Scalar pressure1D = this->
couplingManager().lowDimPriVars(source.id())[Indices::pressureIdx];
237 const Scalar meanDistance = this->
couplingManager().averageDistance(source.id());
238 static const Scalar matrixPerm = getParamFromGroup<Scalar>(
"Matrix",
"SpatialParams.Permeability");
239 static const Scalar rho = getParam<Scalar>(
"Component.LiquidDensity");
240 static const Scalar mu = getParam<Scalar>(
"Component.LiquidKinematicViscosity")*rho;
241 const Scalar sourceValue = rho*(pressure1D - pressure3D)/meanDistance*matrixPerm/mu;
254 {
return PrimaryVariables({1e5}); }
261 for (
const auto& element : elements(this->
gridGeometry().gridView()))
264 fvGeometry.bindElement(element);
266 auto elemVolVars =
localView(gridVars.curGridVolVars());
267 elemVolVars.bindElement(element, fvGeometry, sol);
269 for (
auto&& scv : scvs(fvGeometry))
271 auto pointSources = this->
scvPointSources(element, fvGeometry, elemVolVars, scv);
272 pointSources *= scv.volume()*elemVolVars[scv].extrusionFactor();
277 std::cout <<
"Global integrated source (3D): " <<
source <<
'\n';
282 {
return *couplingManager_; }
286 static constexpr Scalar eps_ = 1.5e-7;
289 std::shared_ptr<CouplingManager> couplingManager_;
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
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.
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 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 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
The matrix problem.
Definition: problem_matrix.hh:108
PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
Evaluates the initial value for a control volume.
Definition: problem_matrix.hh:253
const CouplingManager & couplingManager() const
Get the coupling manager.
Definition: problem_matrix.hh:281
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_matrix.hh:226
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_matrix.hh:180
const std::string & name() const
The problem name.
Definition: problem_matrix.hh:156
void computeSourceIntegral(const SolutionVector &sol, const GridVariables &gridVars)
Definition: problem_matrix.hh:258
MatrixProblem(std::shared_ptr< const GridGeometry > gridGeometry, std::shared_ptr< typename ParentType::SpatialParams > spatialParams, std::shared_ptr< CouplingManager > couplingManager, const std::string ¶mGroup="Matrix")
Definition: problem_matrix.hh:135
void addPointSources(std::vector< PointSource > &pointSources) const
Applies a vector of point sources which are possibly solution dependent.
Definition: problem_matrix.hh:204
Scalar temperature() const
Returns the temperature within the domain [K].
Definition: problem_matrix.hh:164
Definition: problem_matrix.hh:56
std::tuple< OneP, CCTpfaModel > InheritsFrom
Definition: problem_matrix.hh:56
Dune::YaspGrid< 3, Dune::EquidistantOffsetCoordinates< GetPropType< TypeTag, Properties::Scalar >, 3 > > type
Definition: problem_matrix.hh:61
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: problem_matrix.hh:88
Definition of the spatial parameters for the matrix and fracture problem.
Definition: multidomain/embedded/2d3d/1p_1p/spatialparams.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,...
Definition of the spatial parameters for the MaxwellStefan problem.