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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "CentredFitSnGradData.H"
29 #include "surfaceFields.H"
30 #include "SVD.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
35 template<class Polynomial>
37 (
38  const fvMesh& mesh,
39  const extendedCentredCellToFaceStencil& stencil,
40  const scalar linearLimitFactor,
41  const scalar centralWeight
42 )
43 :
44  FitData
45  <
49  >
50  (
51  mesh, stencil, true, linearLimitFactor, centralWeight
52  ),
53  coeffs_(mesh.nFaces())
54 {
55  if (debug)
56  {
58  << "Contructing CentredFitSnGradData<Polynomial>" << endl;
59  }
60 
61  calcFit();
62 
63  if (debug)
64  {
65  Info<< " Finished constructing polynomialFit data" << endl;
66  }
67 }
68 
69 
70 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
71 
72 template<class Polynomial>
74 (
75  scalarList& coeffsi,
76  const List<point>& C,
77  const scalar wLin,
78  const scalar deltaCoeff,
79  const label facei
80 )
81 {
82  vector idir(1,0,0);
83  vector jdir(0,1,0);
84  vector kdir(0,0,1);
85  this->findFaceDirs(idir, jdir, kdir, facei);
86 
87  // Setup the point weights
88  scalarList wts(C.size(), scalar(1));
89  wts[0] = this->centralWeight();
90  wts[1] = this->centralWeight();
91 
92  // Reference point
93  point p0 = this->mesh().faceCentres()[facei];
94 
95  // p0 -> p vector in the face-local coordinate system
96  vector d;
97 
98  // Local coordinate scaling
99  scalar scale = 1;
100 
101  // Matrix of the polynomial components
102  scalarRectangularMatrix B(C.size(), this->minSize(), Zero);
103 
104  forAll(C, ip)
105  {
106  const point& p = C[ip];
107  const vector p0p = p - p0;
108 
109  d.x() = p0p & idir;
110  d.y() = p0p & jdir;
111  d.z() = p0p & kdir;
112 
113  if (ip == 0)
114  {
115  scale = cmptMax(cmptMag((d)));
116  }
117 
118  // Scale the radius vector
119  d /= scale;
120 
121  Polynomial::addCoeffs(B[ip], d, wts[ip], this->dim());
122  }
123 
124  // Additional weighting for constant and linear terms
125  for (label i = 0; i < B.m(); i++)
126  {
127  B(i, 0) *= wts[0];
128  B(i, 1) *= wts[0];
129  }
130 
131  // Set the fit
132  label stencilSize = C.size();
133  coeffsi.setSize(stencilSize);
134 
135  bool goodFit = false;
136  for (int iIt = 0; iIt < 8 && !goodFit; iIt++)
137  {
138  SVD svd(B, SMALL);
139  scalarRectangularMatrix invB(svd.VSinvUt());
140 
141  for (label i=0; i<stencilSize; i++)
142  {
143  coeffsi[i] = wts[1]*wts[i]*invB(1, i)/scale;
144  }
145 
146  goodFit =
147  (
148  mag(wts[0]*wts[0]*invB(0, 0) - wLin)
149  < this->linearLimitFactor()*wLin)
150  && (mag(wts[0]*wts[1]*invB(0, 1) - (1 - wLin)
151  ) < this->linearLimitFactor()*(1 - wLin))
152  && coeffsi[0] < 0 && coeffsi[1] > 0
153  && mag(coeffsi[0] + deltaCoeff) < 0.5*deltaCoeff
154  && mag(coeffsi[1] - deltaCoeff) < 0.5*deltaCoeff;
155 
156  if (!goodFit)
157  {
158  // (not good fit so increase weight in the centre and weight
159  // for constant and linear terms)
160 
162  << "Cannot fit face " << facei << " iteration " << iIt
163  << " with sum of weights " << sum(coeffsi) << nl
164  << " Weights " << coeffsi << nl
165  << " Linear weights " << wLin << " " << 1 - wLin << nl
166  << " deltaCoeff " << deltaCoeff << nl
167  << " sing vals " << svd.S() << nl
168  << "Components of goodFit:\n"
169  << " wts[0]*wts[0]*invB(0, 0) = "
170  << wts[0]*wts[0]*invB(0, 0) << nl
171  << " wts[0]*wts[1]*invB(0, 1) = "
172  << wts[0]*wts[1]*invB(0, 1)
173  << " dim = " << this->dim() << endl;
174 
175  wts[0] *= 10;
176  wts[1] *= 10;
177 
178  for (label j = 0; j < B.n(); j++)
179  {
180  B(0, j) *= 10;
181  B(1, j) *= 10;
182  }
183 
184  for (label i = 0; i < B.m(); i++)
185  {
186  B(i, 0) *= 10;
187  B(i, 1) *= 10;
188  }
189  }
190  }
191 
192  if (goodFit)
193  {
194  // Remove the uncorrected coefficients
195  coeffsi[0] += deltaCoeff;
196  coeffsi[1] -= deltaCoeff;
197  }
198  else
199  {
201  << "Could not fit face " << facei
202  << " Coefficients = " << coeffsi
203  << ", reverting to uncorrected." << endl;
204 
205  coeffsi = 0;
206  }
207 }
208 
209 
210 template<class Polynomial>
212 {
213  const fvMesh& mesh = this->mesh();
214 
215  // Get the cell/face centres in stencil order.
216  // Centred face stencils no good for triangles or tets.
217  // Need bigger stencils
218  List<List<point>> stencilPoints(mesh.nFaces());
219  this->stencil().collectData(mesh.C(), stencilPoints);
220 
221  // find the fit coefficients for every face in the mesh
222 
223  const surfaceScalarField& w = mesh.surfaceInterpolation::weights();
224  const surfaceScalarField& dC = mesh.nonOrthDeltaCoeffs();
225 
226  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
227  {
228  calcFit
229  (
230  coeffs_[facei],
231  stencilPoints[facei],
232  w[facei],
233  dC[facei],
234  facei
235  );
236  }
237 
238  const surfaceScalarField::Boundary& bw = w.boundaryField();
239  const surfaceScalarField::Boundary& bdC = dC.boundaryField();
240 
241  forAll(bw, patchi)
242  {
243  const fvsPatchScalarField& pw = bw[patchi];
244  const fvsPatchScalarField& pdC = bdC[patchi];
245 
246  if (pw.coupled())
247  {
248  label facei = pw.patch().start();
249 
250  forAll(pw, i)
251  {
252  calcFit
253  (
254  coeffs_[facei],
255  stencilPoints[facei],
256  pw[i],
257  pdC[i],
258  facei
259  );
260  facei++;
261  }
262  }
263  }
264 }
265 
266 
267 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::Vector::x
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:78
InfoInFunction
#define InfoInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:316
Foam::FitData
Data for the upwinded and centred polynomial fit interpolation schemes. The linearCorrection_ determi...
Definition: FitData.H:57
p
volScalarField & p
Definition: createFieldRefs.H:8
SVD.H
Foam::extendedCentredCellToFaceStencil
Definition: extendedCentredCellToFaceStencil.H:51
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::fvPatch::start
virtual label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:157
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:68
B
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
surfaceFields.H
Foam::surfaceFields.
Foam::Vector::z
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:90
Foam::CentredFitSnGradData::CentredFitSnGradData
CentredFitSnGradData(const fvMesh &mesh, const extendedCentredCellToFaceStencil &stencil, const scalar linearLimitFactor, const scalar centralWeight)
Construct from components.
Definition: CentredFitSnGradData.C:37
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
CentredFitSnGradData.H
Foam::cmptMag
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:400
Foam::cmptMax
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:253
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::SVD::VSinvUt
scalarRectangularMatrix VSinvUt() const
Return the matrix product V S^(-1) U^T (the pseudo inverse)
Definition: SVD.C:387
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::RectangularMatrix< scalar >
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
extendedCentredCellToFaceStencil.H
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam::CentredFitSnGradData
Data for centred fit snGrad schemes.
Definition: CentredFitSnGradData.H:54
Foam::SVD::S
const scalarDiagonalMatrix & S() const
Return the singular values.
Definition: SVDI.H:52
Foam::Vector::y
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:84
Foam::Polynomial
Polynomial templated on size (order):
Definition: Polynomial.H:67
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::fvsPatchField::patch
const fvPatch & patch() const
Return patch.
Definition: fvsPatchField.H:281
Foam::Vector< scalar >
Foam::List< scalar >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::SVD
Singular value decomposition of a rectangular matrix.
Definition: SVD.H:53
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
p0
const volScalarField & p0
Definition: EEqn.H:36
Foam::C
Graphite solid properties.
Definition: C.H:50
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::GeometricField< scalar, fvsPatchField, surfaceMesh >
Foam::CentredFitSnGradData::calcFit
void calcFit()
Calculate the fit for all the faces.
Definition: CentredFitSnGradData.C:211
Foam::fvsPatchField::coupled
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvsPatchField.H:310
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62