PrimitivePatchInterpolation.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) 2011-2016 OpenFOAM Foundation
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 "faceList.H"
31#include "demandDrivenData.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
36{
37
38// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
39
40template<class Patch>
41const scalarListList&
42PrimitivePatchInterpolation<Patch>::faceToPointWeights() const
43{
44 if (!faceToPointWeightsPtr_)
45 {
46 makeFaceToPointWeights();
47 }
48
49 return *faceToPointWeightsPtr_;
50}
51
52
53template<class Patch>
54void PrimitivePatchInterpolation<Patch>::makeFaceToPointWeights() const
55{
56 if (faceToPointWeightsPtr_)
57 {
59 << "Face-to-edge weights already calculated"
60 << abort(FatalError);
61 }
62
63 const pointField& points = patch_.localPoints();
64 const List<typename Patch::face_type>& faces = patch_.localFaces();
65
66 faceToPointWeightsPtr_ = new scalarListList(points.size());
67 auto& weights = *faceToPointWeightsPtr_;
68
69 // get reference to addressing
70 const labelListList& pointFaces = patch_.pointFaces();
71
72 forAll(pointFaces, pointi)
73 {
74 const labelList& curFaces = pointFaces[pointi];
75
76 scalarList& pw = weights[pointi];
77 pw.setSize(curFaces.size());
78
79 scalar sumw = 0.0;
80
81 forAll(curFaces, facei)
82 {
83 pw[facei] =
84 1.0/mag(faces[curFaces[facei]].centre(points) - points[pointi]);
85 sumw += pw[facei];
86 }
87
88 forAll(curFaces, facei)
89 {
90 pw[facei] /= sumw;
91 }
92 }
93}
94
95
96template<class Patch>
97const scalarList&
98PrimitivePatchInterpolation<Patch>::faceToEdgeWeights() const
99{
100 if (!faceToEdgeWeightsPtr_)
101 {
102 makeFaceToEdgeWeights();
104
105 return *faceToEdgeWeightsPtr_;
106}
108
109template<class Patch>
111{
112 if (faceToEdgeWeightsPtr_)
113 {
115 << "Face-to-edge weights already calculated"
116 << abort(FatalError);
117 }
118
119 const pointField& points = patch_.localPoints();
120 const List<typename Patch::face_type>& faces = patch_.localFaces();
121 const edgeList& edges = patch_.edges();
122 const labelListList& edgeFaces = patch_.edgeFaces();
123
124 faceToEdgeWeightsPtr_ = new scalarList(patch_.nInternalEdges());
125 auto& weights = *faceToEdgeWeightsPtr_;
126
127 forAll(weights, edgei)
128 {
129 vector P = faces[edgeFaces[edgei][0]].centre(points);
130 vector N = faces[edgeFaces[edgei][1]].centre(points);
131 vector S = points[edges[edgei].start()];
132 vector e = edges[edgei].vec(points);
134 scalar alpha =
135 -(((N - P)^(S - P))&((N - P)^e))/(((N - P)^e )&((N - P)^e));
136
137 vector E = S + alpha*e;
138
139 weights[edgei] = mag(N - E)/(mag(N - E) + mag(E - P));
141}
142
143
144template<class Patch>
145void PrimitivePatchInterpolation<Patch>::clearWeights()
147 deleteDemandDrivenData(faceToPointWeightsPtr_);
148 deleteDemandDrivenData(faceToEdgeWeightsPtr_);
149}
150
151
152// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
153
154template<class Patch>
156:
157 patch_(p),
158 faceToPointWeightsPtr_(nullptr),
159 faceToEdgeWeightsPtr_(nullptr)
160{}
161
162
163// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
164
165template<class Patch>
167{
168 clearWeights();
169}
170
171
172// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
173
174template<class Patch>
175template<class Type>
177(
178 const Field<Type>& ff
179) const
180{
181 // Check size of the given field
182 if (ff.size() != patch_.size())
183 {
185 << "given field does not correspond to patch. Patch size: "
186 << patch_.size() << " field size: " << ff.size()
187 << abort(FatalError);
188 }
189
190 tmp<Field<Type>> tresult
191 (
192 new Field<Type>
193 (
194 patch_.nPoints(), Zero
195 )
196 );
197
198 Field<Type>& result = tresult.ref();
199
200 const labelListList& pointFaces = patch_.pointFaces();
201 const scalarListList& weights = faceToPointWeights();
202
203 forAll(pointFaces, pointi)
204 {
205 const labelList& curFaces = pointFaces[pointi];
206 const scalarList& w = weights[pointi];
207
208 forAll(curFaces, facei)
209 {
210 result[pointi] += w[facei]*ff[curFaces[facei]];
211 }
212 }
213
214 return tresult;
215}
216
217
218template<class Patch>
219template<class Type>
221(
222 const tmp<Field<Type>>& tff
223) const
224{
225 tmp<Field<Type>> tint = faceToPointInterpolate(tff());
226 tff.clear();
227 return tint;
228}
229
230
231template<class Patch>
232template<class Type>
234(
235 const Field<Type>& pf
236) const
237{
238 if (pf.size() != patch_.nPoints())
239 {
241 << "given field does not correspond to patch. Patch size: "
242 << patch_.nPoints() << " field size: " << pf.size()
243 << abort(FatalError);
244 }
245
246 tmp<Field<Type>> tresult
247 (
248 new Field<Type>
249 (
250 patch_.size(),
251 Zero
252 )
253 );
254
255 Field<Type>& result = tresult.ref();
256
257 const List<typename Patch::face_type>& localFaces = patch_.localFaces();
258
259 forAll(result, facei)
260 {
261 const labelList& curPoints = localFaces[facei];
262
263 forAll(curPoints, pointi)
264 {
265 result[facei] += pf[curPoints[pointi]];
266 }
267
268 result[facei] /= curPoints.size();
269 }
270
271 return tresult;
272}
273
274
275template<class Patch>
276template<class Type>
278(
279 const tmp<Field<Type>>& tpf
280) const
281{
282 tmp<Field<Type>> tint = pointToFaceInterpolate(tpf());
283 tpf.clear();
284 return tint;
285}
286
287
288template<class Patch>
289template<class Type>
291(
292 const Field<Type>& pf
293) const
294{
295 // Check size of the given field
296 if (pf.size() != patch_.size())
297 {
299 << "given field does not correspond to patch. Patch size: "
300 << patch_.size() << " field size: " << pf.size()
301 << abort(FatalError);
302 }
303
304 auto tresult = tmp<Field<Type>>::New(patch_.nEdges(), Zero);
305 auto& result = tresult.ref();
306
307 const edgeList& edges = patch_.edges();
308 const labelListList& edgeFaces = patch_.edgeFaces();
309
310 const scalarList& weights = faceToEdgeWeights();
311
312 for (label edgei = 0; edgei < patch_.nInternalEdges(); edgei++)
313 {
314 result[edgei] =
315 weights[edgei]*pf[edgeFaces[edgei][0]]
316 + (1.0 - weights[edgei])*pf[edgeFaces[edgei][1]];
317 }
318
319 for (label edgei = patch_.nInternalEdges(); edgei < edges.size(); edgei++)
320 {
321 result[edgei] = pf[edgeFaces[edgei][0]];
322 }
323
324 return tresult;
325}
326
327
328template<class Patch>
329template<class Type>
331(
332 const tmp<Field<Type>>& tpf
333) const
334{
335 tmp<Field<Type>> tint = faceToEdgeInterpolate(tpf());
336 tpf.clear();
337 return tint;
338}
339
340
341template<class Patch>
343{
344 clearWeights();
345
346 return true;
347}
348
349
350// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
351
352} // End namespace Foam
353
354// ************************************************************************* //
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
void setSize(const label n)
Alias for resize()
Definition: List.H:218
Interpolation class within a primitive patch. Allows interpolation from points to faces and vice vers...
bool movePoints()
Do what is necessary if the mesh has moved.
tmp< Field< Type > > faceToPointInterpolate(const Field< Type > &ff) const
Interpolate from faces to points.
tmp< Field< Type > > faceToEdgeInterpolate(const Field< Type > &ff) const
Interpolate from faces to edges.
tmp< Field< Type > > pointToFaceInterpolate(const Field< Type > &pf) const
Interpolate from points to faces.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
volScalarField & p
Template functions to aid in the implementation of demand driven data.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const pointField & points
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition: List.H:66
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:64
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
errorManip< error > abort(error &err)
Definition: errorManip.H:144
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
List< scalarList > scalarListList
A List of scalarList.
Definition: scalarList.H:66
void deleteDemandDrivenData(DataPtr &dataPtr)
volScalarField & alpha
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
const Vector< label > N(dict.get< Vector< label > >("N"))