leastSquaresFaVectors.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) 2016-2017 Wikki Ltd
9 Copyright (C) 2020 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 "edgeFields.H"
31#include "areaFields.H"
32#include "mapPolyMesh.H"
33#include "demandDrivenData.H"
34
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
40}
41
42
43// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44
46:
48 pVectorsPtr_(nullptr),
49 nVectorsPtr_(nullptr)
50{}
51
52
53// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
54
56{
57 deleteDemandDrivenData(pVectorsPtr_);
58 deleteDemandDrivenData(nVectorsPtr_);
59}
60
61
62// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
63
64void Foam::leastSquaresFaVectors::makeLeastSquaresVectors() const
65{
67 << "Constructing finite area least square gradient vectors" << nl;
68
69 pVectorsPtr_ = new edgeVectorField
70 (
72 (
73 "LeastSquaresP",
74 mesh().pointsInstance(),
75 mesh().thisDb(),
78 false
79 ),
80 mesh(),
82 );
83 edgeVectorField& lsP = *pVectorsPtr_;
84
85 nVectorsPtr_ = new edgeVectorField
86 (
88 (
89 "LeastSquaresN",
90 mesh().pointsInstance(),
91 mesh().thisDb(),
94 false
95 ),
96 mesh(),
98 );
99 edgeVectorField& lsN = *nVectorsPtr_;
100
101 // Set local references to mesh data
102 const labelUList& owner = mesh().owner();
103 const labelUList& neighbour = mesh().neighbour();
104
105 const areaVectorField& C = mesh().areaCentres();
106 const edgeScalarField& w = mesh().weights();
107
108
109 // Set up temporary storage for the dd tensor (before inversion)
111
112 forAll(owner, facei)
113 {
114 label own = owner[facei];
115 label nei = neighbour[facei];
116
117 vector d = C[nei] - C[own];
118 symmTensor wdd = (1.0/magSqr(d))*sqr(d);
119
120 dd[own] += wdd;
121 dd[nei] += wdd;
122 }
123
124
125 forAll(lsP.boundaryField(), patchi)
126 {
127 const faePatchScalarField& pw = w.boundaryField()[patchi];
128
129 const faPatch& p = pw.patch();
130 const labelUList& edgeFaces = p.edgeFaces();
131
132 // Build the d-vectors
133 // HJ, reconsider deltas at the boundary, consistent with FVM
134 // Current implementation is good for fixedValue boundaries, but may
135 // cause problems with fixedGradient. HJ, 4/Oct/2010
136 const vectorField pd(p.delta());
137
138 if (p.coupled())
139 {
140 forAll(pd, patchFacei)
141 {
142 const vector& d = pd[patchFacei];
143
144 dd[edgeFaces[patchFacei]] +=
145 (1.0/magSqr(d))*sqr(d);
146 }
147 }
148 else
149 {
150 forAll(pd, patchFacei)
151 {
152 const vector& d = pd[patchFacei];
153
154 dd[edgeFaces[patchFacei]] +=
155 (1.0/magSqr(d))*sqr(d);
156 }
157 }
158 }
159
160
161 // Invert the dd tensor
162 const symmTensorField invDd(inv(dd));
163
164
165 // Revisit all faces and calculate the lsP and lsN vectors
166 forAll(owner, facei)
167 {
168 label own = owner[facei];
169 label nei = neighbour[facei];
170
171 vector d = C[nei] - C[own];
172 scalar magSfByMagSqrd = 1.0/magSqr(d);
173
174 lsP[facei] = magSfByMagSqrd*(invDd[own] & d);
175 lsN[facei] = -magSfByMagSqrd*(invDd[nei] & d);
176 }
177
178 forAll(lsP.boundaryField(), patchi)
179 {
180 faePatchVectorField& patchLsP = lsP.boundaryFieldRef()[patchi];
181
182 const faePatchScalarField& pw = w.boundaryField()[patchi];
183
184 const faPatch& p = pw.patch();
185 const labelUList& edgeFaces = p.edgeFaces();
186
187 // Build the d-vectors
188 const vectorField pd(p.delta());
189
190 if (p.coupled())
191 {
192 forAll(pd, patchFacei)
193 {
194 const vector& d = pd[patchFacei];
195
196 patchLsP[patchFacei] =
197 (1.0/magSqr(d))
198 *(invDd[edgeFaces[patchFacei]] & d);
199 }
200 }
201 else
202 {
203 forAll(pd, patchFacei)
204 {
205 const vector& d = pd[patchFacei];
206
207 patchLsP[patchFacei] =
208 (1.0/magSqr(d))
209 *(invDd[edgeFaces[patchFacei]] & d);
210 }
211 }
212 }
213
215 << "Done constructing finite area least square gradient vectors" << nl;
216}
217
218
220{
221 if (!pVectorsPtr_)
222 {
223 makeLeastSquaresVectors();
224 }
225
226 return *pVectorsPtr_;
227}
228
229
231{
232 if (!nVectorsPtr_)
233 {
234 makeLeastSquaresVectors();
235 }
236
237 return *nVectorsPtr_;
238}
239
240
242{
244 << "Clearing least square data" << nl;
245
246 deleteDemandDrivenData(pVectorsPtr_);
247 deleteDemandDrivenData(nVectorsPtr_);
248
249 return true;
250}
251
252
253// ************************************************************************* //
Graphite solid properties.
Definition: C.H:53
Generic GeometricField class.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:91
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition: faMesh.H:100
const labelUList & owner() const
Internal face owner. Note bypassing virtual mechanism so.
Definition: fvMesh.H:417
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:423
Least-squares gradient scheme vectors for the Finite Area method.
virtual bool movePoints()
Delete the least square vectors when the mesh moves.
virtual ~leastSquaresFaVectors()
Destructor.
const edgeVectorField & pVectors() const
Return reference to owner least square vectors.
const edgeVectorField & nVectors() const
Return reference to neighbour least square vectors.
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
Template functions to aid in the implementation of demand driven data.
const labelList nFaces(UPstream::listGatherValues< label >(aMesh.nFaces()))
#define DebugInFunction
Report an information message using Foam::Info.
Namespace for OpenFOAM.
GeometricField< vector, faePatchField, edgeMesh > edgeVectorField
Definition: edgeFieldsFwd.H:64
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
faePatchField< scalar > faePatchScalarField
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)
faePatchField< vector > faePatchVectorField
UList< label > labelUList
A UList of labels.
Definition: UList.H:85
void deleteDemandDrivenData(DataPtr &dataPtr)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
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