Loading...
Searching...
No Matches
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-2023 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
31// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32
33template<class Patch>
35Foam::PrimitivePatchInterpolation<Patch>::faceToPointWeights() const
36{
37 if (!faceToPointWeightsPtr_)
38 {
39 makeFaceToPointWeights();
40 }
41
42 return *faceToPointWeightsPtr_;
43}
44
45
46template<class Patch>
47void Foam::PrimitivePatchInterpolation<Patch>::makeFaceToPointWeights() const
48{
49 if (faceToPointWeightsPtr_)
50 {
52 << "Face-to-edge weights already calculated"
53 << abort(FatalError);
54 }
55
56 const auto& points = patch_.localPoints();
57 const auto& faces = patch_.localFaces();
58
59 faceToPointWeightsPtr_.reset(new scalarListList(points.size()));
60 auto& weights = *faceToPointWeightsPtr_;
61
62 // get reference to addressing
63 const labelListList& pointFaces = patch_.pointFaces();
64
65 forAll(pointFaces, pointi)
66 {
67 const labelList& curFaces = pointFaces[pointi];
68
69 scalarList& pw = weights[pointi];
70 pw.setSize(curFaces.size());
71
72 scalar sumw = 0.0;
73
74 forAll(curFaces, facei)
75 {
76 pw[facei] =
77 1.0/mag(faces[curFaces[facei]].centre(points) - points[pointi]);
78 sumw += pw[facei];
79 }
80
81 forAll(curFaces, facei)
82 {
83 pw[facei] /= sumw;
84 }
85 }
86}
87
88
89template<class Patch>
91Foam::PrimitivePatchInterpolation<Patch>::faceToEdgeWeights() const
92{
93 if (!faceToEdgeWeightsPtr_)
94 {
95 makeFaceToEdgeWeights();
96 }
97
98 return *faceToEdgeWeightsPtr_;
99}
100
101
102template<class Patch>
103void Foam::PrimitivePatchInterpolation<Patch>::makeFaceToEdgeWeights() const
104{
105 if (faceToEdgeWeightsPtr_)
106 {
108 << "Face-to-edge weights already calculated"
109 << abort(FatalError);
110 }
111
112 const auto& points = patch_.localPoints();
113 const auto& faces = patch_.localFaces();
114 const edgeList& edges = patch_.edges();
115 const labelListList& edgeFaces = patch_.edgeFaces();
116
117 faceToEdgeWeightsPtr_.reset(new scalarList(patch_.nInternalEdges()));
118 auto& weights = *faceToEdgeWeightsPtr_;
119
120 forAll(weights, edgei)
121 {
122 vector P = faces[edgeFaces[edgei][0]].centre(points);
123 vector N = faces[edgeFaces[edgei][1]].centre(points);
124 vector S = points[edges[edgei].start()];
125 vector e = edges[edgei].vec(points);
126
127 scalar alpha =
128 -(((N - P)^(S - P))&((N - P)^e))/(((N - P)^e )&((N - P)^e));
129
130 vector E = S + alpha*e;
131
132 weights[edgei] = mag(N - E)/(mag(N - E) + mag(E - P));
133 }
134}
135
136
137template<class Patch>
138void Foam::PrimitivePatchInterpolation<Patch>::clearWeights()
139{
140 faceToPointWeightsPtr_.reset(nullptr);
141 faceToEdgeWeightsPtr_.reset(nullptr);
142}
143
144
145// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
146
147template<class Patch>
148Foam::PrimitivePatchInterpolation<Patch>::PrimitivePatchInterpolation
149(
150 const Patch& p
151)
152:
153 patch_(p)
154{}
156
157// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
158
159template<class Patch>
160template<class Type>
163(
164 const Field<Type>& ff
165) const
166{
167 // Check size of the given field
168 if (ff.size() != patch_.size())
169 {
171 << "given field does not correspond to patch. Patch size: "
172 << patch_.size() << " field size: " << ff.size()
173 << abort(FatalError);
174 }
175
176 auto tresult = tmp<Field<Type>>::New(patch_.nPoints(), Zero);
177 auto& result = tresult.ref();
178
179 const labelListList& pointFaces = patch_.pointFaces();
180 const scalarListList& weights = faceToPointWeights();
181
182 forAll(pointFaces, pointi)
183 {
184 const labelList& curFaces = pointFaces[pointi];
185 const scalarList& w = weights[pointi];
186
187 forAll(curFaces, facei)
188 {
189 result[pointi] += w[facei]*ff[curFaces[facei]];
190 }
191 }
192
193 return tresult;
194}
195
196
197template<class Patch>
198template<class Type>
201(
202 const tmp<Field<Type>>& tff
203) const
204{
205 tmp<Field<Type>> tint = faceToPointInterpolate(tff());
206 tff.clear();
207 return tint;
208}
209
210
211template<class Patch>
212template<class Type>
215(
216 const Field<Type>& pf
217) const
218{
219 if (pf.size() != patch_.nPoints())
220 {
222 << "given field does not correspond to patch. Patch size: "
223 << patch_.nPoints() << " field size: " << pf.size()
224 << abort(FatalError);
225 }
226
227 auto tresult = tmp<Field<Type>>::New(patch_.size(), Zero);
228 auto& result = tresult.ref();
229
230 const auto& localFaces = patch_.localFaces();
231
232 forAll(result, facei)
233 {
234 const labelList& curPoints = localFaces[facei];
235
236 forAll(curPoints, pointi)
237 {
238 result[facei] += pf[curPoints[pointi]];
239 }
240
241 result[facei] /= curPoints.size();
242 }
243
244 return tresult;
245}
246
247
248template<class Patch>
249template<class Type>
252(
253 const tmp<Field<Type>>& tpf
254) const
255{
256 tmp<Field<Type>> tint = pointToFaceInterpolate(tpf());
257 tpf.clear();
258 return tint;
259}
260
261
262template<class Patch>
263template<class Type>
266(
267 const Field<Type>& pf
268) const
269{
270 // Check size of the given field
271 if (pf.size() != patch_.size())
272 {
274 << "given field does not correspond to patch. Patch size: "
275 << patch_.size() << " field size: " << pf.size()
276 << abort(FatalError);
277 }
278
279 auto tresult = tmp<Field<Type>>::New(patch_.nEdges(), Zero);
280 auto& result = tresult.ref();
281
282 const edgeList& edges = patch_.edges();
283 const labelListList& edgeFaces = patch_.edgeFaces();
284
285 const scalarList& weights = faceToEdgeWeights();
286
287 for (label edgei = 0; edgei < patch_.nInternalEdges(); edgei++)
288 {
289 result[edgei] =
290 weights[edgei]*pf[edgeFaces[edgei][0]]
291 + (1.0 - weights[edgei])*pf[edgeFaces[edgei][1]];
292 }
293
294 for (label edgei = patch_.nInternalEdges(); edgei < edges.size(); edgei++)
295 {
296 result[edgei] = pf[edgeFaces[edgei][0]];
297 }
298
299 return tresult;
300}
301
302
303template<class Patch>
304template<class Type>
307(
308 const tmp<Field<Type>>& tpf
309) const
310{
312 tpf.clear();
313 return tint;
314}
315
316
317template<class Patch>
319{
320 clearWeights();
321
322 return true;
323}
324
325
326// ************************************************************************* //
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
void setSize(label n)
Alias for resize().
Definition List.H:536
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:118
A class for managing temporary objects.
Definition tmp.H:75
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const pointField & points
List< scalarList > scalarListList
List of scalarList.
Definition scalarList.H:35
List< edge > edgeList
List of edge.
Definition edgeList.H:32
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition errorManip.H:139
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...
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
volScalarField & alpha
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
const Vector< label > N(dict.get< Vector< label > >("N"))