3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
problem_tissue.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 *****************************************************************************/
26#ifndef DUMUX_TISSUE_PROBLEM_HH
27#define DUMUX_TISSUE_PROBLEM_HH
28
29#include <dune/grid/yaspgrid.hh>
30
31#include <dune/geometry/quadraturerules.hh>
32#include <dune/geometry/referenceelements.hh>
33#include <dune/localfunctions/lagrange/pqkfactory.hh>
34
35#include <dumux/common/math.hh>
40
44
47
48#include <dumux/multidomain/embedded/couplingmanager1d3d.hh> // for coupling mode
49
51
52namespace Dumux {
53
54template <class TypeTag>
55class TissueProblem;
56
57namespace Properties {
58
59// Create new type tags
60namespace TTag {
61struct Tissue { using InheritsFrom = std::tuple<OneP>; };
62struct TissueCC { using InheritsFrom = std::tuple<Tissue, CCTpfaModel>; };
63struct TissueBox { using InheritsFrom = std::tuple<Tissue, BoxModel>; };
64} // end namespace TTag
65
66// Set the grid type
67template<class TypeTag>
68struct Grid<TypeTag, TTag::Tissue> { using type = Dune::YaspGrid<3, Dune::EquidistantOffsetCoordinates<GetPropType<TypeTag, Properties::Scalar>, 3> >; };
69
70template<class TypeTag>
71struct EnableGridGeometryCache<TypeTag, TTag::Tissue> { static constexpr bool value = true; };
72template<class TypeTag>
73struct EnableGridVolumeVariablesCache<TypeTag, TTag::Tissue> { static constexpr bool value = true; };
74template<class TypeTag>
75struct EnableGridFluxVariablesCache<TypeTag, TTag::Tissue> { static constexpr bool value = true; };
76template<class TypeTag>
77struct SolutionDependentAdvection<TypeTag, TTag::Tissue> { static constexpr bool value = false; };
78template<class TypeTag>
79struct SolutionDependentMolecularDiffusion<TypeTag, TTag::Tissue> { static constexpr bool value = false; };
80template<class TypeTag>
81struct SolutionDependentHeatConduction<TypeTag, TTag::Tissue> { static constexpr bool value = false; };
82
83// Set the problem property
84template<class TypeTag>
85struct Problem<TypeTag, TTag::Tissue> { using type = TissueProblem<TypeTag>; };
86
87// Set the problem property
88template<class TypeTag>
89struct LocalResidual<TypeTag, TTag::Tissue> { using type = OnePIncompressibleLocalResidual<TypeTag>; };
90
91// the fluid system
92template<class TypeTag>
93struct FluidSystem<TypeTag, TTag::Tissue>
94{
97};
98
99// Set the spatial parameters
100template<class TypeTag>
101struct SpatialParams<TypeTag, TTag::Tissue>
102{
105};
106} // end namespace Properties
107
108
114template <class TypeTag>
116{
120 using FVElementGeometry = typename GridGeometry::LocalView;
121 using SubControlVolume = typename GridGeometry::SubControlVolume;
122 using GridView = typename GridGeometry::GridView;
130 using Element = typename GridView::template Codim<0>::Entity;
131 using GlobalPosition = typename GridGeometry::GlobalCoordinate;
132
134
135public:
136 TissueProblem(std::shared_ptr<const GridGeometry> gridGeometry,
137 std::shared_ptr<CouplingManager> couplingManager)
138 : ParentType(gridGeometry, "Tissue")
139 , couplingManager_(couplingManager)
140 {
141 // read parameters from input file
142 name_ = getParam<std::string>("Vtk.OutputName") + "_" + getParamFromGroup<std::string>(this->paramGroup(), "Problem.Name");
143
144 exactPressure_.resize(this->gridGeometry().numDofs());
145 for (const auto& element : elements(this->gridGeometry().gridView()))
146 {
147 auto fvGeometry = localView(this->gridGeometry());
148 fvGeometry.bindElement(element);
149
150 for (auto&& scv : scvs(fvGeometry))
151 exactPressure_[scv.dofIndex()] = exactSolution(scv.dofPosition());
152 }
153 }
154
160 const std::string& name() const
161 { return name_; }
162
168 Scalar temperature() const
169 { return 273.15 + 37; } // in [K]
170
171 // \}
172
176 // \{
177
184 BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
185 {
186 BoundaryTypes values;
187 values.setAllDirichlet();
188 return values;
189 }
190
198 PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
199 {
200 PrimaryVariables values;
201 values = exactSolution(globalPos);
202 return values;
203 }
204
205 // \}
206
210 // \{
211
222 void addPointSources(std::vector<PointSource>& pointSources) const
223 { pointSources = this->couplingManager().bulkPointSources(); }
224
243 template<class ElementVolumeVariables>
244 void pointSource(PointSource& source,
245 const Element &element,
246 const FVElementGeometry& fvGeometry,
247 const ElementVolumeVariables& elemVolVars,
248 const SubControlVolume &scv) const
249 {
250 // compute source at every integration point
251 const Scalar pressure3D = this->couplingManager().bulkPriVars(source.id())[Indices::pressureIdx];
252 const Scalar pressure1D = this->couplingManager().lowDimPriVars(source.id())[Indices::pressureIdx];
253
254 // calculate the source
255 const Scalar radius = this->couplingManager().radius(source.id());
256 const Scalar beta = 2*M_PI/(2*M_PI + std::log(radius));
257 const Scalar sourceValue = beta*(pressure1D - pressure3D);
258 source = sourceValue*source.quadratureWeight()*source.integrationElement();
259 }
260
279 template<class ElementVolumeVariables,
280 bool enable = (CouplingManager::couplingMode == EmbeddedCouplingMode::kernel),
281 std::enable_if_t<enable, int> = 0>
282 NumEqVector source(const Element &element,
283 const FVElementGeometry& fvGeometry,
284 const ElementVolumeVariables& elemVolVars,
285 const SubControlVolume &scv) const
286 {
287 static const Scalar kernelWidth = getParam<Scalar>("MixedDimension.KernelWidth");
288 static const Scalar kernelFactor = std::log(kernelWidth)-9.0/10.0;
289
290 NumEqVector source(0.0);
291
292 const auto eIdx = this->gridGeometry().elementMapper().index(element);
293 const auto& sourceIds = this->couplingManager().bulkSourceIds(eIdx);
294 const auto& sourceWeights = this->couplingManager().bulkSourceWeights(eIdx);
295
296 for (int i = 0; i < sourceIds.size(); ++i)
297 {
298 const auto id = sourceIds[i];
299 const auto weight = sourceWeights[i];
300
301 const Scalar radius = this->couplingManager().radius(id);
302 const Scalar pressure3D = this->couplingManager().bulkPriVars(id)[Indices::pressureIdx];
303 const Scalar pressure1D = this->couplingManager().lowDimPriVars(id)[Indices::pressureIdx];
304
305 // calculate the source
306 const Scalar beta = 2*M_PI/(2*M_PI + std::log(radius));
307 const Scalar sourceValue = beta*(pressure1D - pressure3D / kernelFactor * std::log(radius));
308
309 source[Indices::conti0EqIdx] += sourceValue*weight;
310 }
311
312 const auto volume = scv.volume()*elemVolVars[scv].extrusionFactor();
313 source[Indices::conti0EqIdx] /= volume;
314
315 return source;
316 }
317
319 template<class ElementVolumeVariables,
320 bool enable = (CouplingManager::couplingMode == EmbeddedCouplingMode::kernel),
321 std::enable_if_t<!enable, int> = 0>
322 NumEqVector source(const Element &element,
323 const FVElementGeometry& fvGeometry,
324 const ElementVolumeVariables& elemVolVars,
325 const SubControlVolume &scv) const
326 { return NumEqVector(0.0); }
327
330 template<class MatrixBlock, class VolumeVariables>
331 void addSourceDerivatives(MatrixBlock& block,
332 const Element& element,
333 const FVElementGeometry& fvGeometry,
334 const VolumeVariables& curElemVolVars,
335 const SubControlVolume& scv) const
336 {
337 const auto eIdx = this->gridGeometry().elementMapper().index(element);
338
339 auto key = std::make_pair(eIdx, 0);
340 if (this->pointSourceMap().count(key))
341 {
342 // call the solDependent function. Herein the user might fill/add values to the point sources
343 // we make a copy of the local point sources here
344 auto pointSources = this->pointSourceMap().at(key);
345
346 // add the point source values to the local residual (negative sign is convention for source term)
347 for (const auto& source : pointSources)
348 block[0][0] -= this->couplingManager().pointSourceDerivative(source, Dune::index_constant<0>{}, Dune::index_constant<0>{});
349 }
350 }
351
360 PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
361 { return PrimaryVariables(0.0); }
362
364 Scalar exactSolution(const GlobalPosition &globalPos) const
365 {
366 Dune::FieldVector<double, 2> xy({globalPos[0], globalPos[1]});
367 if (CouplingManager::couplingMode == EmbeddedCouplingMode::cylindersources)
368 {
369 static const auto R = getParam<Scalar>("SpatialParams.Radius");
370 if (xy.two_norm() > R)
371 return -1.0*(1+globalPos[2])/(2*M_PI)*std::log(xy.two_norm());
372 else
373 return -1.0*(1+globalPos[2])/(2*M_PI)*std::log(R);
374
375 }
376 else if (CouplingManager::couplingMode == EmbeddedCouplingMode::kernel)
377 {
378 static const auto rho = getParam<Scalar>("MixedDimension.KernelWidth");
379 const auto& r = xy.two_norm();
380 if (xy.two_norm() >= rho)
381 return -1.0*(1+globalPos[2])/(2*M_PI)*std::log(xy.two_norm());
382 else
383 return -1.0*(1+globalPos[2])/(2*M_PI)*(11.0/15.0*std::pow(r/rho, 5)
384 - 3.0/2.0*std::pow(r/rho, 4)
385 + 5.0/3.0*std::pow(r/rho, 2)
386 - 9.0/10.0 + std::log(rho));
387 }
388 else
389 return -1.0*(1+globalPos[2])/(2*M_PI)*std::log(xy.two_norm());
390 }
391
394 void computeSourceIntegral(const SolutionVector& sol, const GridVariables& gridVars)
395 {
396 PrimaryVariables source(0.0);
397 for (const auto& element : elements(this->gridGeometry().gridView()))
398 {
399 auto fvGeometry = localView(this->gridGeometry());
400 fvGeometry.bindElement(element);
401
402 auto elemVolVars = localView(gridVars.curGridVolVars());
403 elemVolVars.bindElement(element, fvGeometry, sol);
404
405 for (auto&& scv : scvs(fvGeometry))
406 {
407 auto pointSources = this->scvPointSources(element, fvGeometry, elemVolVars, scv);
408 pointSources += this->source(element, fvGeometry, elemVolVars, scv);
409 pointSources *= scv.volume()*elemVolVars[scv].extrusionFactor();
410 source += pointSources;
411 }
412 }
413
414 std::cout << "Global integrated source (3D): " << source << '\n';
415 }
416
422 template<class VtkOutputModule>
424 {
425 vtk.addField(exactPressure_, "exact pressure");
426 }
427
429 const CouplingManager& couplingManager() const
430 { return *couplingManager_; }
431
432private:
433 std::vector<Scalar> exactPressure_;
434
435 static constexpr Scalar eps_ = 1.5e-7;
436 std::string name_;
437
438 std::shared_ptr<CouplingManager> couplingManager_;
439};
440
441} // end namespace Dumux
442
443#endif
Define some often used mathematical functions.
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 tissue 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
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
Definition of a problem, for the 1p2c problem: Component transport of oxygen in interstitial fluid.
Definition: problem_tissue.hh:116
PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
Evaluates the boundary conditions for a Dirichlet boundary segment.
Definition: problem_tissue.hh:198
Scalar temperature() const
Returns the temperature within the domain [K].
Definition: problem_tissue.hh:168
void addSourceDerivatives(MatrixBlock &block, const Element &element, const FVElementGeometry &fvGeometry, const VolumeVariables &curElemVolVars, const SubControlVolume &scv) const
Definition: problem_tissue.hh:331
void addPointSources(std::vector< PointSource > &pointSources) const
Applies a vector of point sources which are possibly solution dependent.
Definition: problem_tissue.hh:222
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_tissue.hh:184
void addVtkOutputFields(VtkOutputModule &vtk) const
Adds additional VTK output data to the VTKWriter.
Definition: problem_tissue.hh:423
const CouplingManager & couplingManager() const
Get the coupling manager.
Definition: problem_tissue.hh:429
const std::string & name() const
The problem name.
Definition: problem_tissue.hh:160
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_tissue.hh:244
Scalar exactSolution(const GlobalPosition &globalPos) const
The exact solution.
Definition: problem_tissue.hh:364
TissueProblem(std::shared_ptr< const GridGeometry > gridGeometry, std::shared_ptr< CouplingManager > couplingManager)
Definition: problem_tissue.hh:136
NumEqVector source(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv) const
Evaluates the source term for all phases within a given sub control volume.
Definition: problem_tissue.hh:282
void computeSourceIntegral(const SolutionVector &sol, const GridVariables &gridVars)
Definition: problem_tissue.hh:394
PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
Evaluates the initial value for a control volume.
Definition: problem_tissue.hh:360
Definition: problem_tissue.hh:61
std::tuple< OneP > InheritsFrom
Definition: problem_tissue.hh:61
Definition: problem_tissue.hh:62
std::tuple< Tissue, CCTpfaModel > InheritsFrom
Definition: problem_tissue.hh:62
Definition: problem_tissue.hh:63
std::tuple< Tissue, BoxModel > InheritsFrom
Definition: problem_tissue.hh:63
Dune::YaspGrid< 3, Dune::EquidistantOffsetCoordinates< GetPropType< TypeTag, Properties::Scalar >, 3 > > type
Definition: problem_tissue.hh:68
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: problem_tissue.hh:95
Definition of the spatial parameters for the tissue problem.
Definition: spatialparams_tissue.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,...