Loading...
Searching...
No Matches
invDistLeastSquaresVectors.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-------------------------------------------------------------------------------
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
30#include "volFields.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
37}
38
39
40// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41
43:
44 MeshObject_type(mesh),
45 pVectors_
46 (
48 (
49 "LeastSquaresP",
50 mesh_.pointsInstance(),
51 mesh_,
52 IOobject::NO_READ,
53 IOobject::NO_WRITE,
54 IOobject::NO_REGISTER
55 ),
56 mesh_,
58 ),
59 nVectors_
60 (
62 (
63 "LeastSquaresN",
64 mesh_.pointsInstance(),
65 mesh_,
66 IOobject::NO_READ,
67 IOobject::NO_WRITE,
68 IOobject::NO_REGISTER
69 ),
70 mesh_,
72 )
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
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 vector d(C[nei] - C[own]);
107 const scalar magSqrDist = d.magSqr();
108
109 if (magSqrDist > ROOTVSMALL)
110 {
111 const symmTensor wdd(sqr(d)/magSqrDist);
112 dd[own] += wdd;
113 dd[nei] += wdd;
114 }
115 }
116
117 auto& blsP = pVectors_.boundaryField();
118
119 forAll(blsP, patchi)
120 {
121 const fvsPatchVectorField& patchLsP = blsP[patchi];
122
123 const fvPatch& p = patchLsP.patch();
124 const labelUList& faceCells = p.patch().faceCells();
125
126 // Build the d-vectors
127 const vectorField pd(p.delta());
128
129 forAll(pd, patchFacei)
130 {
131 const vector& d = pd[patchFacei];
132 const scalar magSqrDist = d.magSqr();
133
134 if (magSqrDist > ROOTVSMALL)
135 {
136 dd[faceCells[patchFacei]] += sqr(d)/magSqrDist;
137 }
138 }
139 }
140
141
142 // Invert the dd tensors - including failsafe checks
143 const symmTensorField invDd(inv(dd));
144
145
146 // Revisit all faces and calculate the pVectors_ and nVectors_ vectors
147 forAll(owner, facei)
148 {
149 const label own = owner[facei];
150 const label nei = neighbour[facei];
151
152 const vector d(C[nei] - C[own]);
153 const scalar magSqrDist = d.magSqr();
154
155 if (magSqrDist > ROOTVSMALL)
156 {
157 pVectors_[facei] = (invDd[own] & d)/magSqrDist;
158 nVectors_[facei] = -(invDd[nei] & d)/magSqrDist;
159 }
160 else
161 {
162 pVectors_[facei] = Zero;
163 nVectors_[facei] = Zero;
164 }
165 }
166
167 forAll(blsP, patchi)
168 {
169 fvsPatchVectorField& patchLsP = blsP[patchi];
170
171 const fvPatch& p = patchLsP.patch();
172 const labelUList& faceCells = p.faceCells();
173
174 // Build the d-vectors
175 const vectorField pd(p.delta());
176
177 forAll(pd, patchFacei)
178 {
179 const vector& d = pd[patchFacei];
180 const scalar magSqrDist = d.magSqr();
181
182 if (magSqrDist > ROOTVSMALL)
183 {
184 patchLsP[patchFacei] =
185 (invDd[faceCells[patchFacei]] & d)/magSqrDist;
186 }
187 else
188 {
189 patchLsP[patchFacei] = Zero;
190 }
192 }
193
194 DebugInfo << "Finished calculating least square gradient vectors" << endl;
195}
196
197
199{
200 calcLeastSquaresVectors();
201 return true;
202}
203
204
205// ************************************************************************* //
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5)
Graphite solid properties.
Definition C.H:49
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
Smooth ATC in cells next to a set of patches supplied by type.
Definition faceCells.H:55
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
Least-squares gradient scheme vectors.
virtual bool movePoints()
Delete the least square vectors when the mesh moves.
virtual ~leastSquaresVectors()
Destructor.
leastSquaresVectors(const fvMesh &)
Construct given an fvMesh.
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:142
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
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
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:127
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
Vector< scalar > vector
Definition vector.H:57
UList< label > labelUList
A UList of labels.
Definition UList.H:75
fvsPatchField< vector > fvsPatchVectorField
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
Definition symmTensor.H:55
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299