Skip to content
GitLab
Explore
Sign in
Register
Primary navigation
Search or go to…
Project
dumux
Manage
Activity
Members
Labels
Plan
Issues
Issue boards
Milestones
Wiki
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Build
Pipelines
Jobs
Pipeline schedules
Artifacts
Deploy
Releases
Container Registry
Model registry
Operate
Environments
Monitor
Incidents
Analyze
Value stream analytics
Contributor analytics
CI/CD analytics
Repository analytics
Model experiments
Help
Help
Support
GitLab documentation
Compare GitLab plans
Community forum
Contribute to GitLab
Provide feedback
Keyboard shortcuts
?
Snippets
Groups
Projects
Show more breadcrumbs
dumux-repositories
dumux
Commits
0d4ecddb
Commit
0d4ecddb
authored
6 years ago
by
Ned Coltman
Committed by
Timo Koch
6 years ago
Browse files
Options
Downloads
Patches
Plain Diff
[box][darcy] Include a class to inherit from, without TypeTags
parent
411d6476
No related branches found
No related tags found
1 merge request
!929
Free the Box Method's Darcy's law from TypeTag
Changes
1
Hide whitespace changes
Inline
Side-by-side
Showing
1 changed file
dumux/discretization/box/darcyslaw.hh
+25
-14
25 additions, 14 deletions
dumux/discretization/box/darcyslaw.hh
with
25 additions
and
14 deletions
dumux/discretization/box/darcyslaw.hh
+
25
−
14
View file @
0d4ecddb
...
...
@@ -30,30 +30,39 @@
#include
<dumux/common/properties.hh>
#include
<dumux/discretization/methods.hh>
namespace
Dumux
{
namespace
Dumux
{
// forward declaration
template
<
class
TypeTag
,
DiscretizationMethod
discMethod
>
class
DarcysLawImplementation
;
// forward declaration
template
<
class
Scalar
,
class
FVGridGeometry
>
class
BoxDarcysLaw
;
/*!
* \ingroup DarcysLaw
* \brief Specialization of Darcy's Law for the box method.
*/
template
<
class
TypeTag
>
template
<
class
TypeTag
>
class
DarcysLawImplementation
<
TypeTag
,
DiscretizationMethod
::
box
>
:
public
BoxDarcysLaw
<
typename
GET_PROP_TYPE
(
TypeTag
,
Scalar
),
typename
GET_PROP_TYPE
(
TypeTag
,
FVGridGeometry
)
>
{
};
/*!
* \ingroup BoxDiscretization
* \brief Darcy's law for box schemes
* \tparam Scalar the scalar type for scalar physical quantities
* \tparam FVGridGeometry the grid geometry
*/
template
<
class
Scalar
,
class
FVGridGeometry
>
class
BoxDarcysLaw
{
using
Problem
=
typename
GET_PROP_TYPE
(
TypeTag
,
Problem
);
using
FVElementGeometry
=
typename
GET_PROP_TYPE
(
TypeTag
,
FVGridGeometry
)
::
LocalView
;
using
FVElementGeometry
=
typename
FVGridGeometry
::
LocalView
;
using
SubControlVolume
=
typename
FVElementGeometry
::
SubControlVolume
;
using
SubControlVolumeFace
=
typename
FVElementGeometry
::
SubControlVolumeFace
;
using
ElemFluxVarCache
=
typename
GET_PROP_TYPE
(
TypeTag
,
GridFluxVariablesCache
)
::
LocalView
;
using
FluxVarCache
=
typename
GET_PROP_TYPE
(
TypeTag
,
FluxVariablesCache
);
using
ElementVolumeVariables
=
typename
GET_PROP_TYPE
(
TypeTag
,
GridVolumeVariables
)
::
LocalView
;
using
GridView
=
typename
GET_PROP_TYPE
(
TypeTag
,
GridView
);
using
Scalar
=
typename
GET_PROP_TYPE
(
TypeTag
,
Scalar
);
using
GridView
=
typename
FVGridGeometry
::
GridView
;
using
Element
=
typename
GridView
::
template
Codim
<
0
>
::
Entity
;
using
IndexType
=
typename
GridView
::
IndexSet
::
IndexType
;
using
CoordScalar
=
typename
GridView
::
ctype
;
enum
{
dim
=
GridView
::
dimension
};
...
...
@@ -64,13 +73,14 @@ class DarcysLawImplementation<TypeTag, DiscretizationMethod::box>
public
:
template
<
class
Problem
,
class
ElementVolumeVariables
,
class
ElementFluxVarsCache
>
static
Scalar
flux
(
const
Problem
&
problem
,
const
Element
&
element
,
const
FVElementGeometry
&
fvGeometry
,
const
ElementVolumeVariables
&
elemVolVars
,
const
SubControlVolumeFace
&
scvf
,
const
IndexType
phaseIdx
,
const
ElemFluxVarCache
&
elemFluxVarCache
)
const
int
phaseIdx
,
const
Elem
ent
FluxVar
s
Cache
&
elemFluxVarCache
)
{
const
auto
&
fluxVarCache
=
elemFluxVarCache
[
scvf
];
const
auto
&
insideScv
=
fvGeometry
.
scv
(
scvf
.
insideScvIdx
());
...
...
@@ -86,7 +96,7 @@ public:
outsideK
*=
outsideVolVars
.
extrusionFactor
();
const
auto
K
=
problem
.
spatialParams
().
harmonicMean
(
insideK
,
outsideK
,
scvf
.
unitOuterNormal
());
static
const
bool
enableGravity
=
getParamFromGroup
<
bool
>
(
GET_PROP_VALUE
(
TypeTag
,
ModelParameter
Group
),
"Problem.EnableGravity"
);
static
const
bool
enableGravity
=
getParamFromGroup
<
bool
>
(
problem
.
param
Group
(
),
"Problem.EnableGravity"
);
const
auto
&
shapeValues
=
fluxVarCache
.
shapeValues
();
...
...
@@ -112,6 +122,7 @@ public:
}
// compute transmissibilities ti for analytical jacobians
template
<
class
Problem
,
class
ElementVolumeVariables
,
class
FluxVarCache
>
static
std
::
vector
<
Scalar
>
calculateTransmissibilities
(
const
Problem
&
problem
,
const
Element
&
element
,
const
FVElementGeometry
&
fvGeometry
,
...
...
This diff is collapsed.
Click to expand it.
Preview
0%
Loading
Try again
or
attach a new file
.
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Save comment
Cancel
Please
register
or
sign in
to comment