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