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 -------------------------------------------------------------------------------
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 "leastSquaresFaVectors.H"
30 #include "edgeFields.H"
31 #include "areaFields.H"
32 #include "mapPolyMesh.H"
33 #include "demandDrivenData.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(leastSquaresFaVectors, 0);
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 
64 void Foam::leastSquaresFaVectors::makeLeastSquaresVectors() const
65 {
67  << "Constructing finite area least square gradient vectors" << nl;
68 
69  pVectorsPtr_ = new edgeVectorField
70  (
71  IOobject
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  (
87  IOobject
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)
110  symmTensorField dd(mesh().nFaces(), Zero);
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 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::SymmTensor< scalar >
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::edgeVectorField
GeometricField< vector, faePatchField, edgeMesh > edgeVectorField
Definition: edgeFieldsFwd.H:57
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
mapPolyMesh.H
Foam::MoveableMeshObject
Definition: MeshObject.H:220
Foam::leastSquaresFaVectors::movePoints
virtual bool movePoints()
Delete the least square vectors when the mesh moves.
Definition: leastSquaresFaVectors.C:241
Foam::faePatchVectorField
faePatchField< vector > faePatchVectorField
Definition: faePatchFieldsFwd.H:47
C
volScalarField & C
Definition: readThermalProperties.H:102
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::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
Foam::Field< symmTensor >
Foam::faePatchScalarField
faePatchField< scalar > faePatchScalarField
Definition: faePatchFieldsFwd.H:44
leastSquaresFaVectors.H
Foam::inv
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:73
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:365
Foam::leastSquaresFaVectors::~leastSquaresFaVectors
virtual ~leastSquaresFaVectors()
Destructor.
Definition: leastSquaresFaVectors.C:55
edgeFields.H
Foam::leastSquaresFaVectors
Least-squares gradient scheme vectors for the Finite Area method.
Definition: leastSquaresFaVectors.H:56
Foam::leastSquaresFaVectors::leastSquaresFaVectors
leastSquaresFaVectors(const faMesh &)
Construct given an faMesh.
Definition: leastSquaresFaVectors.C:45
areaFields.H
Foam::fvMesh::neighbour
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:411
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::fvMesh::owner
const labelUList & owner() const
Internal face owner. Note bypassing virtual mechanism so.
Definition: fvMesh.H:405
Foam::leastSquaresFaVectors::pVectors
const edgeVectorField & pVectors() const
Return reference to owner least square vectors.
Definition: leastSquaresFaVectors.C:219
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:783
Foam::Vector< scalar >
Foam::surfaceInterpolation::weights
virtual const surfaceScalarField & weights() const
Return reference to linear difference weighting factors.
Definition: surfaceInterpolation.C:102
Foam::UList< label >
Foam::faMesh
Finite area mesh. Used for 2-D non-Euclidian finite area method.
Definition: faMesh.H:77
Foam::leastSquaresFaVectors::nVectors
const edgeVectorField & nVectors() const
Return reference to neighbour least square vectors.
Definition: leastSquaresFaVectors.C:230
Foam::MeshObject
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:88
Foam::labelUList
UList< label > labelUList
A UList of labels.
Definition: UList.H:80
Foam::C
Graphite solid properties.
Definition: C.H:50
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:53
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62