invDistLeastSquaresVectors.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2013-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "leastSquaresVectors.H"
27 #include "volFields.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33  defineTypeNameAndDebug(leastSquaresVectors, 0);
34 }
35 
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
40 :
42  pVectors_
43  (
44  IOobject
45  (
46  "LeastSquaresP",
47  mesh_.pointsInstance(),
48  mesh_,
49  IOobject::NO_READ,
50  IOobject::NO_WRITE,
51  false
52  ),
53  mesh_,
55  ),
56  nVectors_
57  (
58  IOobject
59  (
60  "LeastSquaresN",
61  mesh_.pointsInstance(),
62  mesh_,
63  IOobject::NO_READ,
64  IOobject::NO_WRITE,
65  false
66  ),
67  mesh_,
69  )
70 {
71  calcLeastSquaresVectors();
72 }
73 
74 
75 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
76 
78 {}
79 
80 
81 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
82 
83 void Foam::leastSquaresVectors::calcLeastSquaresVectors()
84 {
85  if (debug)
86  {
87  InfoInFunction << "Calculating least square gradient vectors" << endl;
88  }
89 
90  const fvMesh& mesh = mesh_;
91 
92  // Set local references to mesh data
93  const labelUList& owner = mesh_.owner();
94  const labelUList& neighbour = mesh_.neighbour();
95 
96  const volVectorField& C = mesh.C();
97 
98  // Set up temporary storage for the dd tensor (before inversion)
100 
101  forAll(owner, facei)
102  {
103  label own = owner[facei];
104  label nei = neighbour[facei];
105 
106  vector d = C[nei] - C[own];
107  symmTensor wdd = sqr(d)/magSqr(d);
108  dd[own] += wdd;
109  dd[nei] += wdd;
110  }
111 
112 
114  pVectors_.boundaryField();
115 
116  forAll(blsP, patchi)
117  {
118  const fvsPatchVectorField& patchLsP = blsP[patchi];
119 
120  const fvPatch& p = patchLsP.patch();
121  const labelUList& faceCells = p.patch().faceCells();
122 
123  // Build the d-vectors
124  vectorField pd(p.delta());
125 
126  forAll(pd, patchFacei)
127  {
128  const vector& d = pd[patchFacei];
129 
130  dd[faceCells[patchFacei]] += sqr(d)/magSqr(d);
131  }
132  }
133 
134 
135  // Invert the dd tensor
136  const symmTensorField invDd(inv(dd));
137 
138 
139  // Revisit all faces and calculate the pVectors_ and nVectors_ vectors
140  forAll(owner, facei)
141  {
142  label own = owner[facei];
143  label nei = neighbour[facei];
144 
145  vector d = C[nei] - C[own];
146 
147  pVectors_[facei] = (invDd[own] & d)/magSqr(d);
148  nVectors_[facei] = -(invDd[nei] & d)/magSqr(d);
149  }
150 
151  forAll(blsP, patchi)
152  {
153  fvsPatchVectorField& patchLsP = blsP[patchi];
154 
155  const fvPatch& p = patchLsP.patch();
156  const labelUList& faceCells = p.faceCells();
157 
158  // Build the d-vectors
159  vectorField pd(p.delta());
160 
161  forAll(pd, patchFacei)
162  {
163  const vector& d = pd[patchFacei];
164 
165  patchLsP[patchFacei] = (invDd[faceCells[patchFacei]] & d)/magSqr(d);
166  }
167  }
168 
169  if (debug)
170  {
172  <<"Finished calculating least square gradient vectors" << endl;
173  }
174 }
175 
176 
178 {
179  calcLeastSquaresVectors();
180  return true;
181 }
182 
183 
184 // ************************************************************************* //
Graphite solid properties.
Definition: C.H:48
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const Boundary & boundaryField() const
Return const-reference to the boundary field.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
label nCells() const
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
leastSquaresVectors(const fvMesh &)
Construct given an fvMesh.
const dimensionSet dimless
Least-squares gradient scheme vectors.
fvMesh & mesh
const dimensionSet dimLength
virtual bool movePoints()
Delete the least square vectors when the mesh moves.
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:423
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:87
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:340
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:99
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:139
static const zero Zero
Definition: zero.H:97
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
dimensioned< scalar > magSqr(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:417
const fvPatch & patch() const
Return patch.
label patchi
virtual ~leastSquaresVectors()
Destructor.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
const volVectorField & C() const
Return cell centres.
volScalarField & p
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
virtual tmp< vectorField > delta() const
Return cell-centre to face-centre vector.
Definition: fvPatch.C:148
Namespace for OpenFOAM.
#define InfoInFunction
Report an information message using Foam::Info.