25#ifndef DUMUX_TEST_TPFAFACETCOUPLING_ONEP_BULKPROBLEM_HH
26#define DUMUX_TEST_TPFAFACETCOUPLING_ONEP_BULKPROBLEM_HH
28#include <dune/alugrid/grid.hh>
43template<
class TypeTag>
class OnePBulkProblem;
55template<
class TypeTag>
56struct Grid<TypeTag, TTag::OnePBulk> {
using type = Dune::ALUGrid<2, 2, Dune::cube, Dune::nonconforming>; };
58template<
class TypeTag>
61template<
class TypeTag>
69template<
class TypeTag>
85template<
class TypeTag>
91 using PrimaryVariables =
typename GridVariables::PrimaryVariables;
92 using Scalar =
typename GridVariables::Scalar;
94 using GridGeometry =
typename GridVariables::GridGeometry;
95 using FVElementGeometry =
typename GridGeometry::LocalView;
96 using SubControlVolume =
typename GridGeometry::SubControlVolume;
97 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
98 using GridView =
typename GridGeometry::GridView;
99 using Element =
typename GridView::template Codim<0>::Entity;
100 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
108 std::shared_ptr<typename ParentType::SpatialParams>
spatialParams,
109 std::shared_ptr<CouplingManager> couplingManagerPtr,
112 , couplingManagerPtr_(couplingManagerPtr)
113 , lowDimPermeability_(
getParam<Scalar>(
"LowDim.SpatialParams.Permeability"))
114 , aperture_(
getParam<Scalar>(
"Problem.FractureAperture"))
116 problemName_ = getParam<std::string>(
"Vtk.OutputName") +
"_" + getParamFromGroup<std::string>(this->
paramGroup(),
"Problem.Name");
122 const std::string&
name()
const
130 BoundaryTypes values;
131 values.setAllDirichlet();
138 BoundaryTypes values;
139 values.setAllNeumann();
148 Scalar u = (1.0 - lowDimPermeability_)*cos(globalPos[0])*cosh(aperture_/2);
149 return NumEqVector(u);
153 Scalar
exact(
const GlobalPosition& globalPos)
const
157 const auto x = globalPos[0];
158 const auto y = globalPos[1];
159 return lowDimPermeability_*cos(x)*cosh(y) + (1.0 - lowDimPermeability_)*cos(x)*cosh(aperture_/2);
170 const auto x = globalPos[0];
171 const auto y = globalPos[1];
173 GlobalPosition gradU;
174 gradU[0] = -lowDimPermeability_*sin(x)*cosh(y) + (lowDimPermeability_ - 1.0)*sin(x)*cosh(aperture_/2);
175 gradU[1] = lowDimPermeability_*cos(x)*sinh(y);
182 {
return PrimaryVariables(
exact(globalPos)); }
185 template<
class ElementVolumeVariables,
class ElementFluxVarsCache>
187 const FVElementGeometry& fvGeometry,
188 const ElementVolumeVariables& elemVolVars,
189 const ElementFluxVarsCache& elemFluxVarsCache,
190 const SubControlVolumeFace& scvf)
const
192 auto pos = scvf.ipGlobal();
193 const Scalar k = this->
spatialParams().permeabilityAtPos(pos);
195 return NumEqVector( -1.0*k*(gradU*scvf.unitOuterNormal()) );
200 {
return PrimaryVariables(1.0); }
208 {
return *couplingManagerPtr_; }
211 std::shared_ptr<CouplingManager> couplingManagerPtr_;
212 Scalar lowDimPermeability_;
214 std::string problemName_;
A liquid phase consisting of a single component.
Setting constant fluid properties via the input file.
T getParam(Args &&... args)
A free function to get a parameter from the parameter tree singleton.
Definition: parameters.hh:428
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
Class to specify the type of a boundary.
Definition: common/boundarytypes.hh:38
Base class for all finite-volume problems.
Definition: common/fvproblem.hh:50
const std::string & paramGroup() const
The parameter group in which to retrieve runtime parameters.
Definition: common/fvproblem.hh:592
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
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
Definition: multidomain/couplingmanager.hh:46
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 spatial parameters class for the test problem using the 1p cc model.
Definition: multidomain/boundary/stokesdarcy/1p2c_1p2c/spatialparams.hh:41
Test problem for the incompressible one-phase model with coupling across the bulk grid facets.
Definition: problem_1p_bulk.hh:90
PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
Evaluates the Dirichlet boundary condition for a given position.
Definition: analytical/problem_bulk.hh:181
NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
Evaluates the source term at a given position.
Definition: analytical/problem_bulk.hh:144
const std::string & name() const
The problem name.
Definition: analytical/problem_bulk.hh:122
BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
Specifies the type of boundary condition at a given position.
Definition: analytical/problem_bulk.hh:128
BoundaryTypes interiorBoundaryTypes(const Element &element, const SubControlVolumeFace &scvf) const
Specifies the type of interior boundary condition at a given position.
Definition: analytical/problem_bulk.hh:136
NumEqVector neumann(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVarsCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Evaluates the Neumann boundary condition for a boundary segment.
Definition: analytical/problem_bulk.hh:186
Scalar temperature() const
Returns the temperature in in the domain.
Definition: analytical/problem_bulk.hh:203
const CouplingManager & couplingManager() const
Returns reference to the coupling manager.
Definition: analytical/problem_bulk.hh:207
GlobalPosition exactGradient(const GlobalPosition &globalPos) const
Evaluates the exact gradient at a given position.
Definition: analytical/problem_bulk.hh:163
Scalar exact(const GlobalPosition &globalPos) const
Evaluates the exact solution at a given position.
Definition: analytical/problem_bulk.hh:153
OnePBulkProblem(std::shared_ptr< const GridGeometry > gridGeometry, std::shared_ptr< typename ParentType::SpatialParams > spatialParams, std::shared_ptr< CouplingManager > couplingManagerPtr, const std::string ¶mGroup="Bulk")
Definition: analytical/problem_bulk.hh:107
PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
Evaluates the initial conditions.
Definition: analytical/problem_bulk.hh:199
Definition: analytical/problem_bulk.hh:49
std::tuple< OneP > InheritsFrom
Definition: analytical/problem_bulk.hh:49
Definition: analytical/problem_bulk.hh:50
std::tuple< CCTpfaFacetCouplingModel, OnePBulk > InheritsFrom
Definition: analytical/problem_bulk.hh:50
Definition: analytical/problem_bulk.hh:51
std::tuple< BoxFacetCouplingModel, OnePBulk > InheritsFrom
Definition: analytical/problem_bulk.hh:51
Dune::ALUGrid< 2, 2, Dune::cube, Dune::nonconforming > type
Definition: analytical/problem_bulk.hh:56
A single-phase, isothermal flow model using the fully implicit scheme.
Properties (and default properties) for all models using the box scheme together with coupling across...
Properties (and default properties) for all models using cell-centered finite volume scheme with TPFA...
Base class for all porous media problems.
Definition of the spatial parameters for the MaxwellStefan problem.