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