leastSquaresFaVectors.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) 2016-2017 Wikki Ltd
9  Copyright (C) 2020-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
29 #include "leastSquaresFaVectors.H"
30 #include "edgeFields.H"
31 #include "areaFields.H"
32 #include "mapPolyMesh.H"
33 #include "demandDrivenData.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
46 :
48  pVectorsPtr_(nullptr),
49  nVectorsPtr_(nullptr)
50 {}
51 
52 
53 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
54 
56 {
57  deleteDemandDrivenData(pVectorsPtr_);
58  deleteDemandDrivenData(nVectorsPtr_);
59 }
60 
61 
62 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
63 
64 void Foam::leastSquaresFaVectors::makeLeastSquaresVectors() const
65 {
67  << "Constructing finite area least square gradient vectors" << nl;
68 
69  pVectorsPtr_ = new edgeVectorField
70  (
71  IOobject
72  (
73  "LeastSquaresP",
74  mesh().pointsInstance(),
75  mesh().thisDb(),
78  false
79  ),
80  mesh(),
82  );
83  edgeVectorField& lsP = *pVectorsPtr_;
84 
85  nVectorsPtr_ = new edgeVectorField
86  (
87  IOobject
88  (
89  "LeastSquaresN",
90  mesh().pointsInstance(),
91  mesh().thisDb(),
94  false
95  ),
96  mesh(),
98  );
99  edgeVectorField& lsN = *nVectorsPtr_;
100 
101  // Set local references to mesh data
102  const labelUList& owner = mesh().owner();
103  const labelUList& neighbour = mesh().neighbour();
104 
105  const areaVectorField& C = mesh().areaCentres();
106  const edgeScalarField& w = mesh().weights();
107 
108 
109  // Set up temporary storage for the dd tensor (before inversion)
110  symmTensorField dd(mesh().nFaces(), Zero);
111 
112  forAll(owner, facei)
113  {
114  label own = owner[facei];
115  label nei = neighbour[facei];
116 
117  vector d = C[nei] - C[own];
118 
119  // Do not allow any mag(val) < SMALL
120  if (mag(d) < SMALL)
121  {
122  d = vector::uniform(SMALL);
123  }
124 
125  symmTensor wdd = (1.0/magSqr(d))*sqr(d);
126 
127  dd[own] += wdd;
128  dd[nei] += wdd;
129  }
130 
131 
132  forAll(lsP.boundaryField(), patchi)
133  {
134  const faePatchScalarField& pw = w.boundaryField()[patchi];
135 
136  const faPatch& p = pw.patch();
137  const labelUList& edgeFaces = p.edgeFaces();
138 
139  // Build the d-vectors
140  // HJ, reconsider deltas at the boundary, consistent with FVM
141  // Current implementation is good for fixedValue boundaries, but may
142  // cause problems with fixedGradient. HJ, 4/Oct/2010
143  const vectorField pd(p.delta());
144 
145  if (p.coupled())
146  {
147  forAll(pd, patchFacei)
148  {
149  const vector& d = pd[patchFacei];
150 
151  dd[edgeFaces[patchFacei]] +=
152  (1.0/magSqr(d))*sqr(d);
153  }
154  }
155  else
156  {
157  forAll(pd, patchFacei)
158  {
159  const vector& d = pd[patchFacei];
160 
161  dd[edgeFaces[patchFacei]] +=
162  (1.0/magSqr(d))*sqr(d);
163  }
164  }
165  }
166 
167 
168  // Invert the dd tensor
169  const symmTensorField invDd(inv(dd));
170 
171 
172  // Revisit all faces and calculate the lsP and lsN vectors
173  forAll(owner, facei)
174  {
175  label own = owner[facei];
176  label nei = neighbour[facei];
177 
178  vector d = C[nei] - C[own];
179 
180  // Do not allow any mag(val) < SMALL
181  if (mag(d) < SMALL)
182  {
183  d = vector::uniform(SMALL);
184  }
185 
186  scalar magSfByMagSqrd = 1.0/magSqr(d);
187 
188  lsP[facei] = magSfByMagSqrd*(invDd[own] & d);
189  lsN[facei] = -magSfByMagSqrd*(invDd[nei] & d);
190  }
191 
192  forAll(lsP.boundaryField(), patchi)
193  {
194  faePatchVectorField& patchLsP = lsP.boundaryFieldRef()[patchi];
195 
196  const faePatchScalarField& pw = w.boundaryField()[patchi];
197 
198  const faPatch& p = pw.patch();
199  const labelUList& edgeFaces = p.edgeFaces();
200 
201  // Build the d-vectors
202  const vectorField pd(p.delta());
203 
204  if (p.coupled())
205  {
206  forAll(pd, patchFacei)
207  {
208  const vector& d = pd[patchFacei];
209 
210  patchLsP[patchFacei] =
211  (1.0/magSqr(d))
212  *(invDd[edgeFaces[patchFacei]] & d);
213  }
214  }
215  else
216  {
217  forAll(pd, patchFacei)
218  {
219  const vector& d = pd[patchFacei];
220 
221  patchLsP[patchFacei] =
222  (1.0/magSqr(d))
223  *(invDd[edgeFaces[patchFacei]] & d);
224  }
225  }
226  }
227 
229  << "Done constructing finite area least square gradient vectors" << nl;
230 }
231 
232 
234 {
235  if (!pVectorsPtr_)
236  {
237  makeLeastSquaresVectors();
238  }
239 
240  return *pVectorsPtr_;
241 }
242 
243 
245 {
246  if (!nVectorsPtr_)
247  {
248  makeLeastSquaresVectors();
249  }
250 
251  return *nVectorsPtr_;
252 }
253 
254 
256 {
258  << "Clearing least square data" << nl;
259 
260  deleteDemandDrivenData(pVectorsPtr_);
261  deleteDemandDrivenData(nVectorsPtr_);
262 
263  return true;
264 }
265 
266 
267 // ************************************************************************* //
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition: faMesh.H:88
virtual bool movePoints()
Delete the least square vectors when the mesh moves.
Field< symmTensor > symmTensorField
Specialisation of Field<T> for symmTensor.
const edgeVectorField & nVectors() const
Return reference to neighbour least square vectors.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
leastSquaresFaVectors(const faMesh &)
Construct given an faMesh.
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
UList< label > labelUList
A UList of labels.
Definition: UList.H:80
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:545
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:84
virtual const surfaceScalarField & weights() const
Return reference to linear difference weighting factors.
dynamicFvMesh & mesh
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
Definition: symmTensor.H:55
#define DebugInFunction
Report an information message using Foam::Info.
GeometricField< vector, faePatchField, edgeMesh > edgeVectorField
Definition: edgeFieldsFwd.H:58
Vector< scalar > vector
Definition: vector.H:57
GeometricField< scalar, faePatchField, edgeMesh > edgeScalarField
Definition: edgeFieldsFwd.H:46
defineTypeNameAndDebug(combustionModel, 0)
volScalarField & C
faePatchField< scalar > faePatchScalarField
const labelUList & owner() const
Internal face owner. Note bypassing virtual mechanism so.
Definition: fvMesh.H:537
virtual ~leastSquaresFaVectors()
Destructor.
Template functions to aid in the implementation of demand driven data.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Nothing to be read.
const edgeVectorField & pVectors() const
Return reference to owner least square vectors.
faePatchField< vector > faePatchVectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
volScalarField & p
void deleteDemandDrivenData(DataPtr &dataPtr)
Least-squares gradient scheme vectors for the Finite Area method.
GeometricField< vector, faPatchField, areaMesh > areaVectorField
Definition: areaFieldsFwd.H:81
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157