Skip to content
Snippets Groups Projects
Commit f27c0966 authored by Ned Coltman's avatar Ned Coltman
Browse files

[freeflow][higherorder] improve higher order interface

parent 1f7d61e8
No related branches found
No related tags found
1 merge request!1518Feature/improve ff higherorder
......@@ -81,12 +81,13 @@ template<class TypeTag>
struct GridFluxVariablesCache<TypeTag, TTag::StaggeredModel>
{
private:
static constexpr auto enableCache = getPropValue<TypeTag, Properties::EnableGridFluxVariablesCache>();
using Problem = GetPropType<TypeTag, Properties::Problem>;
using FluxVariablesCache = GetPropType<TypeTag, Properties::FluxVariablesCache>;
using FluxVariablesCacheFiller = FluxVariablesCaching::EmptyCacheFiller;
static constexpr auto enableCache = getPropValue<TypeTag, Properties::EnableGridFluxVariablesCache>();
static constexpr auto upwindSchemeOrder = getPropValue<TypeTag, Properties::UpwindSchemeOrder>();
public:
using type = StaggeredGridFluxVariablesCache<Problem, FluxVariablesCache, FluxVariablesCacheFiller, enableCache>;
using type = StaggeredGridFluxVariablesCache<Problem, FluxVariablesCache, FluxVariablesCacheFiller, enableCache, upwindSchemeOrder>;
};
//! Set the face solution type
......
......@@ -38,14 +38,16 @@ namespace Dumux {
* \tparam P The problem type
* \tparam FVC The flux variables cache type
*/
template<class P, class FVC, class FVCF>
template<class P, class FVC, class FVCF, int upwOrder>
struct StaggeredDefaultGridFluxVariablesCacheTraits
{
using Problem = P;
using FluxVariablesCache = FVC;
using FluxVariablesCacheFiller = FVCF;
template<class GridFluxVariablesCache, bool cachingEnabled>
using LocalView = StaggeredElementFluxVariablesCache<GridFluxVariablesCache, cachingEnabled>;
static constexpr int upwindSchemeOrder = upwOrder;
};
/*!
......@@ -56,7 +58,8 @@ template<class Problem,
class FluxVariablesCache,
class FluxVariablesCacheFiller,
bool EnableGridFluxVariablesCache = false,
class Traits = StaggeredDefaultGridFluxVariablesCacheTraits<Problem, FluxVariablesCache, FluxVariablesCacheFiller> >
int upwindSchemeOrder = 1,
class Traits = StaggeredDefaultGridFluxVariablesCacheTraits<Problem, FluxVariablesCache, FluxVariablesCacheFiller, upwindSchemeOrder>>
class StaggeredGridFluxVariablesCache;
/*!
......@@ -64,17 +67,18 @@ class StaggeredGridFluxVariablesCache;
* \brief Flux variables cache class for staggered models.
Specialization in case of storing the flux cache.
*/
template<class P, class FVC, class FVCF, class Traits>
class StaggeredGridFluxVariablesCache<P, FVC, FVCF, true, Traits>
template<class P, class FVC, class FVCF, int upwindSchemeOrder, class Traits>
class StaggeredGridFluxVariablesCache<P, FVC, FVCF, true, upwindSchemeOrder, Traits>
{
using Problem = typename Traits::Problem;
using ThisType = StaggeredGridFluxVariablesCache<P, FVC, FVCF, true, Traits>;
using ThisType = StaggeredGridFluxVariablesCache<P, FVC, FVCF, true, upwindSchemeOrder, Traits>;
public:
//! export the flux variable cache type
using FluxVariablesCache = typename Traits::FluxVariablesCache;
using Scalar = typename FluxVariablesCache::Scalar;
static constexpr bool useHigherOrder = upwindSchemeOrder > 1;
//! export the flux variable cache filler type
using FluxVariablesCacheFiller = typename Traits::FluxVariablesCacheFiller;
......@@ -124,7 +128,7 @@ public:
}
//! Return the StaggeredUpwindMethods
const StaggeredUpwindMethods<Scalar>& staggeredUpwindMethods() const
const StaggeredUpwindMethods<Scalar, upwindSchemeOrder>& staggeredUpwindMethods() const
{
return staggeredUpwindMethods_;
}
......@@ -141,7 +145,7 @@ public:
private:
const Problem* problemPtr_;
StaggeredUpwindMethods<Scalar> staggeredUpwindMethods_;
StaggeredUpwindMethods<Scalar, upwindSchemeOrder> staggeredUpwindMethods_;
std::vector<FluxVariablesCache> fluxVarsCache_;
std::vector<std::size_t> globalScvfIndices_;
......@@ -152,11 +156,11 @@ private:
* \brief Flux variables cache class for staggered models.
Specialization in case of not storing the flux cache.
*/
template<class P, class FVC, class FVCF, class Traits>
class StaggeredGridFluxVariablesCache<P, FVC, FVCF, false, Traits>
template<class P, class FVC, class FVCF, int upwindSchemeOrder, class Traits>
class StaggeredGridFluxVariablesCache<P, FVC, FVCF, false, upwindSchemeOrder, Traits>
{
using Problem = typename Traits::Problem;
using ThisType = StaggeredGridFluxVariablesCache<P, FVC, FVCF, false, Traits>;
using ThisType = StaggeredGridFluxVariablesCache<P, FVC, FVCF, false, upwindSchemeOrder, Traits>;
public:
//! export the flux variable cache type
......@@ -188,14 +192,14 @@ public:
{ return *problemPtr_; }
//! Return the UpwindingMethods
const StaggeredUpwindMethods<Scalar>& staggeredUpwindMethods() const
const StaggeredUpwindMethods<Scalar, upwindSchemeOrder>& staggeredUpwindMethods() const
{
return staggeredUpwindMethods_;
}
private:
const Problem* problemPtr_;
StaggeredUpwindMethods<Scalar> staggeredUpwindMethods_;
StaggeredUpwindMethods<Scalar, upwindSchemeOrder> staggeredUpwindMethods_;
};
......
......@@ -47,7 +47,7 @@ enum class DifferencingScheme
/**
* \brief This file contains different higher order methods for approximating the velocity.
*/
template<class Scalar>
template<class Scalar, int upwindSchemeOrder>
class StaggeredUpwindMethods
{
public:
......@@ -55,11 +55,11 @@ public:
{
upwindWeight_ = getParamFromGroup<Scalar>(paramGroup, "Flux.UpwindWeight");
if (hasParamInGroup(paramGroup, "Flux.TvdApproach"))
if (upwindSchemeOrder > 1)
{
// Read the runtime parameters
tvdApproach_ = tvdApproachFromString(getParamFromGroup<std::string>(paramGroup, "Flux.TvdApproach"));
differencingScheme_ = differencingSchemeFromString(getParamFromGroup<std::string>(paramGroup, "Flux.DifferencingScheme"));
tvdApproach_ = tvdApproachFromString(getParamFromGroup<std::string>(paramGroup, "Flux.TvdApproach", "Uniform"));
differencingScheme_ = differencingSchemeFromString(getParamFromGroup<std::string>(paramGroup, "Flux.DifferencingScheme", "Minmod"));
// Assign the limiter_ depending on the differencing scheme
switch (differencingScheme_)
......@@ -109,6 +109,33 @@ public:
break;
}
}
if (!hasParamInGroup(paramGroup, "Flux.TvdApproach"))
{
std::cout << "No TvdApproach specified. Defaulting to the Uniform method." << "\n";
std::cout << "Other available TVD approaches for uniform (and nonuniform) grids are as follows: \n"
<< " " << tvdApproachToString(TvdApproach::uniform) << ": assumes a Uniform cell size distribution\n"
<< " " << tvdApproachToString(TvdApproach::li) << ": Li's approach for nonuniform cell sizes\n"
<< " " << tvdApproachToString(TvdApproach::hou) << ": Hou's approach for nonuniform cell sizes \n";
std::cout << "Each approach can be specified as written above in the Flux group under the title TvdApproach in your input file. \n";
}
if (!hasParamInGroup(paramGroup, "Flux.DifferencingScheme"))
{
std::cout << "No DifferencingScheme specified. Defaulting to the Minmod scheme." << "\n";
std::cout << "Other available Differencing Schemes are as follows: \n"
<< " " << differencingSchemeToString(DifferencingScheme::vanleer) << ": The Vanleer flux limiter\n"
<< " " << differencingSchemeToString(DifferencingScheme::vanalbada) << ": The Vanalbada flux limiter\n"
<< " " << differencingSchemeToString(DifferencingScheme::minmod) << ": The Minmod flux limiter\n"
<< " " << differencingSchemeToString(DifferencingScheme::superbee) << ": The Superbee flux limiter\n"
<< " " << differencingSchemeToString(DifferencingScheme::umist) << ": The Umist flux limiter\n"
<< " " << differencingSchemeToString(DifferencingScheme::mclimiter) << ": The Mclimiter flux limiter\n"
<< " " << differencingSchemeToString(DifferencingScheme::wahyd) << ": The Wahyd flux limiter";
std::cout << "Each scheme can be specified as written above in the Flux group under the variable DifferencingScheme in your input file. \n";
}
std::cout << "Using the tvdApproach \"" << tvdApproachToString(tvdApproach_)
<< "\" and the differencing Scheme \" " << differencingSchemeToString(differencingScheme_) << "\" \n";
}
else
{
......@@ -116,8 +143,6 @@ public:
tvdApproach_ = TvdApproach::none;
differencingScheme_ = DifferencingScheme::none;
}
std::cout << "Using the tvdApproach " << tvdApproachToString(tvdApproach_)
<< " and the differencing Scheme " << differencingSchemeToString(differencingScheme_) << "\n";
}
/**
......
......@@ -19,7 +19,3 @@ MaxRelativeShift = 1e-5
[Vtk]
WriteFaceData = false
[Flux]
TvdApproach = Hou
DifferencingScheme = Vanleer
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