Loading...
Searching...
No Matches
CentredFitSnGradData.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2013-2016 OpenFOAM Foundation
9 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
31#include "SVD.H"
33
34// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35
36template<class Polynomial>
38(
39 const fvMesh& mesh,
41 const scalar linearLimitFactor,
42 const scalar centralWeight
43)
44:
46 <
50 >
51 (
53 ),
54 coeffs_(mesh.nFaces())
55{
57 << "Constructing CentredFitSnGradData<Polynomial>" << nl;
58
59 calcFit();
60
61 DebugInfo << " Finished constructing polynomialFit data" << endl;
62}
63
64
65// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
66
67template<class Polynomial>
69(
70 scalarList& coeffsi,
71 const List<point>& C,
72 const scalar wLin,
73 const scalar deltaCoeff,
74 const label facei
75)
76{
77 vector idir(1,0,0);
78 vector jdir(0,1,0);
79 vector kdir(0,0,1);
80 this->findFaceDirs(idir, jdir, kdir, facei);
81
82 // Setup the point weights
83 scalarList wts(C.size(), scalar(1));
84 wts[0] = this->centralWeight();
85 wts[1] = this->centralWeight();
86
87 // Reference point
88 point p0 = this->mesh().faceCentres()[facei];
89
90 // p0 -> p vector in the face-local coordinate system
91 vector d;
92
93 // Local coordinate scaling
94 scalar scale = 1;
95
96 // Matrix of the polynomial components
97 scalarRectangularMatrix B(C.size(), this->minSize(), Zero);
98
99 forAll(C, ip)
100 {
101 const point& p = C[ip];
102 const vector p0p = p - p0;
103
104 d.x() = p0p & idir;
105 d.y() = p0p & jdir;
106 d.z() = p0p & kdir;
107
108 if (ip == 0)
109 {
110 scale = cmptMax(cmptMag((d)));
111 }
112
113 // Scale the radius vector
114 d /= scale;
115
116 Polynomial::addCoeffs(B[ip], d, wts[ip], this->dim());
117 }
118
119 // Additional weighting for constant and linear terms
120 for (label i = 0; i < B.m(); i++)
121 {
122 B(i, 0) *= wts[0];
123 B(i, 1) *= wts[0];
124 }
125
126 // Set the fit
127 label stencilSize = C.size();
128 coeffsi.setSize(stencilSize);
129
130 bool goodFit = false;
131 for (int iIt = 0; iIt < 8 && !goodFit; iIt++)
132 {
133 SVD svd(B, SMALL);
134 scalarRectangularMatrix invB(svd.VSinvUt());
135
136 for (label i=0; i<stencilSize; i++)
137 {
138 coeffsi[i] = wts[1]*wts[i]*invB(1, i)/scale;
139 }
140
141 goodFit =
142 (
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;
150
151 if (!goodFit)
152 {
153 // (not good fit so increase weight in the centre and weight
154 // for constant and linear terms)
155
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;
169
170 wts[0] *= 10;
171 wts[1] *= 10;
172
173 for (label j = 0; j < B.n(); j++)
174 {
175 B(0, j) *= 10;
176 B(1, j) *= 10;
177 }
178
179 for (label i = 0; i < B.m(); i++)
180 {
181 B(i, 0) *= 10;
182 B(i, 1) *= 10;
183 }
184 }
185 }
186
187 if (goodFit)
188 {
189 // Remove the uncorrected coefficients
190 coeffsi[0] += deltaCoeff;
191 coeffsi[1] -= deltaCoeff;
192 }
193 else
194 {
196 << "Could not fit face " << facei
197 << " Coefficients = " << coeffsi
198 << ", reverting to uncorrected." << endl;
200 coeffsi = 0;
201 }
202}
203
204
205template<class Polynomial>
207{
208 const fvMesh& mesh = this->mesh();
209
210 // Get the cell/face centres in stencil order.
211 // Centred face stencils no good for triangles or tets.
212 // Need bigger stencils
213 List<List<point>> stencilPoints(mesh.nFaces());
214 this->stencil().collectData(mesh.C(), stencilPoints);
215
216 // find the fit coefficients for every face in the mesh
217
218 const surfaceScalarField& w = mesh.surfaceInterpolation::weights();
219 const surfaceScalarField& dC = mesh.nonOrthDeltaCoeffs();
220
221 for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
222 {
223 calcFit
224 (
225 coeffs_[facei],
226 stencilPoints[facei],
227 w[facei],
228 dC[facei],
229 facei
230 );
231 }
232
233 const surfaceScalarField::Boundary& bw = w.boundaryField();
234 const surfaceScalarField::Boundary& bdC = dC.boundaryField();
235
236 forAll(bw, patchi)
237 {
238 const fvsPatchScalarField& pw = bw[patchi];
239 const fvsPatchScalarField& pdC = bdC[patchi];
240
241 if (pw.coupled())
242 {
243 label facei = pw.patch().start();
244
245 forAll(pw, i)
246 {
247 calcFit
248 (
249 coeffs_[facei],
250 stencilPoints[facei],
251 pw[i],
252 pdC[i],
253 facei
254 );
255 facei++;
256 }
257 }
258 }
259}
260
261
262// ************************************************************************* //
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
Graphite solid properties.
Definition C.H:49
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.
FitData(const fvMesh &mesh, const extendedCentredCellToFaceStencil &stencil, const bool linearCorrection, const scalar linearLimitFactor, const scalar centralWeight)
Definition FitData.C:30
void findFaceDirs(vector &idir, vector &jdir, vector &kdir, const label faci)
Definition FitData.C:65
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...
Definition List.H:72
void setSize(label n)
Alias for resize().
Definition List.H:536
Polynomial templated on size (order):
Definition Polynomial.H:74
Singular value decomposition of a rectangular matrix.
Definition SVD.H:51
scalarRectangularMatrix VSinvUt() const
Return the matrix product V S^(-1) U^T (the pseudo inverse).
Definition SVD.C:479
const scalarDiagonalMatrix & S() const noexcept
Return the singular values.
Definition SVD.H:122
const Cmpt & x() const noexcept
Access to the vector x component.
Definition Vector.H:135
const Cmpt & z() const noexcept
Access to the vector z component.
Definition Vector.H:145
const Cmpt & y() const noexcept
Access to the vector y component.
Definition Vector.H:140
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
label start() const noexcept
The patch start within the polyMesh face list.
Definition fvPatch.H:226
virtual bool coupled() const
True if the patch field is coupled.
const fvPatch & patch() const noexcept
Return the patch.
const vectorField & faceCentres() const
volScalarField & p
const volScalarField & p0
Definition EEqn.H:36
dynamicFvMesh & mesh
#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.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
vector point
Point is a vector.
Definition point.H:37
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).
Definition zero.H:127
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
Vector< scalar > vector
Definition vector.H:57
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
fvsPatchField< scalar > fvsPatchScalarField
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.