DNDSR 0.1.0.dev1+gcd065ad
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
DNDS::CFV::VariationalReconstruction< dim > Class Template Reference

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

#include <VariationalReconstruction.hpp>

Inheritance diagram for DNDS::CFV::VariationalReconstruction< dim >:
[legend]
Collaboration diagram for DNDS::CFV::VariationalReconstruction< dim >:
[legend]

Public Types

using t_base = FiniteVolume
 
using TFTrans = std::function< void(Eigen::Ref< MatrixXR >, Geom::t_index)>
 
using tFGetBoundaryWeight = std::function< real(Geom::t_index, int)>
 
template<int nVarsFixed = 5>
using TFBoundary = std::function< Eigen::Vector< real, nVarsFixed >(const Eigen::Vector< real, nVarsFixed > &, const Eigen::Vector< real, nVarsFixed > &, index, index, int, const Geom::tPoint &, const Geom::tPoint &, Geom::t_index fType)>
 
template<int nVarsFixed = 5>
using TFBoundaryDiff = std::function< Eigen::Vector< real, nVarsFixed >(const Eigen::Vector< real, nVarsFixed > &, const Eigen::Vector< real, nVarsFixed > &, const Eigen::Vector< real, nVarsFixed > &, index, index, int, const Geom::tPoint &, const Geom::tPoint &, Geom::t_index fType)>
 
template<int nVarsFixed>
using TFPost = std::function< void(Eigen::Matrix< real, 1, nVarsFixed > &)>
 
template<int nVarsFixed>
using tLimitBatch = Eigen::Matrix< real, -1, nVarsFixed, 0, CFV::VariationalReconstruction< dim >::maxRecDOFBatch >
 
template<int nVarsFixed>
using tFMEig = std::function< tLimitBatch< nVarsFixed >(const Eigen::Vector< real, nVarsFixed > &uL, const Eigen::Vector< real, nVarsFixed > &uR, const Geom::tPoint &uNorm, const Eigen::Ref< tLimitBatch< nVarsFixed > > &data)>
 
- Public Types inherited from DNDS::CFV::FiniteVolume
template<DeviceBackend B>
using t_deviceView = FiniteVolumeDeviceView< B >
 

Public Member Functions

int getDim () const
 
const VRSettingsgetSettings () const
 
void parseSettings (const VRSettings &s)
 
void parseSettings (VRSettings::json &j)
 Backward-compatible overload accepting JSON (used by Python bindings).
 
 VariationalReconstruction (MPIInfo nMpi, ssp< Geom::UnstructuredMesh > nMesh)
 
void SetPeriodicTransformations (const TFTrans &nFTransPeriodic, const TFTrans &nFTransPeriodicBack)
 
template<size_t pDim>
void SetPeriodicTransformations (std::array< int, pDim > Seq123)
 
void SetPeriodicTransformations ()
 
template<typename... Ts>
void ApplyPeriodicTransform (int if2c, Geom::t_index faceID, Ts &...data) const
 Applies the appropriate periodic transformation to one or more data matrices.
 
void ConstructMetrics ()
 
void ConstructBaseAndWeight (const tFGetBoundaryWeight &id2faceDircWeight=[](Geom::t_index id, int iOrder) { return 1.0;})
 
void ConstructRecCoeff ()
 
auto GetCellRecMatAInv (index iCell)
 
auto GetCellRecMatAInvB (index iCell, int ic2f)
 
autoget_matrixAAInvB ()
 
autoget_vectorAInvB ()
 
template<class TOut >
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 >
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
 
auto GetMatrixSecondary (index iCell, index iFace, index if2c)
 if if2c < 0, then calculated with iCell hen seen as all diffs
 
template<class TDiffIDerived , class TDiffJDerived >
auto FFaceFunctional (const Eigen::MatrixBase< TDiffIDerived > &DiffI, const Eigen::MatrixBase< TDiffJDerived > &DiffJ, index iFace, index iCellL, index iCellR)
 
real GetGreenGauss1WeightOnCell (index iCell)
 
real GetCellAR (index iCell)
 
template<int nVarsFixed = 1>
void MatrixAMult (tURec< nVarsFixed > &uRec, tURec< nVarsFixed > &uRec1)
 
template<int nVarsFixed = 1>
auto DownCastURecOrder (int curOrder, index iCell, tURec< nVarsFixed > &uRec, int downCastMethod)
 
template<int nVarsFixed>
void BuildURec (tURec< nVarsFixed > &u, int nVars, bool buildSon=true, bool buildTrans=true)
 
template<int nVarsFixed>
void BuildUGrad (tUGrad< nVarsFixed, dim > &u, int nVars, bool buildSon=true, bool buildTrans=true)
 
void BuildScalar (tScalarPair &u, bool buildSon=true, bool buildTrans=true)
 
template<int nVarsFixed = 1>
void BuildUDofNode (tUDof< nVarsFixed > &u, int nVars, bool buildSon=true, bool buildTrans=true)
 
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
 
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
 
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
 
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.
 
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
 
template<int nVarsFixed, int nVarsSee = 2>
void DoCalculateSmoothIndicator (tScalarPair &si, tURec< nVarsFixed > &uRec, tUDof< nVarsFixed > &u, const std::array< int, nVarsSee > &varsSee)
 
template<int nVarsFixed>
void DoCalculateSmoothIndicatorV1 (tScalarPair &si, tURec< nVarsFixed > &uRec, tUDof< nVarsFixed > &u, const Eigen::Vector< real, nVarsFixed > &varsSee, const TFPost< nVarsFixed > &FPost)
 
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.
 
void WriteSerializeRecMatrix (const Serializer::SerializerBaseSSP &serializerP)
 
- Public Member Functions inherited from DNDS::CFV::FiniteVolume
const FiniteVolumeSettingsgetSettings () const
 
void parseSettings (FiniteVolumeSettings::json &j)
 
int GetAxisSymmetric () const
 
int SetAxisSymmetric (int v)
 
 FiniteVolume (MPIInfo nMpi, ssp< Geom::UnstructuredMesh > nMesh)
 
int getDim () const
 
template<class TArrayPair , class... TOthers>
void MakePairDefaultOnCell (TArrayPair &aPair, const std::string &name, TOthers... others)
 make pair with default MPI type, match cell layout
 
template<class TArrayPair , class... TOthers>
void MakePairDefaultOnFace (TArrayPair &aPair, const std::string &name, TOthers... others)
 make pair with default MPI type, match face layout
 
void SetCellAtrBasic ()
 
void ConstructCellVolume ()
 
void ConstructCellBary ()
 
void ConstructCellCent ()
 
void ConstructCellIntJacobiDet ()
 
void ConstructCellIntPPhysics ()
 
void ConstructCellAlignedHBox ()
 
void ConstructCellMajorHBoxCoordInertia ()
 
void SetFaceAtrBasic ()
 
void ConstructFaceArea ()
 
void ConstructFaceCent ()
 
void ConstructFaceIntJacobiDet ()
 
void ConstructFaceIntPPhysics ()
 
void ConstructFaceUnitNorm ()
 
void ConstructFaceMeanNorm ()
 
void ConstructCellSmoothScale ()
 
template<int nVarsFixed = 1>
void BuildUDof (tUDof< nVarsFixed > &u, int nVars, bool buildSon=true, bool buildTrans=true, Geom::MeshLoc varloc=Geom::MeshLoc::Cell)
 
template<int nVarsFixed, int dim>
void BuildUGradD (tUGrad< nVarsFixed, dim > &u, int nVars, bool buildSon=true, bool buildTrans=true, Geom::MeshLoc varloc=Geom::MeshLoc::Cell)
 
RecAtrGetCellAtr (index iCell)
 
int GetCellOrder (index iCell)
 
RecAtrGetFaceAtr (index iFace)
 
int maxNDOF ()
 
real GetCellVol (index iCell) const
 
real GetFaceArea (index iFace) const
 
real GetCellSmoothScaleRatio (index iCell) const
 
real GetGlobalVol () const
 
real GetCellJacobiDet (index iCell, rowsize iG)
 
real GetFaceJacobiDet (index iFace, rowsize iG)
 
Geom::Elem::Quadrature GetFaceQuad (index iFace) const
 
Geom::Elem::Quadrature GetFaceQuadO1 (index iFace) const
 
real GetFaceParamArea (index iFace) const
 
Geom::Elem::Quadrature GetCellQuad (index iCell) const
 
Geom::Elem::Quadrature GetCellQuadO1 (index iCell) const
 
real GetCellParamVol (index iCell) const
 
Geom::tPoint GetCellBary (index iCell)
 
bool CellIsFaceBack (index iCell, index iFace) const
 
index CellFaceOther (index iCell, index iFace) const
 
Geom::tPoint GetFaceNorm (index iFace, int iG) const
 
Geom::tPoint GetFaceNormFromCell (index iFace, index iCell, rowsize if2c, int iG)
 
Geom::tPoint GetFaceQuadraturePPhys (index iFace, int iG)
 
Geom::tPoint GetFaceQuadraturePPhysFromCell (index iFace, index iCell, rowsize if2c, int iG)
 
Geom::tPoint GetFacePointFromCell (index iFace, index iCell, rowsize if2c, const Geom::tPoint &pnt)
 
Geom::tPoint GetOtherCellBaryFromCell (index iCell, index iCellOther, index iFace)
 
Geom::tPoint GetOtherCellPointFromCell (index iCell, index iCellOther, index iFace, const Geom::tPoint &pnt)
 
Geom::tGPoint GetOtherCellInertiaFromCell (index iCell, index iCellOther, index iFace)
 
Geom::tPoint GetCellQuadraturePPhys (index iCell, int iG)
 
real GetCellMaxLenScale (index iCell)
 
real GetCellNodeMinLenScale (index iCell)
 
index getArrayBytes ()
 
template<DeviceBackend B>
t_deviceView< B > deviceView ()
 
- Public Member Functions inherited from DNDS::DeviceTransferable< FiniteVolume >
void to_device (DeviceBackend B)
 Mirror every registered ArrayPair to the given device backend.
 
void to_host ()
 Pull every registered pair back to host memory.
 
DeviceBackend device ()
 Consistent device backend across all registered pairs.
 
index getDeviceArrayBytes ()
 Total footprint of every registered father+son pair in bytes.
 

Static Public Attributes

static const int maxRecDOFBatch = dim == 2 ? 4 : 10
 
static const int maxRecDOF = dim == 2 ? 9 : 19
 
static const int maxNDiff = dim == 2 ? 10 : 20
 
static const int maxNeighbour = 7
 

Additional Inherited Members

- Public Attributes inherited from DNDS::CFV::FiniteVolume
MPI_int mRank {0}
 
MPIInfo mpi
 
ssp< Geom::UnstructuredMeshmesh
 
- Protected Member Functions inherited from DNDS::CFV::FiniteVolume
auto device_array_list ()
 constructed using ConstructMetrics()
 
template<typename F >
void for_each_device_member (F &&f)
 
- Protected Attributes inherited from DNDS::CFV::FiniteVolume
FiniteVolumeSettings settings
 
real sumVolume {0}
 
real minVolume {veryLargeReal}
 
real maxVolume {0}
 
real volGlobal {0}
 
tScalarPair volumeLocal
 constructed using ConstructMetrics()
 
tScalarPair faceArea
 constructed using ConstructMetrics()
 
tRecAtrPair cellAtr
 constructed using ConstructMetrics()
 
tRecAtrPair faceAtr
 constructed using ConstructMetrics()
 
tCoeffPair cellIntJacobiDet
 constructed using ConstructMetrics()
 
tCoeffPair faceIntJacobiDet
 constructed using ConstructMetrics()
 
t3VecsPair faceUnitNorm
 constructed using ConstructMetrics()
 
t3VecPair faceMeanNorm
 constructed using ConstructMetrics()
 
t3VecPair cellBary
 constructed using ConstructMetrics()
 
t3VecPair faceCent
 constructed using ConstructMetrics()
 
t3VecPair cellCent
 constructed using ConstructMetrics()
 
t3VecsPair cellIntPPhysics
 constructed using ConstructMetrics()
 
t3VecsPair faceIntPPhysics
 constructed using ConstructMetrics()
 
t3VecPair cellAlignedHBox
 constructed using ConstructMetrics()
 
t3VecPair cellMajorHBox
 constructed using ConstructMetrics()
 
t3MatPair cellMajorCoord
 constructed using ConstructMetrics()
 
t3MatPair cellInertia
 constructed using ConstructMetrics()
 
tScalarPair cellSmoothScale
 constructed using ConstructMetrics()
 
Geom::Base::CFVPeriodicity periodicity
 
int axisSymmetric = 0
 
std::set< indexaxisFaces
 

Detailed Description

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

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

Definition at line 109 of file VariationalReconstruction.hpp.

Member Typedef Documentation

◆ t_base

template<int dim = 2>
using DNDS::CFV::VariationalReconstruction< dim >::t_base = FiniteVolume

Definition at line 113 of file VariationalReconstruction.hpp.

◆ TFBoundary

template<int dim = 2>
template<int nVarsFixed = 5>
using DNDS::CFV::VariationalReconstruction< dim >::TFBoundary = std::function<Eigen::Vector<real, nVarsFixed>( const Eigen::Vector<real, nVarsFixed> &, const Eigen::Vector<real, nVarsFixed> &, index, index, int, const Geom::tPoint &, const Geom::tPoint &, Geom::t_index fType )>

Definition at line 894 of file VariationalReconstruction.hpp.

◆ TFBoundaryDiff

template<int dim = 2>
template<int nVarsFixed = 5>
using DNDS::CFV::VariationalReconstruction< dim >::TFBoundaryDiff = std::function<Eigen::Vector<real, nVarsFixed>( const Eigen::Vector<real, nVarsFixed> &, const Eigen::Vector<real, nVarsFixed> &, const Eigen::Vector<real, nVarsFixed> &, index, index, int, const Geom::tPoint &, const Geom::tPoint &, Geom::t_index fType )>

Definition at line 904 of file VariationalReconstruction.hpp.

◆ tFGetBoundaryWeight

template<int dim = 2>
using DNDS::CFV::VariationalReconstruction< dim >::tFGetBoundaryWeight = std::function<real(Geom::t_index, int)>

Definition at line 284 of file VariationalReconstruction.hpp.

◆ tFMEig

template<int dim = 2>
template<int nVarsFixed>
using DNDS::CFV::VariationalReconstruction< dim >::tFMEig = std::function<tLimitBatch<nVarsFixed>( const Eigen::Vector<real, nVarsFixed> &uL, const Eigen::Vector<real, nVarsFixed> &uR, const Geom::tPoint &uNorm, const Eigen::Ref<tLimitBatch<nVarsFixed> > &data)>

Definition at line 1060 of file VariationalReconstruction.hpp.

◆ TFPost

template<int dim = 2>
template<int nVarsFixed>
using DNDS::CFV::VariationalReconstruction< dim >::TFPost = std::function<void(Eigen::Matrix<real, 1, nVarsFixed> &)>

Definition at line 1043 of file VariationalReconstruction.hpp.

◆ TFTrans

template<int dim = 2>
using DNDS::CFV::VariationalReconstruction< dim >::TFTrans = std::function<void(Eigen::Ref<MatrixXR>, Geom::t_index)>

Definition at line 143 of file VariationalReconstruction.hpp.

◆ tLimitBatch

template<int dim = 2>
template<int nVarsFixed>
using DNDS::CFV::VariationalReconstruction< dim >::tLimitBatch = Eigen::Matrix<real, -1, nVarsFixed, 0, CFV::VariationalReconstruction<dim>::maxRecDOFBatch>

Definition at line 1057 of file VariationalReconstruction.hpp.

Constructor & Destructor Documentation

◆ VariationalReconstruction()

template<int dim = 2>
DNDS::CFV::VariationalReconstruction< dim >::VariationalReconstruction ( MPIInfo  nMpi,
ssp< Geom::UnstructuredMesh nMesh 
)
inline

Definition at line 222 of file VariationalReconstruction.hpp.

Member Function Documentation

◆ ApplyPeriodicTransform()

template<int dim = 2>
template<typename... Ts>
void DNDS::CFV::VariationalReconstruction< dim >::ApplyPeriodicTransform ( int  if2c,
Geom::t_index  faceID,
Ts &...  data 
) const
inline

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
if2cFace-to-cell local index (0 = back, 1 = front)
faceIDBoundary zone ID of the face
dataOne or more Eigen matrix references to transform

Definition at line 269 of file VariationalReconstruction.hpp.

Here is the call graph for this function:

◆ BuildScalar()

template<int dim = 2>
void DNDS::CFV::VariationalReconstruction< dim >::BuildScalar ( tScalarPair u,
bool  buildSon = true,
bool  buildTrans = true 
)
inline

Definition at line 852 of file VariationalReconstruction.hpp.

Here is the call graph for this function:

◆ BuildUDofNode()

template<int dim = 2>
template<int nVarsFixed = 1>
void DNDS::CFV::VariationalReconstruction< dim >::BuildUDofNode ( tUDof< nVarsFixed > &  u,
int  nVars,
bool  buildSon = true,
bool  buildTrans = true 
)
inline

Definition at line 873 of file VariationalReconstruction.hpp.

◆ BuildUGrad()

template<int dim = 2>
template<int nVarsFixed>
void DNDS::CFV::VariationalReconstruction< dim >::BuildUGrad ( tUGrad< nVarsFixed, dim > &  u,
int  nVars,
bool  buildSon = true,
bool  buildTrans = true 
)
inline

Definition at line 831 of file VariationalReconstruction.hpp.

◆ BuildURec()

template<int dim = 2>
template<int nVarsFixed>
void DNDS::CFV::VariationalReconstruction< dim >::BuildURec ( tURec< nVarsFixed > &  u,
int  nVars,
bool  buildSon = true,
bool  buildTrans = true 
)
inline

Definition at line 808 of file VariationalReconstruction.hpp.

Here is the call graph for this function:

◆ ConstructBaseAndWeight()

template<int dim>
template void DNDS::CFV::VariationalReconstruction< dim >::ConstructBaseAndWeight ( const tFGetBoundaryWeight id2faceDircWeight = [](Geom::t_index id, int iOrder) { return 1.0;})

note: if use cell-2-cell distance, note periodic, use wrapped call here

Definition at line 55 of file VariationalReconstruction.cpp.

Here is the call graph for this function:

◆ ConstructMetrics()

template<int dim>
template void DNDS::CFV::VariationalReconstruction< dim >::ConstructMetrics ( )

Definition at line 9 of file VariationalReconstruction.cpp.

Here is the call graph for this function:

◆ ConstructRecCoeff()

template<int dim>
template void DNDS::CFV::VariationalReconstruction< dim >::ConstructRecCoeff ( )

Definition at line 509 of file VariationalReconstruction.cpp.

Here is the call graph for this function:

◆ DoCalculateSmoothIndicator()

template<int dim>
template<int nVarsFixed, int nVarsSee>
void DNDS::CFV::VariationalReconstruction< dim >::DoCalculateSmoothIndicator ( tScalarPair si,
tURec< nVarsFixed > &  uRec,
tUDof< nVarsFixed > &  u,
const std::array< int, nVarsSee > &  varsSee 
)

Definition at line 42 of file VariationalReconstruction_LimiterProcedure.hxx.

Here is the call graph for this function:

◆ DoCalculateSmoothIndicatorV1()

template<int dim>
template<int nVarsFixed>
void DNDS::CFV::VariationalReconstruction< dim >::DoCalculateSmoothIndicatorV1 ( tScalarPair si,
tURec< nVarsFixed > &  uRec,
tUDof< nVarsFixed > &  u,
const Eigen::Vector< real, nVarsFixed > &  varsSee,
const TFPost< nVarsFixed > &  FPost 
)

Definition at line 112 of file VariationalReconstruction_LimiterProcedure.hxx.

Here is the call graph for this function:

◆ DoLimiterWBAP_3()

template<int dim>
template<int nVarsFixed>
void DNDS::CFV::VariationalReconstruction< dim >::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.

Definition at line 338 of file VariationalReconstruction_LimiterProcedure.hxx.

Here is the call graph for this function:

◆ DoLimiterWBAP_C()

template<int dim>
template<int nVarsFixed>
void DNDS::CFV::VariationalReconstruction< dim >::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.

no lim need to copy !!!!

Definition at line 190 of file VariationalReconstruction_LimiterProcedure.hxx.

Here is the call graph for this function:

◆ DoReconstruction2nd()

template<int dim = 2>
template<int nVarsFixed = 5>
void DNDS::CFV::VariationalReconstruction< dim >::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
uRecoutput, reconstructed gradients
uinput, mean values
FBoundarysee TFBoundary
methodcurrently 1==2nd-order-GaussGreen

◆ DoReconstruction2ndGrad()

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

fallback reconstruction method, explicit 2nd order FV reconstruction

Parameters
uRecoutput, reconstructed gradients
uinput, mean values
FBoundarysee TFBoundary
methodcurrently 1==2nd-order-GaussGreen

◆ DoReconstructionIter()

template<int dim>
template<int nVarsFixed>
void DNDS::CFV::VariationalReconstruction< dim >::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

Parameters
uRecinput/ouput, reconstructed coefficients, is output if putIntoNew==false, otherwise is used as medium value
uRecNewinput/ouput, reconstructed coefficients, is output if putIntoNew==true, otherwise is used as medium value
uinput, mean values
FBoundarysee TFBoundary
putIntoNewput valid output into uRecNew or not
recordIncif true, uRecNew holds the incremental value to this iteration
Warning
mind that uRec could be overwritten

Definition at line 487 of file VariationalReconstruction_Reconstruction.hxx.

◆ DoReconstructionIterDiff()

template<int dim>
template<int nVarsFixed>
void DNDS::CFV::VariationalReconstruction< dim >::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
uRecinput, the base value
uRecDiffinput, the diff value, $x$ in $y=Jx$, or the disturbance value
uRecNewoutput, the result, $y$, or the response value
uinput, mean value
FBoundaryVec 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

Definition at line 619 of file VariationalReconstruction_Reconstruction.hxx.

◆ DoReconstructionIterSOR()

template<int dim>
template<int nVarsFixed>
void DNDS::CFV::VariationalReconstruction< dim >::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

Definition at line 674 of file VariationalReconstruction_Reconstruction.hxx.

◆ DownCastURecOrder()

template<int dim = 2>
template<int nVarsFixed = 1>
auto DNDS::CFV::VariationalReconstruction< dim >::DownCastURecOrder ( int  curOrder,
index  iCell,
tURec< nVarsFixed > &  uRec,
int  downCastMethod 
)
inline

Definition at line 726 of file VariationalReconstruction.hpp.

◆ FDiffBaseValue()

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

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

Definition at line 318 of file VariationalReconstruction.hpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ FFaceFunctional()

template<int dim = 2>
auto DNDS::CFV::VariationalReconstruction< dim >::FFaceFunctional ( const Eigen::MatrixBase< TDiffIDerived > &  DiffI,
const Eigen::MatrixBase< TDiffJDerived > &  DiffJ,
index  iFace,
index  iCellL,
index  iCellR 
)
inline

Definition at line 460 of file VariationalReconstruction.hpp.

Here is the call graph for this function:

◆ get_matrixAAInvB()

template<int dim = 2>
auto & DNDS::CFV::VariationalReconstruction< dim >::get_matrixAAInvB ( )
inline

Definition at line 301 of file VariationalReconstruction.hpp.

◆ get_vectorAInvB()

template<int dim = 2>
auto & DNDS::CFV::VariationalReconstruction< dim >::get_vectorAInvB ( )
inline

Definition at line 306 of file VariationalReconstruction.hpp.

◆ GetCellAR()

template<int dim = 2>
real DNDS::CFV::VariationalReconstruction< dim >::GetCellAR ( index  iCell)
inline

Definition at line 710 of file VariationalReconstruction.hpp.

Here is the caller graph for this function:

◆ GetCellRecMatAInv()

template<int dim = 2>
auto DNDS::CFV::VariationalReconstruction< dim >::GetCellRecMatAInv ( index  iCell)
inline

Definition at line 291 of file VariationalReconstruction.hpp.

◆ GetCellRecMatAInvB()

template<int dim = 2>
auto DNDS::CFV::VariationalReconstruction< dim >::GetCellRecMatAInvB ( index  iCell,
int  ic2f 
)
inline

Definition at line 296 of file VariationalReconstruction.hpp.

◆ getDim()

template<int dim = 2>
int DNDS::CFV::VariationalReconstruction< dim >::getDim ( ) const
inline

Definition at line 112 of file VariationalReconstruction.hpp.

◆ GetGreenGauss1WeightOnCell()

template<int dim = 2>
real DNDS::CFV::VariationalReconstruction< dim >::GetGreenGauss1WeightOnCell ( index  iCell)
inline

Definition at line 695 of file VariationalReconstruction.hpp.

Here is the call graph for this function:

◆ GetIntPointDiffBaseValue()

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

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
Warning
maxDiff is max(diffList) + 1 not len(difflist)
Todo:
: //TODO add support for rotational periodic boundary!
Todo:
//!!!!TODO: take care of periodic case

Definition at line 380 of file VariationalReconstruction.hpp.

Here is the call graph for this function:

◆ GetMatrixSecondary()

template<int dim = 2>
auto DNDS::CFV::VariationalReconstruction< dim >::GetMatrixSecondary ( index  iCell,
index  iFace,
index  if2c 
)
inline

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

Definition at line 448 of file VariationalReconstruction.hpp.

Here is the call graph for this function:

◆ getSettings()

template<int dim = 2>
const VRSettings & DNDS::CFV::VariationalReconstruction< dim >::getSettings ( ) const
inline

Definition at line 119 of file VariationalReconstruction.hpp.

◆ MatrixAMult()

template<int dim = 2>
template<int nVarsFixed = 1>
void DNDS::CFV::VariationalReconstruction< dim >::MatrixAMult ( tURec< nVarsFixed > &  uRec,
tURec< nVarsFixed > &  uRec1 
)
inline

Definition at line 718 of file VariationalReconstruction.hpp.

◆ parseSettings() [1/2]

template<int dim = 2>
void DNDS::CFV::VariationalReconstruction< dim >::parseSettings ( const VRSettings s)
inline

Definition at line 124 of file VariationalReconstruction.hpp.

Here is the caller graph for this function:

◆ parseSettings() [2/2]

template<int dim = 2>
void DNDS::CFV::VariationalReconstruction< dim >::parseSettings ( VRSettings::json j)
inline

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

Definition at line 135 of file VariationalReconstruction.hpp.

Here is the call graph for this function:

◆ SetPeriodicTransformations() [1/3]

template<int dim = 2>
void DNDS::CFV::VariationalReconstruction< dim >::SetPeriodicTransformations ( )
inline

Definition at line 250 of file VariationalReconstruction.hpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetPeriodicTransformations() [2/3]

template<int dim = 2>
void DNDS::CFV::VariationalReconstruction< dim >::SetPeriodicTransformations ( const TFTrans nFTransPeriodic,
const TFTrans nFTransPeriodicBack 
)
inline

Definition at line 228 of file VariationalReconstruction.hpp.

◆ SetPeriodicTransformations() [3/3]

template<int dim = 2>
template<size_t pDim>
void DNDS::CFV::VariationalReconstruction< dim >::SetPeriodicTransformations ( std::array< int, pDim Seq123)
inline

Definition at line 237 of file VariationalReconstruction.hpp.

Here is the call graph for this function:

◆ WriteSerializeRecMatrix()

template<int dim = 2>
void DNDS::CFV::VariationalReconstruction< dim >::WriteSerializeRecMatrix ( const Serializer::SerializerBaseSSP serializerP)
inline

Definition at line 1096 of file VariationalReconstruction.hpp.

Member Data Documentation

◆ maxNDiff

template<int dim = 2>
const int DNDS::CFV::VariationalReconstruction< dim >::maxNDiff = dim == 2 ? 10 : 20
static

Definition at line 1053 of file VariationalReconstruction.hpp.

◆ maxNeighbour

template<int dim = 2>
const int DNDS::CFV::VariationalReconstruction< dim >::maxNeighbour = 7
static

Definition at line 1054 of file VariationalReconstruction.hpp.

◆ maxRecDOF

template<int dim = 2>
const int DNDS::CFV::VariationalReconstruction< dim >::maxRecDOF = dim == 2 ? 9 : 19
static

Definition at line 1052 of file VariationalReconstruction.hpp.

◆ maxRecDOFBatch

template<int dim = 2>
const int DNDS::CFV::VariationalReconstruction< dim >::maxRecDOFBatch = dim == 2 ? 4 : 10
static

Definition at line 1051 of file VariationalReconstruction.hpp.


The documentation for this class was generated from the following files: