36template<
class Polynomial>
54 coeffs_(
mesh.nFaces())
57 <<
"Constructing CentredFitSnGradData<Polynomial>" <<
nl;
67template<
class Polynomial>
73 const scalar deltaCoeff,
80 this->findFaceDirs(idir, jdir, kdir, facei);
84 wts[0] = this->centralWeight();
85 wts[1] = this->centralWeight();
116 Polynomial::addCoeffs(
B[ip], d, wts[ip], this->dim());
120 for (label i = 0; i <
B.m(); i++)
127 label stencilSize =
C.size();
130 bool goodFit =
false;
131 for (
int iIt = 0; iIt < 8 && !goodFit; iIt++)
136 for (label i=0; i<stencilSize; i++)
138 coeffsi[i] = wts[1]*wts[i]*invB(1, i)/scale;
143 mag(wts[0]*wts[0]*invB(0, 0) - wLin)
144 < this->linearLimitFactor()*wLin)
145 && (
mag(wts[0]*wts[1]*invB(0, 1) - (1 - wLin)
146 ) < this->linearLimitFactor()*(1 - wLin))
147 && coeffsi[0] < 0 && coeffsi[1] > 0
148 &&
mag(coeffsi[0] + deltaCoeff) < 0.5*deltaCoeff
149 &&
mag(coeffsi[1] - deltaCoeff) < 0.5*deltaCoeff;
157 <<
"Cannot fit face " << facei <<
" iteration " << iIt
158 <<
" with sum of weights " <<
sum(coeffsi) <<
nl
159 <<
" Weights " << coeffsi <<
nl
160 <<
" Linear weights " << wLin <<
" " << 1 - wLin <<
nl
161 <<
" deltaCoeff " << deltaCoeff <<
nl
162 <<
" sing vals " << svd.S() <<
nl
163 <<
"Components of goodFit:\n"
164 <<
" wts[0]*wts[0]*invB(0, 0) = "
165 << wts[0]*wts[0]*invB(0, 0) <<
nl
166 <<
" wts[0]*wts[1]*invB(0, 1) = "
167 << wts[0]*wts[1]*invB(0, 1)
168 <<
" dim = " << this->dim() <<
endl;
173 for (label j = 0; j <
B.n(); j++)
179 for (label i = 0; i <
B.m(); i++)
190 coeffsi[0] += deltaCoeff;
191 coeffsi[1] -= deltaCoeff;
196 <<
"Could not fit face " << facei
197 <<
" Coefficients = " << coeffsi
198 <<
", reverting to uncorrected." <<
endl;
205template<
class Polynomial>
213 List<List<point>> stencilPoints(
mesh.nFaces());
214 this->stencil().collectData(
mesh.C(), stencilPoints);
221 for (label facei = 0; facei <
mesh.nInternalFaces(); facei++)
226 stencilPoints[facei],
243 label facei = pw.patch().start();
250 stencilPoints[facei],
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
Graphite solid properties.
void calcFit(scalarList &coeffsi, const List< point > &, const scalar wLin, const scalar deltaCoeff, const label faci)
Calculate the fit for the specified face and set the coefficients.
CentredFitSnGradData(const fvMesh &mesh, const extendedCentredCellToFaceStencil &stencil, const scalar linearLimitFactor, const scalar centralWeight)
Construct from components.
scalar linearLimitFactor() const
const extendedCentredCellToFaceStencil & stencil() const
FitData(const fvMesh &mesh, const extendedCentredCellToFaceStencil &stencil, const bool linearCorrection, const scalar linearLimitFactor, const scalar centralWeight)
scalar centralWeight() const
void findFaceDirs(vector &idir, vector &jdir, vector &kdir, const label faci)
GeometricBoundaryField< scalar, fvsPatchField, surfaceMesh > Boundary
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
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 scalarDiagonalMatrix & S() const noexcept
Return the singular values.
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.
Mesh data needed to do the Finite Volume discretisation.
label start() const noexcept
The patch start within the polyMesh face list.
virtual bool coupled() const
True if the patch field is coupled.
const fvPatch & patch() const noexcept
Return the patch.
const vectorField & faceCentres() const
const volScalarField & p0
#define DebugInfo
Report an information message using Foam::Info.
#define WarningInFunction
Report a warning using Foam::Warning.
#define DebugInFunction
Report an information message using Foam::Info.
RectangularMatrix< scalar > scalarRectangularMatrix
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
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)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
List< scalar > scalarList
List of scalar.
fvsPatchField< scalar > fvsPatchScalarField
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.