leastSquaresVectors.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) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2020-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
34namespace 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
86void 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 const surfaceScalarField& w = mesh.weights();
98 const surfaceScalarField& magSf = mesh.magSf();
99
100
101 // Set up temporary storage for the dd tensor (before inversion)
102 symmTensorField dd(mesh_.nCells(), Zero);
103
104 forAll(owner, facei)
105 {
106 const label own = owner[facei];
107 const label nei = neighbour[facei];
108
109 const vector d(C[nei] - C[own]);
110 const symmTensor wdd((magSf[facei]/magSqr(d))*sqr(d));
111
112 dd[own] += (1.0 - w[facei])*wdd;
113 dd[nei] += w[facei]*wdd;
114 }
115
116
117 surfaceVectorField::Boundary& pVectorsBf =
118 pVectors_.boundaryFieldRef();
119
120 forAll(pVectorsBf, patchi)
121 {
122 const fvsPatchScalarField& pw = w.boundaryField()[patchi];
123 const fvsPatchScalarField& pMagSf = magSf.boundaryField()[patchi];
124
125 const fvPatch& p = pw.patch();
126 const labelUList& faceCells = p.patch().faceCells();
127
128 // Build the d-vectors
129 const vectorField pd(p.delta());
130
131 if (pw.coupled())
132 {
133 forAll(pd, patchFacei)
134 {
135 const vector& d = pd[patchFacei];
136
137 dd[faceCells[patchFacei]] +=
138 ((1 - pw[patchFacei])*pMagSf[patchFacei]/magSqr(d))*sqr(d);
139 }
140 }
141 else
142 {
143 forAll(pd, patchFacei)
144 {
145 const vector& d = pd[patchFacei];
146
147 dd[faceCells[patchFacei]] +=
148 (pMagSf[patchFacei]/magSqr(d))*sqr(d);
149 }
150 }
151 }
152
153
154 // Invert the dd tensor
155 const symmTensorField invDd(inv(dd));
156
157
158 // Revisit all faces and calculate the pVectors_ and nVectors_ 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 const scalar magSfByMagSqrd = magSf[facei]/magSqr(d);
166
167 pVectors_[facei] = (1.0 - w[facei])*magSfByMagSqrd*(invDd[own] & d);
168 nVectors_[facei] = -w[facei]*magSfByMagSqrd*(invDd[nei] & d);
169 }
170
171 forAll(pVectorsBf, patchi)
172 {
173 fvsPatchVectorField& patchLsP = pVectorsBf[patchi];
174
175 const fvsPatchScalarField& pw = w.boundaryField()[patchi];
176 const fvsPatchScalarField& pMagSf = magSf.boundaryField()[patchi];
177
178 const fvPatch& p = pw.patch();
179 const labelUList& faceCells = p.faceCells();
180
181 // Build the d-vectors
182 const vectorField pd(p.delta());
183
184 if (pw.coupled())
185 {
186 forAll(pd, patchFacei)
187 {
188 const vector& d = pd[patchFacei];
189
190 patchLsP[patchFacei] =
191 ((1.0 - pw[patchFacei])*pMagSf[patchFacei]/magSqr(d))
192 *(invDd[faceCells[patchFacei]] & d);
193 }
194 }
195 else
196 {
197 forAll(pd, patchFacei)
198 {
199 const vector& d = pd[patchFacei];
200
201 patchLsP[patchFacei] =
202 pMagSf[patchFacei]*(1.0/magSqr(d))
203 *(invDd[faceCells[patchFacei]] & d);
204 }
205 }
206 }
207
208 DebugInfo << "Finished calculating least square gradient vectors" << nl;
209}
210
211
213{
214 calcLeastSquaresVectors();
215 return true;
216}
217
218
219// ************************************************************************* //
GeometricBoundaryField< vector, fvsPatchField, surfaceMesh > Boundary
Type of boundary fields.
const volVectorField & C() const
Return cell centres as volVectorField.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
Least-squares gradient scheme vectors.
virtual bool movePoints()
Delete the least square vectors when the mesh moves.
virtual ~leastSquaresVectors()
Destructor.
virtual const surfaceScalarField & weights() const
Return reference to linear difference weighting factors.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
dynamicFvMesh & mesh
#define DebugInfo
Report an information message using Foam::Info.
#define DebugInFunction
Report an information message using Foam::Info.
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:83
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
Definition: symmTensor.H:59
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Field< vector > vectorField
Specialisation of Field<T> for vector.
Field< symmTensor > symmTensorField
Specialisation of Field<T> for symmTensor.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
UList< label > labelUList
A UList of labels.
Definition: UList.H:85
fvsPatchField< vector > fvsPatchVectorField
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
fvsPatchField< scalar > fvsPatchScalarField
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
volScalarField & C
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333