unweightedLeastSquaresVectors.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) 2013-2016 OpenFOAM Foundation
9  Copyright (C) 2020-2021 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 "leastSquaresVectors.H"
30 #include "volFields.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(leastSquaresVectors, 0);
37 }
38 
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
43 :
44  MeshObject<fvMesh, Foam::MoveableMeshObject, leastSquaresVectors>(mesh),
45  pVectors_
46  (
47  IOobject
48  (
49  "LeastSquaresP",
50  mesh_.pointsInstance(),
51  mesh_,
52  IOobject::NO_READ,
53  IOobject::NO_WRITE,
54  false
55  ),
56  mesh_,
58  ),
59  nVectors_
60  (
61  IOobject
62  (
63  "LeastSquaresN",
64  mesh_.pointsInstance(),
65  mesh_,
66  IOobject::NO_READ,
67  IOobject::NO_WRITE,
68  false
69  ),
70  mesh_,
72  )
73 {
74  calcLeastSquaresVectors();
75 }
76 
77 
78 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
79 
81 {}
82 
83 
84 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
85 
86 void Foam::leastSquaresVectors::calcLeastSquaresVectors()
87 {
88  DebugInFunction << "Calculating least square gradient vectors" << nl;
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)
99  symmTensorField dd(mesh_.nCells(), Zero);
100 
101  forAll(owner, facei)
102  {
103  const label own = owner[facei];
104  const label nei = neighbour[facei];
105 
106  const 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  const vectorField pdSqr(sqr(p.delta()));
124 
125  forAll(pd, patchFacei)
126  {
127  dd[faceCells[patchFacei]] += pdSqr[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  const label own = owner[facei];
140  const label nei = neighbour[facei];
141 
142  const 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  const vectorField pd(p.delta());
157 
158  forAll(pd, patchFacei)
159  {
160  patchLsP[patchFacei] =
161  (invDd[faceCells[patchFacei]] & pd[patchFacei]);
162  }
163  }
164 
165  DebugInfo
166  << "Finished calculating least square gradient vectors" << endl;
167 }
168 
169 
171 {
172  calcLeastSquaresVectors();
173  return true;
174 }
175 
176 
177 // ************************************************************************* //
volFields.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::leastSquaresVectors::~leastSquaresVectors
virtual ~leastSquaresVectors()
Destructor.
Definition: invDistLeastSquaresVectors.C:80
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
leastSquaresVectors.H
C
volScalarField & C
Definition: readThermalProperties.H:102
Foam::leastSquaresVectors::movePoints
virtual bool movePoints()
Delete the least square vectors when the mesh moves.
Definition: invDistLeastSquaresVectors.C:173
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
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::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
Foam::fvsPatchVectorField
fvsPatchField< vector > fvsPatchVectorField
Definition: fvsPatchFieldsFwd.H:48
Foam::leastSquaresVectors::leastSquaresVectors
leastSquaresVectors(const fvMesh &)
Construct given an fvMesh.
Definition: invDistLeastSquaresVectors.C:42
Foam::inv
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:73
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:388
Foam::symmTensor
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
Definition: symmTensor.H:59
Foam::fvMesh::C
const volVectorField & C() const
Return cell centres as volVectorField.
Definition: fvMeshGeometry.C:341
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::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:62
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:382
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::labelUList
UList< label > labelUList
A UList of labels.
Definition: UList.H:85
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189