CFV Module

Compact Finite Volume schemes and variational reconstruction.

Namespace Overview

namespace CFV

Typedefs

using tPy_FiniteVolume = py_class_ssp<FiniteVolume>
using tPy_ModelEvaluator = py_class_ssp<ModelEvaluator>
template<int dim = 2>
using tPy_VariationalReconstruction = py_class_ssp<VariationalReconstruction<dim>>
using tCoeffPair = DNDS::ArrayPair<DNDS::ArrayEigenVector<NonUniformSize>>
using tCoeff = decltype(tCoeffPair::father)
using t3VecsPair = DNDS::ArrayPair<DNDS::ArrayEigenUniMatrixBatch<3, 1>>
using t3Vecs = decltype(t3VecsPair::father)
using t3VecPair = Geom::tCoordPair
using t3Vec = Geom::tCoord
using t3MatPair = DNDS::ArrayPair<DNDS::ArrayEigenMatrix<3, 3>>
using t3Mat = decltype(t3MatPair::father)
using tVVecPair = ::DNDS::ArrayPair<DNDS::ArrayEigenVector<DynamicSize>>
using tVVec = decltype(tVVecPair::father)
using tMatsPair = DNDS::ArrayPair<DNDS::ArrayEigenUniMatrixBatch<DynamicSize, DynamicSize>>
using tMats = decltype(tMatsPair::father)
using tVecsPair = DNDS::ArrayPair<DNDS::ArrayEigenUniMatrixBatch<DynamicSize, 1>>
using tVecs = decltype(tVecsPair::father)
using tVMatPair = DNDS::ArrayPair<DNDS::ArrayEigenMatrix<DynamicSize, DynamicSize>>
using tVMat = decltype(tVMatPair::father)
using tRecAtrPair = DNDS::ArrayPair<DNDS::ParArray<RecAtr, 1>>
using tRecAtr = decltype(tRecAtrPair::father)
template<int nVarsFixed>
using tURec = ArrayDof<DynamicSize, nVarsFixed>
template<int nVarsFixed>
using tUDof = ArrayDof<nVarsFixed, 1>
template<int nVarsFixed, int gDim>
using tUGrad = ArrayDof<gDim, nVarsFixed>
using tScalarPair = DNDS::ArrayPair<DNDS::ArrayEigenVector<1>>
using tScalar = decltype(tScalarPair::father)
using tPy_RecAtr = py::classh<RecAtr>

Functions

template<>
void finiteVolumeCellOpTest_run<B>(FiniteVolume::t_deviceView<B> &fv, tUDof<DynamicSize>::t_deviceView<B> &u, tUGrad<DynamicSize, 3>::t_deviceView<B> &u_grad, int nIter, const t_jsonconfig &settings)
template<DeviceBackend B, bool iVarOne = false> DNDS_DEVICE_CALLABLE void finiteVolumeCellOpTest (FiniteVolume::t_deviceView< B > &fv, tUDof< DynamicSize >::t_deviceView< B > &u, tUGrad< DynamicSize, 3 >::t_deviceView< B > &u_grad, index iCell, real *local_buf, int iVar=UnInitRowsize)
template<DeviceBackend B>
void finiteVolumeCellOpTest_run(FiniteVolume::t_deviceView<B> &fv, tUDof<DynamicSize>::t_deviceView<B> &u, tUGrad<DynamicSize, 3>::t_deviceView<B> &u_grad, int nIter, const t_jsonconfig &settings)
template<DeviceBackend B>
void finiteVolumeCellOpTest_main(FiniteVolume &fv, tUDof<DynamicSize> &u, tUGrad<DynamicSize, 3> &u_grad, int nIter, const t_jsonconfig &settings)
template<DeviceBackend B, int nVarsFixed, bool iVarOne = false> DNDS_DEVICE_CALLABLE void finiteVolumeCellOpTest_Fixed (FiniteVolume::t_deviceView< B > &fv, typename tUDof< nVarsFixed >::template t_deviceView< B > &u, typename tUGrad< nVarsFixed, 3 >::template t_deviceView< B > &u_grad, index iCell, int iVar=UnInitRowsize)
template<DeviceBackend B, int nVarsFixed>
void finiteVolumeCellOpTest_Fixed_main(FiniteVolume &fv, tUDof<nVarsFixed> &u, tUGrad<nVarsFixed, 3> &u_grad, int nIter, const t_jsonconfig &settings)
template<DeviceBackend B, int nVarsFixed, bool iVarOne = false> DNDS_DEVICE_CALLABLE void finiteVolumeCellOpTest_SOA_ver0 (FiniteVolume::t_deviceView< B > &fv, std::array< tUDof< 1 >::t_deviceView< B >, nVarsFixed > &u, std::array< tUGrad< 1, 3 >::t_deviceView< B >, nVarsFixed > &u_grad, index iCell, int iVar=UnInitRowsize)
template<DeviceBackend B, int nVarsFixed>
void finiteVolumeCellOpTest_SOA_ver0_main(FiniteVolume &fv, std::array<tUDof<1>, nVarsFixed> &u, std::array<tUGrad<1, 3>, nVarsFixed> &u_grad, int nIter, const t_jsonconfig &settings)
template<DeviceBackend B>
void pybind11_BenchmarkFiniteVolume_define(py::module_ &m)
template<DeviceBackend B, int nVarsFixed>
void pybind11_BenchmarkFiniteVolume_define_Fixed(py::module_ &m)
template<DeviceBackend B, int nVarsFixed>
void pybind11_BenchmarkFiniteVolume_define_SOA_ver0(py::module_ &m)
template<int nVarsFixed = 1>
void BuildUDofOnMesh(tUDof<nVarsFixed> &u, const std::string &name, const MPIInfo &mpi, const ssp<Geom::UnstructuredMesh> &mesh, int nVars, bool buildSon = true, bool buildTrans = true, Geom::MeshLoc varloc = Geom::MeshLoc::Cell)
template<int nVarsFixed, int dim>
void BuildUGradDOnMesh(tUGrad<nVarsFixed, dim> &u, const std::string &name, const MPIInfo &mpi, const ssp<Geom::UnstructuredMesh> &mesh, int nVars, bool buildSon = true, bool buildTrans = true, Geom::MeshLoc varloc = Geom::MeshLoc::Cell)
void pybind11_FiniteVolume_define(py::module_ &m)
template<int dim, typename ThetaArray>
inline Eigen::ArrayXd PolynomialSquaredNorm(const ThetaArray &theta)
template<int dim, typename ThetaArray1, typename ThetaArray2>
inline Eigen::ArrayXd PolynomialDotProduct(const ThetaArray1 &theta1, const ThetaArray2 &theta2)
template<typename TinOthers, typename Tout>
static inline void FWBAP_L2_Multiway_Polynomial2D(const TinOthers &uOthers, int Nother, Tout &uOut, real n1 = 1)
template<typename Tcenter, typename TinOthers, typename Tout>
static inline void FMEMM_Multiway_Polynomial2D(const Tcenter &u, const TinOthers &uOthers, int Nother, Tout &uOut, real n1 = 1)
template<typename TinOthers, typename Tout>
static inline void FWBAP_L2_Multiway_PolynomialOrth(const TinOthers &uOthers, int Nother, Tout &uOut, real n1 = 1)
template<typename TinOthers, typename Tout>
inline void FWBAP_L2_Multiway(const TinOthers &uOthers, int Nother, Tout &uOut, real n1 = 1.0)
template<typename Tin1, typename Tin2, typename Tout>
inline void FWBAP_L2_Biway(const Tin1 &u1, const Tin2 &u2, Tout &uOut, real n)
template<typename Tin1, typename Tin2, typename Tout>
inline void FWBAP_L2_Cut_Biway(const Tin1 &u1, const Tin2 &u2, Tout &uOut, real n)
template<typename Tin1, typename Tin2, typename Tout>
inline void FMINMOD_Biway(const Tin1 &u1, const Tin2 &u2, Tout &uOut, real n)
template<typename Tin1, typename Tin2, typename Tout>
inline void FVanLeer_Biway(const Tin1 &u1, const Tin2 &u2, Tout &uOut, real n)
template<int dim, int nVarsFixed, typename Tin1, typename Tin2, typename Tout>
inline void FWBAP_L2_Biway_PolynomialNorm(const Tin1 &u1, const Tin2 &u2, Tout &uOut, real n)
template<int dim, int nVarsFixed, typename Tin1, typename Tin2, typename Tout>
inline void FMEMM_Biway_PolynomialNorm(const Tin1 &u1, const Tin2 &u2, Tout &uOut, real n)
template<typename Tin1, typename Tin2, typename Tout>
inline void FWBAP_L2_Biway_PolynomialOrth(const Tin1 &u1, const Tin2 &u2, Tout &uOut, real n)
void pybind11_ModelEvaluator_define(py::module_ &m)
template<int dim, int powV = 1, class TDiffI, class TDiffJ>
inline void AccumulateDiffOrderContributions(const Eigen::MatrixBase<TDiffI> &DiffI, const Eigen::MatrixBase<TDiffJ> &DiffJ, MatrixXR &Conj, const Eigen::Vector<real, Eigen::Dynamic> &wgd, int cnDiffs, real faceL)
template<int dim = 2>
void pybind11_VariationalReconstruction_define(py::module_ &m)
template<int dim, int nVarsFixed, typename Tin1, typename Tin2, typename Tout>
inline void DispatchBiwayLimiter(int limiterBiwayAlter, const Tin1 &uThis, const Tin2 &uOther, Tout &uOut, real n)
DNDS_SWITCH_INTELLISENSE (template< int dim > template< int nVarsFixed >, template<> template<>) auto VariationalReconstruction< dim >
template<int dim>
inline std::pair<int, int> GetRecDOFRange(int pOrder)
template<int nVarsFixed>
inline void pybind11_define_tURec_alias(py::module_ &m, py::module_ &m_dnds)
template<int nVarsFixed>
inline void pybind11_define_tUDof_alias(py::module_ &m, py::module_ &m_dnds)
template<int nVarsFixed>
inline void pybind11_define_tUGrad_alias(py::module_ &m, py::module_ &m_dnds)
inline void pybind11_define_tVVecPair_alias(py::module_ &m, py::module_ &m_dnds)
inline void pybind11_define_tMatsPair_alias(py::module_ &m, py::module_ &m_dnds)
inline void pybind11_define_RecAtr(py::module_ &m)
inline void pybind11_define_RecAtrArrayPair_and_alias(py::module_ &m)
inline void pybind11_VRDefines_define(py::module_ &m, py::module_ &m_dnds)
DNDS_DEFINE_ENUM_JSON (VRSettings::FunctionalSettings::ScaleType, {{VRSettings::FunctionalSettings::ScaleType::UnknownScale, nullptr}, {VRSettings::FunctionalSettings::ScaleType::MeanAACBB, "MeanAACBB"}, {VRSettings::FunctionalSettings::ScaleType::BaryDiff, "BaryDiff"}, {VRSettings::FunctionalSettings::ScaleType::CellMax, "CellMax"}}) DNDS_DEFINE_ENUM_JSON(VRSettings
DNDS_DEFINE_ENUM_JSON (VRSettings::FunctionalSettings::GeomWeightScheme, {{VRSettings::FunctionalSettings::GeomWeightScheme::UnknownGeomWeight, nullptr}, {VRSettings::FunctionalSettings::GeomWeightScheme::GWNone, "GWNone"}, {VRSettings::FunctionalSettings::GeomWeightScheme::HQM_SD, "HQM_SD"}, {VRSettings::FunctionalSettings::GeomWeightScheme::SD_Power, "SD_Power"}}) DNDS_DEFINE_ENUM_JSON(VRSettings

Variables

static const DeviceBackend B = DeviceBackend::Host
static const int dim = 3
static const int nVarsFixed = 6
class FiniteVolume : public DNDS::DeviceTransferable<FiniteVolume>

Subclassed by DNDS::CFV::VariationalReconstruction< dim >

Public Functions

template<class TArrayPair, class ...TOthers>
inline void MakePairDefaultOnCell(TArrayPair &aPair, const std::string &name, TOthers... others)
template<class TArrayPair, class ...TOthers>
inline void MakePairDefaultOnFace(TArrayPair &aPair, const std::string &name, TOthers... others)
template<DeviceBackend B, int nVarsFixed>
struct finiteVolumeCellOpTest_Fixed_entry
template<int nVarsFixed>
struct finiteVolumeCellOpTest_Fixed_entry<DeviceBackend::Host, nVarsFixed>
template<DeviceBackend B, int nVarsFixed>
struct finiteVolumeCellOpTest_SOA_ver0_entry
template<int nVarsFixed>
struct finiteVolumeCellOpTest_SOA_ver0_entry<DeviceBackend::Host, nVarsFixed>
template<DeviceBackend B>
class FiniteVolumeDeviceView
struct FiniteVolumeSettings

Subclassed by DNDS::CFV::VRSettings

Public Functions

inline DNDS_HOST void WriteIntoJson (json &jsonSetting) const
inline DNDS_HOST void ParseFromJson (const json &jsonSetting)

Public Members

int intOrder = {1}
bool ignoreMeshGeometryDeficiency = false
class ModelEvaluator

Public Functions

inline TU generateBoundaryValue(TU &ULxy, const TU &ULMeanXy, index iCell, index iFace, int iG, const TVec &uNorm, const TMat &normBase, const Geom::tPoint &pPhysics, real t, Geom::t_index btype, bool fixUL = false, int geomMode = 0, int linMode = 0)
struct EvaluateRHSOptions
struct ModelSettings
struct RecAtr
template<int dim = 2>
class VariationalReconstruction : public DNDS::CFV::FiniteVolume

Public Functions

inline void parseSettings(VRSettings::json &j)
template<typename ...Ts>
inline void ApplyPeriodicTransform(int if2c, Geom::t_index faceID, Ts&... data) const
template<class TOut>
inline void FDiffBaseValue(TOut &DiBj, const Geom::tPoint &pPhy, index iCell, index iFace, int iG, int flag = 0)
template<class TList>
inline MatrixXR GetIntPointDiffBaseValue(index iCell, index iFace, rowsize if2c, int iG, TList &&diffList = EigenAll, uint8_t maxDiff = UINT8_MAX)
inline auto GetMatrixSecondary(index iCell, index iFace, index if2c)
template<int nVarsFixed = 5>
void DoReconstruction2ndGrad(tUGrad<nVarsFixed, dim> &uRec, tUDof<nVarsFixed> &u, const TFBoundary<nVarsFixed> &FBoundary, int method)
template<int nVarsFixed = 5>
void DoReconstruction2nd(tURec<nVarsFixed> &uRec, tUDof<nVarsFixed> &u, const TFBoundary<nVarsFixed> &FBoundary, int method, const std::vector<int> &mask)
template<int nVarsFixed = 5>
void DoReconstructionIter(tURec<nVarsFixed> &uRec, tURec<nVarsFixed> &uRecNew, tUDof<nVarsFixed> &u, const TFBoundary<nVarsFixed> &FBoundary, bool putIntoNew = false, bool recordInc = false, bool uRecIsZero = false)
template<int nVarsFixed = 5>
void DoReconstructionIterDiff(tURec<nVarsFixed> &uRec, tURec<nVarsFixed> &uRecDiff, tURec<nVarsFixed> &uRecNew, tUDof<nVarsFixed> &u, const TFBoundaryDiff<nVarsFixed> &FBoundaryDiff)
template<int nVarsFixed = 5>
void DoReconstructionIterSOR(tURec<nVarsFixed> &uRec, tURec<nVarsFixed> &uRecInc, tURec<nVarsFixed> &uRecNew, tUDof<nVarsFixed> &u, const TFBoundaryDiff<nVarsFixed> &FBoundaryDiff, bool reverse = false)
template<int nVarsFixed>
void DoLimiterWBAP_C(tUDof<nVarsFixed> &u, tURec<nVarsFixed> &uRec, tURec<nVarsFixed> &uRecNew, tURec<nVarsFixed> &uRecBuf, tScalarPair &si, bool ifAll, const tFMEig<nVarsFixed> &FM, const tFMEig<nVarsFixed> &FMI, bool putIntoNew = false)
template<int nVarsFixed>
void DoLimiterWBAP_3(tUDof<nVarsFixed> &u, tURec<nVarsFixed> &uRec, tURec<nVarsFixed> &uRecNew, tURec<nVarsFixed> &uRecBuf, tScalarPair &si, bool ifAll, const tFMEig<nVarsFixed> &FM, const tFMEig<nVarsFixed> &FMI, bool putIntoNew = false)
struct VRSettings : public DNDS::CFV::FiniteVolumeSettings

Public Functions

inline void WriteIntoJson(json &jsonSetting) const
inline void ParseFromJson(const json &jsonSetting)

Public Members

int intOrderVRBC = {-2}
bool cacheDiffBase = false
uint8_t cacheDiffBaseSize = UINT8_MAX
bool SORInstead = true
real smoothThreshold = 0.01
real WBAP_nStd = 10
bool normWBAP = false
int limiterBiwayAlter = 0
int subs2ndOrder = 0
int subs2ndOrderGGScheme = 1
real svdTolerance = 0
real bcWeight = 1
struct BaseSettings
struct FunctionalSettings
namespace _Limiters_Internal

Functions

inline real W12n1(real u1, real u2)
inline real W12center(real *u, const int J, real n)

Key Classes

template<int dim = 2>
class VariationalReconstruction : public DNDS::CFV::FiniteVolume

Inheritance diagram for DNDS::CFV::VariationalReconstruction:

digraph {
    graph [bgcolor="#00000000"]
    node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2]
    edge [color="#1414CE"]
    "2" [label="DNDS::CFV::FiniteVolume" tooltip="DNDS::CFV::FiniteVolume"]
    "1" [label="DNDS::CFV::VariationalReconstruction< dim >" tooltip="DNDS::CFV::VariationalReconstruction< dim >" fillcolor="#BFBFBF"]
    "3" [label="DNDS::DeviceTransferable< FiniteVolume >" tooltip="DNDS::DeviceTransferable< FiniteVolume >"]
    "2" -> "3" [dir=forward tooltip="public-inheritance"]
    "1" -> "2" [dir=forward tooltip="public-inheritance"]
}

Collaboration diagram for DNDS::CFV::VariationalReconstruction:

digraph {
    graph [bgcolor="#00000000"]
    node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2]
    edge [color="#1414CE"]
    "16" [label="DNDS::ArrayPair< DNDS::ArrayEigenMatrix< 3, 3 > >" tooltip="DNDS::ArrayPair< DNDS::ArrayEigenMatrix< 3, 3 > >"]
    "12" [label="DNDS::ArrayPair< DNDS::ArrayEigenUniMatrixBatch< 3, 1 > >" tooltip="DNDS::ArrayPair< DNDS::ArrayEigenUniMatrixBatch< 3, 1 > >"]
    "6" [label="DNDS::ArrayPair< DNDS::ArrayEigenVector< 1 > >" tooltip="DNDS::ArrayPair< DNDS::ArrayEigenVector< 1 > >"]
    "10" [label="DNDS::ArrayPair< DNDS::ArrayEigenVector< NonUniformSize > >" tooltip="DNDS::ArrayPair< DNDS::ArrayEigenVector< NonUniformSize > >"]
    "8" [label="DNDS::ArrayPair< DNDS::ParArray< RecAtr, 1 > >" tooltip="DNDS::ArrayPair< DNDS::ParArray< RecAtr, 1 > >"]
    "17" [label="DNDS::ArrayTransformerType< DNDS::ArrayEigenMatrix< 3, 3 > >" tooltip="DNDS::ArrayTransformerType< DNDS::ArrayEigenMatrix< 3, 3 > >"]
    "13" [label="DNDS::ArrayTransformerType< DNDS::ArrayEigenUniMatrixBatch< 3, 1 > >" tooltip="DNDS::ArrayTransformerType< DNDS::ArrayEigenUniMatrixBatch< 3, 1 > >"]
    "7" [label="DNDS::ArrayTransformerType< DNDS::ArrayEigenVector< 1 > >" tooltip="DNDS::ArrayTransformerType< DNDS::ArrayEigenVector< 1 > >"]
    "11" [label="DNDS::ArrayTransformerType< DNDS::ArrayEigenVector< NonUniformSize > >" tooltip="DNDS::ArrayTransformerType< DNDS::ArrayEigenVector< NonUniformSize > >"]
    "9" [label="DNDS::ArrayTransformerType< DNDS::ParArray< RecAtr, 1 > >" tooltip="DNDS::ArrayTransformerType< DNDS::ParArray< RecAtr, 1 > >"]
    "15" [label="DNDS::ArrayTransformerType< ParArray< real, 1 > >" tooltip="DNDS::ArrayTransformerType< ParArray< real, 1 > >"]
    "14" [label="DNDS::ArrayPair< TArray >" tooltip="DNDS::ArrayPair< TArray >"]
    "2" [label="DNDS::CFV::FiniteVolume" tooltip="DNDS::CFV::FiniteVolume"]
    "5" [label="DNDS::CFV::FiniteVolumeSettings" tooltip="DNDS::CFV::FiniteVolumeSettings"]
    "1" [label="DNDS::CFV::VariationalReconstruction< dim >" tooltip="DNDS::CFV::VariationalReconstruction< dim >" fillcolor="#BFBFBF"]
    "18" [label="DNDS::Geom::Base::CFVPeriodicity" tooltip="DNDS::Geom::Base::CFVPeriodicity"]
    "19" [label="DNDS::Geom::Periodicity" tooltip="DNDS::Geom::Periodicity"]
    "21" [label="DNDS::Geom::tGPointPortable" tooltip="DNDS::Geom::tGPointPortable"]
    "24" [label="DNDS::Geom::tPointPortable" tooltip="DNDS::Geom::tPointPortable"]
    "4" [label="DNDS::MPIInfo" tooltip="DNDS::MPIInfo"]
    "3" [label="DNDS::DeviceTransferable< FiniteVolume >" tooltip="DNDS::DeviceTransferable< FiniteVolume >"]
    "20" [label="std::array< DNDS::Geom::tGPointPortable, 4 >" tooltip="std::array< DNDS::Geom::tGPointPortable, 4 >"]
    "23" [label="std::array< DNDS::Geom::tPointPortable, 4 >" tooltip="std::array< DNDS::Geom::tPointPortable, 4 >"]
    "25" [label="std::array< real, 3 >" tooltip="std::array< real, 3 >"]
    "22" [label="std::array< real, 9 >" tooltip="std::array< real, 9 >"]
    "26" [label="std::set< index >" tooltip="std::set< index >"]
    "16" -> "17" [dir=forward tooltip="usage"]
    "12" -> "13" [dir=forward tooltip="usage"]
    "6" -> "7" [dir=forward tooltip="usage"]
    "10" -> "11" [dir=forward tooltip="usage"]
    "8" -> "9" [dir=forward tooltip="usage"]
    "14" -> "15" [dir=forward tooltip="usage"]
    "2" -> "3" [dir=forward tooltip="public-inheritance"]
    "2" -> "4" [dir=forward tooltip="usage"]
    "2" -> "5" [dir=forward tooltip="usage"]
    "2" -> "6" [dir=forward tooltip="usage"]
    "2" -> "8" [dir=forward tooltip="usage"]
    "2" -> "10" [dir=forward tooltip="usage"]
    "2" -> "12" [dir=forward tooltip="usage"]
    "2" -> "14" [dir=forward tooltip="usage"]
    "2" -> "16" [dir=forward tooltip="usage"]
    "2" -> "18" [dir=forward tooltip="usage"]
    "2" -> "26" [dir=forward tooltip="usage"]
    "1" -> "2" [dir=forward tooltip="public-inheritance"]
    "18" -> "19" [dir=forward tooltip="public-inheritance"]
    "19" -> "20" [dir=forward tooltip="usage"]
    "19" -> "23" [dir=forward tooltip="usage"]
    "21" -> "22" [dir=forward tooltip="usage"]
    "24" -> "25" [dir=forward tooltip="usage"]
    "20" -> "21" [dir=forward tooltip="usage"]
    "23" -> "24" [dir=forward tooltip="usage"]
}

The VR class that provides any information needed in high-order CFV.

VR holds a primitive mesh and any needed derived information about geometry and general reconstruction coefficients

The (differentiable) input of VR is merely: geometry (mesh), and the weight inputs if using distributed weight; output are all values derived using construct* method

Public Functions

inline void parseSettings(VRSettings::json &j)

Backward-compatible overload accepting JSON (used by Python bindings).

template<typename ...Ts>
inline void ApplyPeriodicTransform(int if2c, Geom::t_index faceID, Ts&... data) const

Applies the appropriate periodic transformation to one or more data matrices.

Determines from if2c and faceID whether the current cell is the donor or main side of the periodic pair, and calls FTransPeriodic or FTransPeriodicBack accordingly. No-op when the mesh is not periodic.

Parameters:
  • if2c – Face-to-cell local index (0 = back, 1 = front)

  • faceID – Boundary zone ID of the face

  • data – One or more Eigen matrix references to transform

template<class TOut>
inline void FDiffBaseValue(TOut &DiBj, const Geom::tPoint &pPhy, index iCell, index iFace, int iG, int flag = 0)

flag = 0 means use moment data, or else use no moment (as 0) pPhy must be relative to cell if iFace < 0, means anywhere if iFace > 0, iG == -1, means center; iG < -1, then anywhere

template<class TList>
inline MatrixXR GetIntPointDiffBaseValue(index iCell, index iFace, rowsize if2c, int iG, TList &&diffList = EigenAll, uint8_t maxDiff = UINT8_MAX)

if if2c < 0, then calculated, if maxDiff == 255, then seen as all diffs if iFace < 0, then seen as cell int points; if iG < 1, then seen as center

Todo:

: divide GetIntPointDiffBaseValue into different calls

Todo:

: //TODO add support for rotational periodic boundary!

Warning

maxDiff is max(diffList) + 1 not len(difflist)

inline auto GetMatrixSecondary(index iCell, index iFace, index if2c)

if if2c < 0, then calculated with iCell hen seen as all diffs

template<int nVarsFixed = 5>
void DoReconstruction2ndGrad(tUGrad<nVarsFixed, dim> &uRec, tUDof<nVarsFixed> &u, const TFBoundary<nVarsFixed> &FBoundary, int method)

fallback reconstruction method, explicit 2nd order FV reconstruction

Parameters:
  • uRec – output, reconstructed gradients

  • u – input, mean values

  • FBoundary – see TFBoundary

  • method – currently 1==2nd-order-GaussGreen

template<int nVarsFixed = 5>
void DoReconstruction2nd(tURec<nVarsFixed> &uRec, tUDof<nVarsFixed> &u, const TFBoundary<nVarsFixed> &FBoundary, int method, const std::vector<int> &mask)

fallback reconstruction method, explicit 2nd order FV reconstruction

Parameters:
  • uRec – output, reconstructed gradients

  • u – input, mean values

  • FBoundary – see TFBoundary

  • method – currently 1==2nd-order-GaussGreen

template<int nVarsFixed = 5>
void DoReconstructionIter(tURec<nVarsFixed> &uRec, tURec<nVarsFixed> &uRecNew, tUDof<nVarsFixed> &u, const TFBoundary<nVarsFixed> &FBoundary, bool putIntoNew = false, bool recordInc = false, bool uRecIsZero = false)

iterative variational reconstruction method

DoReconstructionIter updates uRec (could locate into uRecNew) $$ ur_i = A_{i}^{-1}B_{ij}ur_j + A_{i}^{-1}b_{i} $$ which is a fixed-point solver (using SOR/Jacobi sweeping)

if recordInc, value in the output array is actually defined as : $$ -(A_{i}^{-1}B_{ij}ur_j + A_{i}^{-1}b_{i}) + ur_i $$ which is the RHS of Block-Jacobi preconditioned system

Warning

mind that uRec could be overwritten

Parameters:
  • uRec – input/ouput, reconstructed coefficients, is output if putIntoNew==false, otherwise is used as medium value

  • uRecNew – input/ouput, reconstructed coefficients, is output if putIntoNew==true, otherwise is used as medium value

  • u – input, mean values

  • FBoundary – see TFBoundary

  • putIntoNew – put valid output into uRecNew or not

  • recordInc – if true, uRecNew holds the incremental value to this iteration

template<int nVarsFixed = 5>
void DoReconstructionIterDiff(tURec<nVarsFixed> &uRec, tURec<nVarsFixed> &uRecDiff, tURec<nVarsFixed> &uRecNew, tUDof<nVarsFixed> &u, const TFBoundaryDiff<nVarsFixed> &FBoundaryDiff)

puts into uRecNew with Mat * uRecDiff; uses the Block-jacobi preconditioned reconstruction system as Mat.

the matrix operation can be viewed as $$ vr_i = ur_j - A_{i}^{-1}B_{ij}ur_j $$ which is the Jacobian of the DoReconstructionIter’s RHS output. Note that we account for nonlinear effects from BC

Parameters:
  • uRec – input, the base value

  • uRecDiff – input, the diff value, $x$ in $y=Jx$, or the disturbance value

  • uRecNew – output, the result, $y$, or the response value

  • u – input, mean value

  • FBoundary – Vec F(const Vec &uL, const Vec& dUL, const tPoint &unitNorm, const tPoint &p, t_index faceID), with Vec == Eigen::Vector<real, nVarsFixed> uRecDiff should be untouched

template<int nVarsFixed = 5>
void DoReconstructionIterSOR(tURec<nVarsFixed> &uRec, tURec<nVarsFixed> &uRecInc, tURec<nVarsFixed> &uRecNew, tUDof<nVarsFixed> &u, const TFBoundaryDiff<nVarsFixed> &FBoundaryDiff, bool reverse = false)

do a SOR iteration from uRecNew, with uRecInc as the RHSterm of Block-Jacobi preconditioned system

//TODO

template<int nVarsFixed>
void DoLimiterWBAP_C(tUDof<nVarsFixed> &u, tURec<nVarsFixed> &uRec, tURec<nVarsFixed> &uRecNew, tURec<nVarsFixed> &uRecBuf, tScalarPair &si, bool ifAll, const tFMEig<nVarsFixed> &FM, const tFMEig<nVarsFixed> &FMI, bool putIntoNew = false)

FM(uLeft,uRight,norm) gives vsize * vsize mat of Left Eigen Vectors.

template<int nVarsFixed>
void DoLimiterWBAP_3(tUDof<nVarsFixed> &u, tURec<nVarsFixed> &uRec, tURec<nVarsFixed> &uRecNew, tURec<nVarsFixed> &uRecBuf, tScalarPair &si, bool ifAll, const tFMEig<nVarsFixed> &FM, const tFMEig<nVarsFixed> &FMI, bool putIntoNew = false)

FM(uLeft,uRight,norm) gives vsize * vsize mat of Left Eigen Vectors.