Skip to content
Snippets Groups Projects
Commit d56ead20 authored by Timo Koch's avatar Timo Koch
Browse files

[cvfe] Introduce Navier-Stokes CVFEVolumeVariables to avoid code duplication

parent 2f5363a4
No related branches found
No related tags found
1 merge request!3349Feature/cvfe improvements
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*****************************************************************************
* See the file COPYING for full copying permissions. *
* *
* This program is free software: you can redistribute it and/or modify *
* it under the terms of the GNU General Public License as published by *
* the Free Software Foundation, either version 3 of the License, or *
* (at your option) any later version. *
* *
* This program is distributed in the hope that it will be useful, *
* but WITHOUT ANY WARRANTY; without even the implied warranty of *
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
* GNU General Public License for more details. *
* *
* You should have received a copy of the GNU General Public License *
* along with this program. If not, see <http://www.gnu.org/licenses/>. *
*****************************************************************************/
/*!
* \file
* \ingroup NavierStokesModel
*
* \copydoc Dumux::NavierStokesVolumeVariables
*/
#ifndef DUMUX_NAVIERSTOKES_MOMENTUM_CVFE_VOLUME_VARIABLES_HH
#define DUMUX_NAVIERSTOKES_MOMENTUM_CVFE_VOLUME_VARIABLES_HH
namespace Dumux {
/*!
* \ingroup NavierStokesModel
* \brief Volume variables for the single-phase Navier-Stokes model.
*/
template <class Traits>
class NavierStokesMomentumCVFEVolumeVariables
{
using Scalar = typename Traits::PrimaryVariables::value_type;
static_assert(Traits::PrimaryVariables::dimension == Traits::ModelTraits::dim());
public:
//! export the type used for the primary variables
using PrimaryVariables = typename Traits::PrimaryVariables;
//! export the indices type
using Indices = typename Traits::ModelTraits::Indices;
/*!
* \brief Update all quantities for a given control volume
*
* \param elemSol A vector containing all primary variables connected to the element
* \param problem The object specifying the problem which ought to
* be simulated
* \param element An element which contains part of the control volume
* \param scv The sub-control volume
*/
template<class ElementSolution, class Problem, class Element, class SubControlVolume>
void update(const ElementSolution& elemSol,
const Problem& problem,
const Element& element,
const SubControlVolume& scv)
{
priVars_ = elemSol[scv.indexInElement()];
extrusionFactor_ = problem.spatialParams().extrusionFactor(element, scv, elemSol);
}
/*!
* \brief Return how much the sub-control volume is extruded.
*
* This means the factor by which a lower-dimensional (1D or 2D)
* entity needs to be expanded to get a full dimensional cell. The
* default is 1.0 which means that 1D problems are actually
* thought as pipes with a cross section of 1 m^2 and 2D problems
* are assumed to extend 1 m to the back.
*/
Scalar extrusionFactor() const
{ return extrusionFactor_; }
PrimaryVariables velocity() const
{ return priVars_; }
Scalar velocity(const int dirIdx) const
{ return priVars_[dirIdx]; }
/*!
* \brief Return a component of primary variable vector
* \param pvIdx The index of the primary variable of interest
*/
Scalar priVar(const int pvIdx) const
{ return priVars_[pvIdx]; }
/*!
* \brief Return the primary variable vector
*/
const PrimaryVariables& priVars() const
{ return priVars_; }
private:
PrimaryVariables priVars_;
Scalar extrusionFactor_;
};
} // end namespace Dumux
#endif
......@@ -57,8 +57,8 @@
#include <dumux/flux/fluxvariablescaching.hh>
#include <dumux/freeflow/navierstokes/momentum/cvfe/localresidual.hh>
#include <dumux/freeflow/navierstokes/momentum/cvfe/volumevariables.hh>
#include "volumevariables.hh"
#include "indices.hh"
namespace Dumux {
......@@ -191,7 +191,7 @@ private:
using Traits = NavierStokesMomentumDiamondVolumeVariablesTraits<PV, FSY, FST, MT>;
public:
using type = NavierStokesMomentumDiamondVolumeVariables<Traits>;
using type = NavierStokesMomentumCVFEVolumeVariables<Traits>;
};
// This is the default (model not coupled with a mass (pressure) discretization)
......
......@@ -25,81 +25,14 @@
#ifndef DUMUX_NAVIERSTOKES_MOMENTUM_DIAMOND_VOLUME_VARIABLES_HH
#define DUMUX_NAVIERSTOKES_MOMENTUM_DIAMOND_VOLUME_VARIABLES_HH
#warning "This header is deprecated and will be removed after 3.7"
namespace Dumux {
/*!
* \ingroup NavierStokesModel
* \brief Volume variables for the single-phase Navier-Stokes model.
*/
template <class Traits>
class NavierStokesMomentumDiamondVolumeVariables
{
using Scalar = typename Traits::PrimaryVariables::value_type;
static_assert(Traits::PrimaryVariables::dimension == Traits::ModelTraits::dim());
public:
//! export the type used for the primary variables
using PrimaryVariables = typename Traits::PrimaryVariables;
#include <dumux/freeflow/navierstokes/momentum/cvfe/volumevariables.hh>
//! export the indices type
using Indices = typename Traits::ModelTraits::Indices;
/*!
* \brief Update all quantities for a given control volume
*
* \param elemSol A vector containing all primary variables connected to the element
* \param problem The object specifying the problem which ought to
* be simulated
* \param element An element which contains part of the control volume
* \param scv The sub-control volume
*/
template<class ElementSolution, class Problem, class Element, class SubControlVolume>
void update(const ElementSolution& elemSol,
const Problem& problem,
const Element& element,
const SubControlVolume& scv)
{
priVars_ = elemSol[scv.indexInElement()];
extrusionFactor_ = problem.spatialParams().extrusionFactor(element, scv, elemSol);
}
/*!
* \brief Return how much the sub-control volume is extruded.
*
* This means the factor by which a lower-dimensional (1D or 2D)
* entity needs to be expanded to get a full dimensional cell. The
* default is 1.0 which means that 1D problems are actually
* thought as pipes with a cross section of 1 m^2 and 2D problems
* are assumed to extend 1 m to the back.
*/
Scalar extrusionFactor() const
{ return extrusionFactor_; }
PrimaryVariables velocity() const
{ return priVars_; }
Scalar velocity(const int dirIdx) const
{ return priVars_[dirIdx]; }
/*!
* \brief Return a component of primary variable vector
* \param pvIdx The index of the primary variable of interest
*/
Scalar priVar(const int pvIdx) const
{ return priVars_[pvIdx]; }
/*!
* \brief Return the primary variable vector
*/
const PrimaryVariables& priVars() const
{ return priVars_; }
namespace Dumux {
private:
PrimaryVariables priVars_;
Scalar extrusionFactor_;
};
template<class TypeTag>
using NavierStokesMomentumDiamondVolumeVariables [[deprecated("Use NavierStokesMomentumCVFEVolumeVariables")]] = NavierStokesMomentumCVFEVolumeVariables<TypeTag>;
} // end namespace Dumux
......
......@@ -57,8 +57,8 @@
#include <dumux/flux/fluxvariablescaching.hh>
#include <dumux/freeflow/navierstokes/momentum/cvfe/localresidual.hh>
#include <dumux/freeflow/navierstokes/momentum/cvfe/volumevariables.hh>
#include "volumevariables.hh"
#include "indices.hh"
namespace Dumux {
......@@ -191,7 +191,7 @@ private:
using Traits = NavierStokesMomentumPQ1BubbleVolumeVariablesTraits<PV, FSY, FST, MT>;
public:
using type = NavierStokesMomentumPQ1BubbleVolumeVariables<Traits>;
using type = NavierStokesMomentumCVFEVolumeVariables<Traits>;
};
// This is the default (model not coupled with a mass (pressure) discretization)
......
......@@ -25,7 +25,9 @@
#ifndef DUMUX_NAVIERSTOKES_MOMENTUM_PQ1BUBBLE_VOLUME_VARIABLES_HH
#define DUMUX_NAVIERSTOKES_MOMENTUM_PQ1BUBBLE_VOLUME_VARIABLES_HH
#include <dumux/freeflow/navierstokes/momentum/diamond/volumevariables.hh>
#warning "This header is deprecated and will be removed after 3.7"
#include <dumux/freeflow/navierstokes/momentum/cvfe/volumevariables.hh>
namespace Dumux {
......@@ -34,7 +36,7 @@ namespace Dumux {
* \brief Volume variables for the single-phase Navier-Stokes model.
*/
template <class Traits>
using NavierStokesMomentumPQ1BubbleVolumeVariables = NavierStokesMomentumDiamondVolumeVariables<Traits>;
using NavierStokesMomentumPQ1BubbleVolumeVariables [[deprecated("Use NavierStokesMomentumCVFEVolumeVariables")]] = NavierStokesMomentumCVFEVolumeVariables<Traits>;
} // end namespace Dumux
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment