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 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 :
45  MeshObject_type(mesh)
46 {}
47 
48 
49 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
50 
51 void Foam::leastSquaresFaVectors::makeLeastSquaresVectors() const
52 {
54  << "Constructing finite area (invDist) least square gradient vectors"
55  << nl;
56 
57  pVectorsPtr_ = std::make_unique<edgeVectorField>
58  (
59  IOobject
60  (
61  "LeastSquaresP",
62  mesh().pointsInstance(),
63  mesh().thisDb(),
67  ),
68  mesh(),
70  );
71  auto& lsP = *pVectorsPtr_;
72 
73  nVectorsPtr_ = std::make_unique<edgeVectorField>
74  (
75  IOobject
76  (
77  "LeastSquaresN",
78  mesh().pointsInstance(),
79  mesh().thisDb(),
83  ),
84  mesh(),
86  );
87  auto& lsN = *nVectorsPtr_;
88 
89  // Set local references to mesh data
90  const labelUList& owner = mesh().owner();
91  const labelUList& neighbour = mesh().neighbour();
92 
93  const areaVectorField& C = mesh().areaCentres();
94 
95  // Set up temporary storage for the dd tensor (before inversion)
96  symmTensorField dd(mesh().nFaces(), Zero);
97 
98  // No contribution when mag(val) < SMALL
99  const scalar minLenSqr(SMALL*SMALL);
100 
101  forAll(owner, facei)
102  {
103  const label own = owner[facei];
104  const label nei = neighbour[facei];
105 
106  const vector d(C[nei] - C[own]);
107 
108  // No contribution when mag(val) < SMALL
109  const scalar magSqrDist = d.magSqr();
110  if (magSqrDist >= minLenSqr)
111  {
112  const symmTensor wdd(sqr(d)/magSqrDist);
113  dd[own] += wdd;
114  dd[nei] += wdd;
115  }
116  }
117 
118 
119  for (const auto& patchLsP : lsP.boundaryField())
120  {
121  const faPatch& p = patchLsP.patch();
122  const labelUList& edgeFaces = p.edgeFaces();
123 
124  // Build the d-vectors
125  const vectorField pd(p.delta());
126 
127  forAll(pd, patchFacei)
128  {
129  const vector& d = pd[patchFacei];
130 
131  // No contribution when mag(val) < SMALL
132  const scalar magSqrDist = d.magSqr();
133  if (magSqrDist >= minLenSqr)
134  {
135  dd[edgeFaces[patchFacei]] += sqr(d)/magSqrDist;
136  }
137  }
138  }
139 
140 
141  // Invert the dd tensors - including failsafe checks
142  const symmTensorField invDd(inv(dd));
143 
144 
145  // Revisit all faces and calculate the lsP and lsN vectors
146  forAll(owner, facei)
147  {
148  const label own = owner[facei];
149  const label nei = neighbour[facei];
150 
151  const vector d(C[nei] - C[own]);
152 
153  // No contribution when mag(val) < SMALL
154  const scalar magSqrDist = d.magSqr();
155  if (magSqrDist >= minLenSqr)
156  {
157  lsP[facei] = (invDd[own] & d)/magSqrDist;
158  lsN[facei] = -(invDd[nei] & d)/magSqrDist;
159  }
160  // ... already zero from initialisation
161  // else
162  // {
163  // lsP[facei] = Zero;
164  // lsN[facei] = Zero;
165  // }
166  }
167 
168  for (auto& patchLsP : lsP.boundaryFieldRef())
169  {
170  const faPatch& p = patchLsP.patch();
171  const labelUList& edgeFaces = p.edgeFaces();
172 
173  // Build the d-vectors
174  const vectorField pd(p.delta());
175 
176  forAll(pd, patchFacei)
177  {
178  const vector& d = pd[patchFacei];
179 
180  // No contribution when mag(val) < SMALL
181  const scalar magSqrDist = d.magSqr();
182  if (magSqrDist >= minLenSqr)
183  {
184  patchLsP[patchFacei] =
185  (invDd[edgeFaces[patchFacei]] & d)/magSqrDist;
186  }
187  // ... already zero from initialisation
188  // else
189  // {
190  // patchLsP[patchFacei] = Zero;
191  // }
192  }
193  }
194 
195  DebugInfo
196  << "Done constructing finite area least square gradient vectors" << nl;
197 }
198 
199 
201 {
202  if (!pVectorsPtr_)
203  {
204  makeLeastSquaresVectors();
205  }
206 
207  return *pVectorsPtr_;
208 }
209 
210 
212 {
213  if (!nVectorsPtr_)
214  {
215  makeLeastSquaresVectors();
216  }
217 
218  return *nVectorsPtr_;
219 }
220 
221 
223 {
225  << "Clearing least square data" << nl;
226 
227  pVectorsPtr_.reset(nullptr);
228  nVectorsPtr_.reset(nullptr);
229 
230  return true;
231 }
232 
233 
234 // ************************************************************************* //
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition: faMesh.H:133
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.
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
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.
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
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
Least-squares gradient scheme vectors for the Finite Area method.
GeometricField< vector, faPatchField, areaMesh > areaVectorField
Definition: areaFieldsFwd.H:76
Do not request registration (bool: false)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127