35template<
class Form,
class ExtendedStencil,
class Polynomial>
45 MeshObject_type(
mesh),
50 #ifdef SPHERICAL_GEOMETRY
53 dim_(
mesh.nGeometricD()),
62 <<
" should be between zero and 3"
70template<
class FitDataType,
class ExtendedStencil,
class Polynomial>
83 #ifndef SPHERICAL_GEOMETRY
84 if (
mesh.nGeometricD() <= 2)
86 if (
mesh.geometricD()[0] == -1)
112 kdir -= (idir & kdir)*idir;
114 scalar magk =
mag(kdir);
131template<
class FitDataType,
class ExtendedStencil,
class Polynomial>
135 const List<point>&
C,
143 findFaceDirs(idir, jdir, kdir, facei);
147 wts[0] = centralWeight_;
148 if (linearCorrection_)
150 wts[1] = centralWeight_;
172 d.x() = (
p -
p0)&idir;
173 d.y() = (
p -
p0)&jdir;
174 #ifndef SPHERICAL_GEOMETRY
175 d.z() = (
p -
p0)&kdir;
188 Polynomial::addCoeffs
198 for (label i = 0; i <
B.m(); i++)
205 label stencilSize =
C.size();
206 coeffsi.setSize(stencilSize);
208 bool goodFit =
false;
209 for (
int iIt = 0; iIt < 8 && !goodFit; iIt++)
217 for (label i=0; i<stencilSize; i++)
219 coeffsi[i] = wts[0]*wts[i]*invB(0, i);
220 if (
mag(coeffsi[i]) > maxCoeff)
222 maxCoeff =
mag(coeffsi[i]);
227 if (linearCorrection_)
230 (
mag(coeffsi[0] - wLin) < linearLimitFactor_*wLin)
231 && (
mag(coeffsi[1] - (1 - wLin)) < linearLimitFactor_*(1 - wLin))
238 (
mag(coeffsi[0] - 1.0) < linearLimitFactor_*1.0)
267 if (linearCorrection_)
272 for (label j = 0; j <
B.n(); j++)
278 for (label i = 0; i <
B.m(); i++)
288 if (linearCorrection_)
292 coeffsi[1] -= 1 - wLin;
305 <<
"Could not fit face " << facei
306 <<
" Weights = " << coeffsi
307 <<
", reverting to linear." <<
nl
308 <<
" Linear weights " << wLin <<
" " << 1 - wLin <<
endl;
316template<
class FitDataType,
class ExtendedStencil,
class Polynomial>
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5)
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
Graphite solid properties.
bool movePoints()
Recalculate weights (but not stencil) when the mesh moves.
scalar linearLimitFactor() const
Factor the fit is allowed to deviate from the base scheme.
void calcFit(scalarList &coeffsi, const List< point > &, const scalar wLin, const label faci)
Calculate the fit for the specified face and set the coefficients.
const ExtendedStencil & stencil() const
Return reference to the stencil.
bool linearCorrection() const
FitData(const fvMesh &mesh, const ExtendedStencil &stencil, const bool linearCorrection, const scalar linearLimitFactor, const scalar centralWeight)
Construct from components.
scalar centralWeight() const
Return weight for central stencil.
virtual void calcFit()=0
Calculate the fit for all the faces.
void findFaceDirs(vector &idir, vector &jdir, vector &kdir, const label faci)
Find the normal direction (i) and j and k directions for face faci.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
void setSize(label n)
Alias for resize().
const fvMesh & mesh() const noexcept
Polynomial templated on size (order):
Singular value decomposition of a rectangular matrix.
scalarRectangularMatrix VSinvUt() const
Return the matrix product V S^(-1) U^T (the pseudo inverse).
const Cmpt & x() const noexcept
Access to the vector x component.
const Cmpt & z() const noexcept
Access to the vector z component.
const Cmpt & y() const noexcept
Access to the vector y component.
A face is a list of labels corresponding to mesh vertices.
Mesh data needed to do the Finite Volume discretisation.
virtual const faceList & faces() const
Return raw faces.
virtual const pointField & points() const
Return raw points.
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
const vectorField & faceCentres() const
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const volScalarField & p0
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
#define WarningInFunction
Report a warning using Foam::Warning.
RectangularMatrix< scalar > scalarRectangularMatrix
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
vector point
Point is a vector.
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
errorManipArg< error, int > exit(error &err, const int errNo=1)
List< scalar > scalarList
List of scalar.
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.