unweightedLeastSquaresVectors.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 :
41  MeshObject<fvMesh, Foam::MoveableMeshObject, leastSquaresVectors>(mesh),
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  symmTensor wdd = sqr(C[nei] - C[own]);
107  dd[own] += wdd;
108  dd[nei] += wdd;
109  }
110 
111 
112  surfaceVectorField::Boundary& blsP =
113  pVectors_.boundaryField();
114 
115  forAll(blsP, patchi)
116  {
117  const fvsPatchVectorField& patchLsP = blsP[patchi];
118 
119  const fvPatch& p = patchLsP.patch();
120  const labelUList& faceCells = p.patch().faceCells();
121 
122  // Build the d-vectors
123  vectorField pd(p.delta());
124 
125  forAll(pd, patchFacei)
126  {
127  dd[faceCells[patchFacei]] += sqr(pd[patchFacei]);
128  }
129  }
130 
131 
132  // Invert the dd tensor
133  const symmTensorField invDd(inv(dd));
134 
135 
136  // Revisit all faces and calculate the pVectors_ and nVectors_ vectors
137  forAll(owner, facei)
138  {
139  label own = owner[facei];
140  label nei = neighbour[facei];
141 
142  vector d = C[nei] - C[own];
143 
144  pVectors_[facei] = (invDd[own] & d);
145  nVectors_[facei] = -(invDd[nei] & d);
146  }
147 
148  forAll(blsP, patchi)
149  {
150  fvsPatchVectorField& patchLsP = blsP[patchi];
151 
152  const fvPatch& p = patchLsP.patch();
153  const labelUList& faceCells = p.faceCells();
154 
155  // Build the d-vectors
156  vectorField pd(p.delta());
157 
158  forAll(pd, patchFacei)
159  {
160  patchLsP[patchFacei] =
161  (invDd[faceCells[patchFacei]] & pd[patchFacei]);
162  }
163  }
164 
165  if (debug)
166  {
168  << "Finished calculating least square gradient vectors" << endl;
169  }
170 }
171 
172 
174 {
175  calcLeastSquaresVectors();
176  return true;
177 }
178 
179 
180 // ************************************************************************* //
fvsPatchField< vector > fvsPatchVectorField
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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:59
Field< symmTensor > symmTensorField
Specialisation of Field<T> for symmTensor.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
label nCells() const
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
leastSquaresVectors(const fvMesh &)
Construct given an fvMesh.
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
const dimensionSet dimless
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
UList< label > labelUList
Definition: UList.H:65
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:284
SymmTensor< scalar > symmTensor
SymmTensor of scalars.
Definition: symmTensor.H:48
static const zero Zero
Definition: zero.H:97
defineTypeNameAndDebug(combustionModel, 0)
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:278
label patchi
virtual ~leastSquaresVectors()
Destructor.
Field< vector > vectorField
Specialisation of Field<T> for vector.
Namespace for OpenFOAM.
#define InfoInFunction
Report an information message using Foam::Info.