DNDSR 0.1.0.dev1+gcd065ad
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
VRSettings.hpp
Go to the documentation of this file.
1#pragma once
2
3// #ifndef __DNDS_REALLY_COMPILING__
4// #define __DNDS_REALLY_COMPILING__
5// #define __DNDS_REALLY_COMPILING__HEADER_ON__
6// #endif
7#include "DNDS/Defines.hpp"
8#include "DNDS/MPI.hpp"
9#include "Geom/Quadrature.hpp"
10#include "Geom/Mesh.hpp"
14#include "Geom/BaseFunction.hpp"
15#include "Limiters.hpp"
16#include "DNDS/JsonUtil.hpp"
17#include "DNDS/ConfigParam.hpp"
18#include "DNDS/ConfigEnum.hpp"
20// #ifdef __DNDS_REALLY_COMPILING__HEADER_ON__
21// #undef __DNDS_REALLY_COMPILING__
22// #endif
23
24namespace DNDS::CFV
25{
26 /**
27 * @brief
28 * A means to translate nlohmann json into c++ primitive data types and back;
29 * and stores then during computation.
30 *
31 */
33 {
34 using json = nlohmann::ordered_json;
36
37 int intOrderVR{-1}; /// @brief integration degree for VR matrices, <0 means using intOrder
38 int intOrderVRBC{-2}; /// @brief integration degree for VR matrices on BC faces, -1 means using int Order, < -1 means using intOrderVRValue() (the same as intOrderVR)
39 bool cacheDiffBase = false; /// @brief if cache the base function values on each of the quadrature points
40 uint8_t cacheDiffBaseSize = UINT8_MAX;
41
42 real jacobiRelax = 1.0; /// @brief VR SOR/Jacobi iteration relaxation factor
43 bool SORInstead = true; /// @brief use SOR instead of relaxed Jacobi iteration
44
45 real smoothThreshold = 0.01; /// @brief limiter's smooth indicator threshold
46 real WBAP_nStd = 10; /// @brief n used in WBAP limiters
47 bool normWBAP = false; /// @brief if switch to normWBAP
48 int limiterBiwayAlter = 0; /// @brief 0=wbap-L2-biway, 1=minmod-biway
49 int subs2ndOrder = 0; /// @brief 0: vfv; 1: gauss rule; 2: least square; 11: GGMP;
50 int subs2ndOrderGGScheme = 1; /// @brief 1: gauss rule using distance for interpolation; 0: no interpolation
51 real svdTolerance = 0; /// @brief tolerance used in svd
52
54
56 {
57 bool localOrientation = false;
58 bool anisotropicLengths = false;
59
61
63 {
64 DNDS_FIELD(localOrientation, "Use local orientation for basis");
65 DNDS_FIELD(anisotropicLengths, "Use anisotropic length scales");
66 }
68
70 {
71 enum class ScaleType
72 {
73 UnknownScale = -1,
74 MeanAACBB = 0,
75 BaryDiff = 1,
76 CellMax = 2,
78
80
89
91
92 std::array<real, 5> manualDirWeights{{1, 1, 0.5, 1. / 6, 1. / 24}};
93
101
106
109
121
123
127 int greenGaussSpacial = 0; // 1 for uniform weight
128
130
132 {
133 DNDS_FIELD(scaleType, "Functional scale type");
134 DNDS_FIELD(scaleMultiplier, "Functional scale multiplier",
136 DNDS_FIELD(dirWeightScheme, "Directional weight scheme");
137 DNDS_FIELD(dirWeightCombPowV, "Directional weight combination power");
138 DNDS_FIELD(manualDirWeights, "Manual directional weights vector");
139 DNDS_FIELD(geomWeightScheme, "Geometric weight scheme");
140 DNDS_FIELD(geomWeightPower, "Geometric weight power");
141 DNDS_FIELD(geomWeightPower1, "Geometric weight power 1");
142 DNDS_FIELD(geomWeightPower2, "Geometric weight power 2");
143 DNDS_FIELD(useAnisotropicFunctional, "Use anisotropic functional");
144 DNDS_FIELD(tanWeightScale, "Tangential weight scale");
145 DNDS_FIELD(anisotropicType, "Anisotropic functional type");
146 DNDS_FIELD(inertiaWeightPower, "Inertia weight power");
147 DNDS_FIELD(geomWeightBias, "Geometric weight bias");
148 DNDS_FIELD(greenGauss1Weight, "Green-Gauss type-1 weight");
149 DNDS_FIELD(greenGauss1Bias, "Green-Gauss type-1 bias");
150 DNDS_FIELD(greenGauss1Penalty, "Green-Gauss type-1 penalty");
151 DNDS_FIELD(greenGaussSpacial, "Green-Gauss spatial mode: 0=default, 1=uniform");
152 }
155
156 // VRSettings()
157 // {
158 // }
159
162 {
163 cacheDiffBaseSize = uint8_t(dim + 1);
164 }
165
167 {
168 // Base class fields (FiniteVolumeSettings) — flattened into the same JSON object.
169 // Cast base-class pointer-to-member to derived type for template deduction.
170 config.field(static_cast<int T::*>(&T::maxOrder), "maxOrder",
171 "Polynomial degree of reconstruction",
173 config.field(static_cast<int T::*>(&T::intOrder), "intOrder",
174 "Global integration degree",
176 config.field(static_cast<bool T::*>(&T::ignoreMeshGeometryDeficiency), "ignoreMeshGeometryDeficiency",
177 "Ignore mesh geometry deficiency warnings");
178 config.field(static_cast<int T::*>(&T::nIterCellSmoothScale), "nIterCellSmoothScale",
179 "Cell smooth scale iterations",
181
182 // VRSettings own fields
183 DNDS_FIELD(intOrderVR, "VR integration degree (<0 = use intOrder)");
184 DNDS_FIELD(intOrderVRBC, "VR BC integration degree (-1=intOrder, <-1=intOrderVR)");
185 DNDS_FIELD(cacheDiffBase, "Cache base function values at quadrature points");
186 DNDS_FIELD(cacheDiffBaseSize, "Cached diff base size (dim+1)");
187 DNDS_FIELD(jacobiRelax, "VR SOR/Jacobi relaxation factor",
189 DNDS_FIELD(SORInstead, "Use SOR instead of relaxed Jacobi");
190 DNDS_FIELD(smoothThreshold, "Smooth indicator threshold",
192 DNDS_FIELD(WBAP_nStd, "WBAP limiter n parameter",
194 DNDS_FIELD(normWBAP, "Use normWBAP limiter variant");
195 DNDS_FIELD(limiterBiwayAlter, "Limiter biway alter: 0=wbap-L2, 1=minmod");
196 DNDS_FIELD(subs2ndOrder, "2nd order substitution: 0=vfv, 1=gauss, 2=LS, 11=GGMP");
197 DNDS_FIELD(subs2ndOrderGGScheme, "2nd order GG scheme: 0=no interp, 1=distance interp");
198 DNDS_FIELD(svdTolerance, "SVD tolerance for reconstruction",
200 config.field_section(&T::baseSettings, "baseSettings",
201 "Basis function settings");
202 config.field_section(&T::functionalSettings, "functionalSettings",
203 "Functional/weight settings");
204 DNDS_FIELD(bcWeight, "Boundary condition weight",
206 }
207
208 /// @brief Backward-compatible write (used by Python bindings).
209 void WriteIntoJson(json &jsonSetting) const
210 {
211 to_json(jsonSetting, *this);
212 }
213
214 /// @brief Backward-compatible read (used by Python bindings).
215 void ParseFromJson(const json &jsonSetting)
216 {
217 from_json(jsonSetting, *this);
218 }
219
220 [[nodiscard]] bool intOrderVRIsSame() const { return intOrderVR == intOrder || intOrderVR < 0; }
221 [[nodiscard]] int intOrderVRValue() const { return intOrderVR < 0 ? intOrder : intOrderVR; }
222
223 [[nodiscard]] bool intOrderVRBCIsSame() const
224 {
225 if (intOrderVRBC >= 0)
226 return intOrderVRBC == intOrder;
227 if (intOrderVRBC == -1)
228 return true;
229 return intOrderVRIsSame();
230 }
231 [[nodiscard]] int intOrderVRBCValue() const
232 {
233 if (intOrderVRBC >= 0)
234 return intOrderVRBC;
235 if (intOrderVRBC == -1)
236 return intOrder;
237 return intOrderVRValue();
238 }
239 };
240
255
275}
Eigen-matrix array: each row is an Eigen::Map<Matrix> over contiguous real storage.
Batch of uniform-sized Eigen matrices per row, with variable batch count.
Eigen-vector array: each row is an Eigen::Map over contiguous real storage.
Extended enum-to-JSON serialization macro that also exposes allowed string values for JSON Schema gen...
#define DNDS_DEFINE_ENUM_JSON(EnumType_,...)
Define JSON serialization for an enum AND expose its allowed string values.
pybind11-style configuration registration with macro-based field declaration and namespace-scoped tag...
#define DNDS_FIELD(name_, desc_,...)
Register a field inside a DNDS_DECLARE_CONFIG body.
#define DNDS_DECLARE_CONFIG(Type_)
Open a config section registration body.
Core type aliases, constants, and metaprogramming utilities for the DNDS framework.
#define DNDS_DEVICE_TRIVIAL_COPY_DEFINE(T, T_Self)
Definition Defines.hpp:83
#define DNDS_DEVICE_CALLABLE
Definition Defines.hpp:76
#define DNDS_HOST
Definition Defines.hpp:78
#define DNDS_DEVICE_TRIVIAL_COPY_DEFINE_NO_EMPTY_CTOR(T, T_Self)
Definition Defines.hpp:91
JSON-to-Eigen conversion utilities and nlohmann_json helper macros.
MPI wrappers: MPIInfo, collective operations, type mapping, CommStrategy.
RangeTag range(double min)
Create a minimum-only range constraint.
void from_json(const nlohmann::ordered_json &j, host_device_vector< real > &v)
Definition JsonUtil.hpp:183
double real
Canonical floating-point scalar used throughout DNDSR (double precision).
Definition Defines.hpp:105
void to_json(nlohmann::ordered_json &j, const host_device_vector< real > &v)
Definition JsonUtil.hpp:177
A means to translate nlohmann json into c++ primitive data types and back; and stores then during com...
int intOrder
polynomial degree of reconstruction
DNDS_DEVICE_CALLABLE FiniteVolumeSettings()=default
enum DNDS::CFV::VRSettings::FunctionalSettings::DirWeightScheme dirWeightScheme
enum DNDS::CFV::VRSettings::FunctionalSettings::ScaleType scaleType
enum DNDS::CFV::VRSettings::FunctionalSettings::AnisotropicType anisotropicType
DNDS_DEVICE_CALLABLE FunctionalSettings()=default
enum DNDS::CFV::VRSettings::FunctionalSettings::GeomWeightScheme geomWeightScheme
A means to translate nlohmann json into c++ primitive data types and back; and stores then during com...
real WBAP_nStd
limiter's smooth indicator threshold
DNDS_DECLARE_CONFIG(VRSettings)
int limiterBiwayAlter
if switch to normWBAP
bool SORInstead
VR SOR/Jacobi iteration relaxation factor.
struct DNDS::CFV::VRSettings::BaseSettings baseSettings
struct DNDS::CFV::VRSettings::FunctionalSettings functionalSettings
int intOrderVRBC
integration degree for VR matrices, <0 means using intOrder
bool normWBAP
n used in WBAP limiters
bool intOrderVRBCIsSame() const
real smoothThreshold
use SOR instead of relaxed Jacobi iteration
real svdTolerance
1: gauss rule using distance for interpolation; 0: no interpolation
int intOrderVRValue() const
nlohmann::ordered_json json
bool cacheDiffBase
integration degree for VR matrices on BC faces, -1 means using int Order, < -1 means using intOrderVR...
void WriteIntoJson(json &jsonSetting) const
Backward-compatible write (used by Python bindings).
bool intOrderVRIsSame() const
void ParseFromJson(const json &jsonSetting)
Backward-compatible read (used by Python bindings).
uint8_t cacheDiffBaseSize
if cache the base function values on each of the quadrature points
int intOrderVRBCValue() const
DNDS_HOST VRSettings()=default
real bcWeight
tolerance used in svd
int subs2ndOrderGGScheme
0: vfv; 1: gauss rule; 2: least square; 11: GGMP;
DNDS_HOST VRSettings(int dim)
int subs2ndOrder
0=wbap-L2-biway, 1=minmod-biway