pointVolInterpolation.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) Wikki Ltd
9  Copyright (C) 2019 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 "pointVolInterpolation.H"
30 #include "fvMesh.H"
31 #include "pointFields.H"
32 #include "volFields.H"
33 #include "emptyFvPatch.H"
34 #include "SubField.H"
35 #include "demandDrivenData.H"
36 #include "globalMeshData.H"
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42  defineTypeNameAndDebug(pointVolInterpolation, 0);
43 }
44 
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
48 void Foam::pointVolInterpolation::makeWeights() const
49 {
50  if (volWeightsPtr_)
51  {
53  << "weighting factors already calculated"
54  << abort(FatalError);
55  }
56 
57  if (debug)
58  {
59  Info<< "pointVolInterpolation::makeWeights() : "
60  << "constructing weighting factors"
61  << endl;
62  }
63 
64  const pointField& points = vMesh().points();
65  const labelListList& cellPoints = vMesh().cellPoints();
66  const vectorField& cellCentres = vMesh().cellCentres();
67 
68  // Allocate storage for weighting factors
69  volWeightsPtr_ = new FieldField<Field, scalar>(cellCentres.size());
70  FieldField<Field, scalar>& weightingFactors = *volWeightsPtr_;
71 
72  forAll(weightingFactors, pointi)
73  {
74  weightingFactors.set
75  (
76  pointi,
77  new scalarField(cellPoints[pointi].size())
78  );
79  }
80 
81 
82  // Calculate inverse distances between points and cell centres
83  // and store in weighting factor array
84  forAll(cellCentres, cellI)
85  {
86  const labelList& curCellPoints = cellPoints[cellI];
87 
88  forAll(curCellPoints, cellPointI)
89  {
90  weightingFactors[cellI][cellPointI] = 1.0/
91  mag
92  (
93  cellCentres[cellI] - points[curCellPoints[cellPointI]]
94  );
95  }
96  }
97 
98  scalarField pointVolSumWeights(cellCentres.size(), Zero);
99 
100  // Calculate weighting factors using inverse distance weighting
101  forAll(cellCentres, cellI)
102  {
103  const labelList& curCellPoints = cellPoints[cellI];
104 
105  forAll(curCellPoints, cellPointI)
106  {
107  pointVolSumWeights[cellI] += weightingFactors[cellI][cellPointI];
108  }
109  }
110 
111  forAll(cellCentres, cellI)
112  {
113  const labelList& curCellPoints = cellPoints[cellI];
114 
115  forAll(curCellPoints, cellPointI)
116  {
117  weightingFactors[cellI][cellPointI] /= pointVolSumWeights[cellI];
118  }
119  }
120 
121  if (debug)
122  {
123  Info<< "pointVolInterpolation::makeWeights() : "
124  << "finished constructing weighting factors"
125  << endl;
126  }
127 }
128 
129 
130 // Do what is necessary if the mesh has changed topology
131 void Foam::pointVolInterpolation::clearAddressing() const
132 {
133  deleteDemandDrivenData(patchInterpolatorsPtr_);
134 }
135 
136 
137 // Do what is necessary if the mesh has moved
138 void Foam::pointVolInterpolation::clearGeom() const
139 {
140  deleteDemandDrivenData(volWeightsPtr_);
141 }
142 
143 
145 Foam::pointVolInterpolation::patchInterpolators() const
146 {
147  if (!patchInterpolatorsPtr_)
148  {
149  const fvBoundaryMesh& bdry = vMesh().boundary();
150 
151  patchInterpolatorsPtr_ =
152  new PtrList<primitivePatchInterpolation>(bdry.size());
153 
154  forAll(bdry, patchI)
155  {
156  patchInterpolatorsPtr_->set
157  (
158  patchI,
159  new primitivePatchInterpolation(bdry[patchI].patch())
160  );
161  }
162  }
163 
164  return *patchInterpolatorsPtr_;
165 }
166 
167 
168 // * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
169 
171 (
172  const pointMesh& pm,
173  const fvMesh& vm
174 )
175 :
176  pointMesh_(pm),
177  fvMesh_(vm),
178  volWeightsPtr_(nullptr),
179  patchInterpolatorsPtr_(nullptr)
180 {}
181 
182 
183 // * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * //
184 
186 {
187  clearAddressing();
188  clearGeom();
189 }
190 
191 
192 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
193 
194 // Return point weights
197 {
198  // If weighting factor array not present then construct
199  if (!volWeightsPtr_)
200  {
201  makeWeights();
202  }
203 
204  return *volWeightsPtr_;
205 }
206 
207 
208 // Do what is necessary if the mesh has moved
210 {
211  clearAddressing();
212  clearGeom();
213 }
214 
215 
216 // Do what is necessary if the mesh has moved
218 {
219  clearGeom();
220 
221  return true;
222 }
223 
224 
225 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
volFields.H
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1069
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
SubField.H
Foam::FieldField
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:53
Foam::pointVolInterpolation::updateTopology
void updateTopology()
Update mesh topology using the morph engine.
Definition: pointVolInterpolation.C:209
Foam::primitiveMesh::cellPoints
const labelListList & cellPoints() const
Definition: primitiveMeshCellPoints.C:34
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
globalMeshData.H
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
pointVolInterpolation.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::pointVolInterpolation::vMesh
const fvMesh & vMesh() const
Definition: pointVolInterpolation.H:98
Foam::primitivePatchInterpolation
PrimitivePatchInterpolation< primitivePatch > primitivePatchInterpolation
Foam::primitivePatchInterpolation.
Definition: primitivePatchInterpolation.H:46
Foam::pointVolInterpolation::movePoints
bool movePoints()
Correct weighting factors for moving mesh.
Definition: pointVolInterpolation.C:217
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::pointVolInterpolation::volWeights
const FieldField< Field, scalar > & volWeights() const
Return reference to weights arrays.
Definition: pointVolInterpolation.C:196
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::pointVolInterpolation::~pointVolInterpolation
~pointVolInterpolation()
Definition: pointVolInterpolation.C:185
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
Foam::FatalError
error FatalError
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::pointMesh
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:51
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
emptyFvPatch.H
Foam::pointVolInterpolation::pointVolInterpolation
pointVolInterpolation(const pointMesh &, const fvMesh &)
Constructor given pointMesh and fvMesh.
Definition: pointVolInterpolation.C:171
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:84
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
pointFields.H