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-2023 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 (invDist) least square gradient vectors"
68  << nl;
69 
70  pVectorsPtr_ = new edgeVectorField
71  (
72  IOobject
73  (
74  "LeastSquaresP",
75  mesh().pointsInstance(),
76  mesh().thisDb(),
80  ),
81  mesh(),
83  );
84  auto& lsP = *pVectorsPtr_;
85 
86  nVectorsPtr_ = new edgeVectorField
87  (
88  IOobject
89  (
90  "LeastSquaresN",
91  mesh().pointsInstance(),
92  mesh().thisDb(),
96  ),
97  mesh(),
99  );
100  auto& lsN = *nVectorsPtr_;
101 
102  // Set local references to mesh data
103  const labelUList& owner = mesh().owner();
104  const labelUList& neighbour = mesh().neighbour();
105 
106  const areaVectorField& C = mesh().areaCentres();
107 
108  // Set up temporary storage for the dd tensor (before inversion)
109  symmTensorField dd(mesh().nFaces(), Zero);
110 
111  // No contribution when mag(val) < SMALL
112  const scalar minLenSqr(SMALL*SMALL);
113 
114  forAll(owner, facei)
115  {
116  const label own = owner[facei];
117  const label nei = neighbour[facei];
118 
119  const vector d(C[nei] - C[own]);
120 
121  // No contribution when mag(val) < SMALL
122  const scalar magSqrDist = d.magSqr();
123  if (magSqrDist >= minLenSqr)
124  {
125  const symmTensor wdd(sqr(d)/magSqrDist);
126  dd[own] += wdd;
127  dd[nei] += wdd;
128  }
129  }
130 
131 
132  for (const auto& patchLsP : lsP.boundaryField())
133  {
134  const faPatch& p = patchLsP.patch();
135  const labelUList& edgeFaces = p.edgeFaces();
136 
137  // Build the d-vectors
138  const vectorField pd(p.delta());
139 
140  forAll(pd, patchFacei)
141  {
142  const vector& d = pd[patchFacei];
143 
144  // No contribution when mag(val) < SMALL
145  const scalar magSqrDist = d.magSqr();
146  if (magSqrDist >= minLenSqr)
147  {
148  dd[edgeFaces[patchFacei]] += sqr(d)/magSqrDist;
149  }
150  }
151  }
152 
153 
154  // Invert the dd tensors - including failsafe checks
155  const symmTensorField invDd(inv(dd));
156 
157 
158  // Revisit all faces and calculate the lsP and lsN vectors
159  forAll(owner, facei)
160  {
161  const label own = owner[facei];
162  const label nei = neighbour[facei];
163 
164  const vector d(C[nei] - C[own]);
165 
166  // No contribution when mag(val) < SMALL
167  const scalar magSqrDist = d.magSqr();
168  if (magSqrDist >= minLenSqr)
169  {
170  lsP[facei] = (invDd[own] & d)/magSqrDist;
171  lsN[facei] = -(invDd[nei] & d)/magSqrDist;
172  }
173  // ... already zero from initialisation
174  // else
175  // {
176  // lsP[facei] = Zero;
177  // lsN[facei] = Zero;
178  // }
179  }
180 
181  for (auto& patchLsP : lsP.boundaryFieldRef())
182  {
183  const faPatch& p = patchLsP.patch();
184  const labelUList& edgeFaces = p.edgeFaces();
185 
186  // Build the d-vectors
187  const vectorField pd(p.delta());
188 
189  forAll(pd, patchFacei)
190  {
191  const vector& d = pd[patchFacei];
192 
193  // No contribution when mag(val) < SMALL
194  const scalar magSqrDist = d.magSqr();
195  if (magSqrDist >= minLenSqr)
196  {
197  patchLsP[patchFacei] =
198  (invDd[edgeFaces[patchFacei]] & d)/magSqrDist;
199  }
200  // ... already zero from initialisation
201  // else
202  // {
203  // patchLsP[patchFacei] = Zero;
204  // }
205  }
206  }
207 
208  DebugInfo
209  << "Done constructing finite area least square gradient vectors" << nl;
210 }
211 
212 
214 {
215  if (!pVectorsPtr_)
216  {
217  makeLeastSquaresVectors();
218  }
219 
220  return *pVectorsPtr_;
221 }
222 
223 
225 {
226  if (!nVectorsPtr_)
227  {
228  makeLeastSquaresVectors();
229  }
230 
231  return *nVectorsPtr_;
232 }
233 
234 
236 {
238  << "Clearing least square data" << nl;
239 
240  deleteDemandDrivenData(pVectorsPtr_);
241  deleteDemandDrivenData(nVectorsPtr_);
242 
243  return true;
244 }
245 
246 
247 // ************************************************************************* //
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition: faMesh.H:87
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.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
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:78
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:580
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:85
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
scalar magSqr() const
The length (L2-norm) squared of the vector.
Definition: VectorI.H:76
Vector< scalar > vector
Definition: vector.H:57
#define DebugInfo
Report an information message using Foam::Info.
defineTypeNameAndDebug(combustionModel, 0)
volScalarField & C
const labelUList & owner() const
Internal face owner. Note bypassing virtual mechanism so.
Definition: fvMesh.H:572
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.
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
Do not request registration (bool: false)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127