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 -------------------------------------------------------------------------------
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 
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 template<class Patch>
36 {
37  if (!faceToPointWeightsPtr_)
38  {
39  makeFaceToPointWeights();
40  }
41 
42  return *faceToPointWeightsPtr_;
43 }
44 
45 
46 template<class Patch>
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 
89 template<class Patch>
90 const Foam::scalarList&
92 {
93  if (!faceToEdgeWeightsPtr_)
94  {
95  makeFaceToEdgeWeights();
96  }
97 
98  return *faceToEdgeWeightsPtr_;
99 }
100 
101 
102 template<class Patch>
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 
137 template<class Patch>
139 {
140  faceToPointWeightsPtr_.reset(nullptr);
141  faceToEdgeWeightsPtr_.reset(nullptr);
142 }
143 
144 
145 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
146 
147 template<class Patch>
149 (
150  const Patch& p
151 )
152 :
153  patch_(p),
154  faceToPointWeightsPtr_(nullptr),
155  faceToEdgeWeightsPtr_(nullptr)
156 {}
157 
159 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
160 
161 template<class Patch>
162 template<class Type>
165 (
166  const Field<Type>& ff
167 ) const
168 {
169  // Check size of the given field
170  if (ff.size() != patch_.size())
171  {
173  << "given field does not correspond to patch. Patch size: "
174  << patch_.size() << " field size: " << ff.size()
175  << abort(FatalError);
176  }
177 
178  auto tresult = tmp<Field<Type>>::New(patch_.nPoints(), Zero);
179  auto& result = tresult.ref();
180 
181  const labelListList& pointFaces = patch_.pointFaces();
182  const scalarListList& weights = faceToPointWeights();
183 
184  forAll(pointFaces, pointi)
185  {
186  const labelList& curFaces = pointFaces[pointi];
187  const scalarList& w = weights[pointi];
188 
189  forAll(curFaces, facei)
190  {
191  result[pointi] += w[facei]*ff[curFaces[facei]];
192  }
193  }
194 
195  return tresult;
196 }
197 
198 
199 template<class Patch>
200 template<class Type>
203 (
204  const tmp<Field<Type>>& tff
205 ) const
206 {
207  tmp<Field<Type>> tint = faceToPointInterpolate(tff());
208  tff.clear();
209  return tint;
210 }
211 
212 
213 template<class Patch>
214 template<class Type>
217 (
218  const Field<Type>& pf
219 ) const
220 {
221  if (pf.size() != patch_.nPoints())
222  {
224  << "given field does not correspond to patch. Patch size: "
225  << patch_.nPoints() << " field size: " << pf.size()
226  << abort(FatalError);
227  }
228 
229  auto tresult = tmp<Field<Type>>::New(patch_.size(), Zero);
230  auto& result = tresult.ref();
231 
232  const auto& localFaces = patch_.localFaces();
233 
234  forAll(result, facei)
235  {
236  const labelList& curPoints = localFaces[facei];
237 
238  forAll(curPoints, pointi)
239  {
240  result[facei] += pf[curPoints[pointi]];
241  }
242 
243  result[facei] /= curPoints.size();
244  }
245 
246  return tresult;
247 }
248 
249 
250 template<class Patch>
251 template<class Type>
254 (
255  const tmp<Field<Type>>& tpf
256 ) const
257 {
258  tmp<Field<Type>> tint = pointToFaceInterpolate(tpf());
259  tpf.clear();
260  return tint;
261 }
262 
263 
264 template<class Patch>
265 template<class Type>
268 (
269  const Field<Type>& pf
270 ) const
271 {
272  // Check size of the given field
273  if (pf.size() != patch_.size())
274  {
276  << "given field does not correspond to patch. Patch size: "
277  << patch_.size() << " field size: " << pf.size()
278  << abort(FatalError);
279  }
280 
281  auto tresult = tmp<Field<Type>>::New(patch_.nEdges(), Zero);
282  auto& result = tresult.ref();
283 
284  const edgeList& edges = patch_.edges();
285  const labelListList& edgeFaces = patch_.edgeFaces();
286 
287  const scalarList& weights = faceToEdgeWeights();
288 
289  for (label edgei = 0; edgei < patch_.nInternalEdges(); edgei++)
290  {
291  result[edgei] =
292  weights[edgei]*pf[edgeFaces[edgei][0]]
293  + (1.0 - weights[edgei])*pf[edgeFaces[edgei][1]];
294  }
295 
296  for (label edgei = patch_.nInternalEdges(); edgei < edges.size(); edgei++)
297  {
298  result[edgei] = pf[edgeFaces[edgei][0]];
299  }
300 
301  return tresult;
302 }
303 
304 
305 template<class Patch>
306 template<class Type>
309 (
310  const tmp<Field<Type>>& tpf
311 ) const
312 {
313  tmp<Field<Type>> tint = faceToEdgeInterpolate(tpf());
314  tpf.clear();
315  return tint;
316 }
317 
318 
319 template<class Patch>
321 {
322  clearWeights();
323 
324  return true;
325 }
326 
327 
328 // ************************************************************************* //
List< scalar > scalarList
List of scalar.
Definition: scalarList.H:32
Interpolation class within a primitive patch. Allows interpolation from points to faces and vice vers...
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
tmp< Field< Type > > pointToFaceInterpolate(const Field< Type > &pf) const
Interpolate from points to faces.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
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.
tmp< Field< Type > > faceToEdgeInterpolate(const Field< Type > &ff) const
Interpolate from faces to edges.
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
void setSize(const label n)
Alias for resize()
Definition: List.H:316
const pointField & points
Generic templated field type.
Definition: Field.H:62
List< scalarList > scalarListList
List of scalarList.
Definition: scalarList.H:35
bool movePoints()
Do what is necessary if the mesh has moved.
Vector< scalar > vector
Definition: vector.H:57
errorManip< error > abort(error &err)
Definition: errorManip.H:139
const Vector< Cmpt > & centre(const Foam::UList< Vector< Cmpt >> &) const noexcept
Return this (for point which is a typedef to Vector<scalar>)
Definition: VectorI.H:67
const Vector< label > N(dict.get< Vector< label >>("N"))
tmp< Field< Type > > faceToPointInterpolate(const Field< Type > &ff) const
Interpolate from faces to points.
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
const FieldField< fvPatchField, Type > & ff(const FieldField< fvPatchField, Type > &bf)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127