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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
29 #include "faceList.H"
30 #include "demandDrivenData.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
38 
39 template<class Patch>
40 const scalarListList&
41 PrimitivePatchInterpolation<Patch>::faceToPointWeights() const
42 {
43  if (!faceToPointWeightsPtr_)
44  {
45  makeFaceToPointWeights();
46  }
47 
48  return *faceToPointWeightsPtr_;
49 }
50 
51 
52 template<class Patch>
53 void PrimitivePatchInterpolation<Patch>::makeFaceToPointWeights() const
54 {
55  if (faceToPointWeightsPtr_)
56  {
58  << "Face-to-edge weights already calculated"
59  << abort(FatalError);
60  }
61 
62  const pointField& points = patch_.localPoints();
63  const List<typename Patch::FaceType>& faces = patch_.localFaces();
64 
65  faceToPointWeightsPtr_ = new scalarListList(points.size());
66  scalarListList& weights = *faceToPointWeightsPtr_;
67 
68  // get reference to addressing
69  const labelListList& pointFaces = patch_.pointFaces();
70 
71  forAll(pointFaces, pointi)
72  {
73  const labelList& curFaces = pointFaces[pointi];
74 
75  scalarList& pw = weights[pointi];
76  pw.setSize(curFaces.size());
77 
78  scalar sumw = 0.0;
79 
80  forAll(curFaces, facei)
81  {
82  pw[facei] =
83  1.0/mag(faces[curFaces[facei]].centre(points) - points[pointi]);
84  sumw += pw[facei];
85  }
86 
87  forAll(curFaces, facei)
88  {
89  pw[facei] /= sumw;
90  }
91  }
92 }
93 
94 
95 template<class Patch>
96 const scalarList&
97 PrimitivePatchInterpolation<Patch>::faceToEdgeWeights() const
98 {
99  if (!faceToEdgeWeightsPtr_)
100  {
101  makeFaceToEdgeWeights();
102  }
103 
104  return *faceToEdgeWeightsPtr_;
105 }
106 
107 
108 template<class Patch>
109 void PrimitivePatchInterpolation<Patch>::makeFaceToEdgeWeights() const
110 {
111  if (faceToEdgeWeightsPtr_)
112  {
114  << "Face-to-edge weights already calculated"
115  << abort(FatalError);
116  }
117 
118  const pointField& points = patch_.localPoints();
119  const List<typename Patch::FaceType>& faces = patch_.localFaces();
120  const edgeList& edges = patch_.edges();
121  const labelListList& edgeFaces = patch_.edgeFaces();
122 
123  faceToEdgeWeightsPtr_ = new scalarList(patch_.nInternalEdges());
124  scalarList& weights = *faceToEdgeWeightsPtr_;
125 
126  forAll(weights, edgei)
127  {
128  vector P = faces[edgeFaces[edgei][0]].centre(points);
129  vector N = faces[edgeFaces[edgei][1]].centre(points);
130  vector S = points[edges[edgei].start()];
131  vector e = edges[edgei].vec(points);
132 
133  scalar alpha =
134  -(((N - P)^(S - P))&((N - P)^e))/(((N - P)^e )&((N - P)^e));
135 
136  vector E = S + alpha*e;
137 
138  weights[edgei] = mag(N - E)/(mag(N - E) + mag(E - P));
139  }
140 }
141 
142 
143 template<class Patch>
144 void PrimitivePatchInterpolation<Patch>::clearWeights()
145 {
146  deleteDemandDrivenData(faceToPointWeightsPtr_);
147  deleteDemandDrivenData(faceToEdgeWeightsPtr_);
148 }
149 
150 
151 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
152 
153 template<class Patch>
155 :
156  patch_(p),
157  faceToPointWeightsPtr_(nullptr),
158  faceToEdgeWeightsPtr_(nullptr)
159 {}
160 
161 
162 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
163 
164 template<class Patch>
166 {
167  clearWeights();
168 }
169 
170 
171 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
172 
173 template<class Patch>
174 template<class Type>
176 (
177  const Field<Type>& ff
178 ) const
179 {
180  // Check size of the given field
181  if (ff.size() != patch_.size())
182  {
184  << "given field does not correspond to patch. Patch size: "
185  << patch_.size() << " field size: " << ff.size()
186  << abort(FatalError);
187  }
188 
189  tmp<Field<Type>> tresult
190  (
191  new Field<Type>
192  (
193  patch_.nPoints(), Zero
194  )
195  );
196 
197  Field<Type>& result = tresult.ref();
198 
199  const labelListList& pointFaces = patch_.pointFaces();
200  const scalarListList& weights = faceToPointWeights();
201 
202  forAll(pointFaces, pointi)
203  {
204  const labelList& curFaces = pointFaces[pointi];
205  const scalarList& w = weights[pointi];
206 
207  forAll(curFaces, facei)
208  {
209  result[pointi] += w[facei]*ff[curFaces[facei]];
210  }
211  }
212 
213  return tresult;
214 }
215 
216 
217 template<class Patch>
218 template<class Type>
220 (
221  const tmp<Field<Type>>& tff
222 ) const
223 {
224  tmp<Field<Type>> tint = faceToPointInterpolate(tff());
225  tff.clear();
226  return tint;
227 }
228 
229 
230 template<class Patch>
231 template<class Type>
233 (
234  const Field<Type>& pf
235 ) const
236 {
237  if (pf.size() != patch_.nPoints())
238  {
240  << "given field does not correspond to patch. Patch size: "
241  << patch_.nPoints() << " field size: " << pf.size()
242  << abort(FatalError);
243  }
244 
245  tmp<Field<Type>> tresult
246  (
247  new Field<Type>
248  (
249  patch_.size(),
250  Zero
251  )
252  );
253 
254  Field<Type>& result = tresult.ref();
255 
256  const List<typename Patch::FaceType>& localFaces = patch_.localFaces();
257 
258  forAll(result, facei)
259  {
260  const labelList& curPoints = localFaces[facei];
261 
262  forAll(curPoints, pointi)
263  {
264  result[facei] += pf[curPoints[pointi]];
265  }
266 
267  result[facei] /= curPoints.size();
268  }
269 
270  return tresult;
271 }
272 
273 
274 template<class Patch>
275 template<class Type>
277 (
278  const tmp<Field<Type>>& tpf
279 ) const
280 {
281  tmp<Field<Type>> tint = pointToFaceInterpolate(tpf());
282  tpf.clear();
283  return tint;
284 }
285 
286 
287 template<class Patch>
288 template<class Type>
290 (
291  const Field<Type>& pf
292 ) const
293 {
294  // Check size of the given field
295  if (pf.size() != patch_.size())
296  {
298  << "given field does not correspond to patch. Patch size: "
299  << patch_.size() << " field size: " << pf.size()
300  << abort(FatalError);
301  }
302 
303  tmp<Field<Type>> tresult
304  (
305  new Field<Type>(patch_.nEdges(), Zero)
306  );
307 
308  Field<Type>& result = tresult.ref();
309 
310  const edgeList& edges = patch_.edges();
311  const labelListList& edgeFaces = patch_.edgeFaces();
312 
313  const scalarList& weights = faceToEdgeWeights();
314 
315  for (label edgei = 0; edgei < patch_.nInternalEdges(); edgei++)
316  {
317  result[edgei] =
318  weights[edgei]*pf[edgeFaces[edgei][0]]
319  + (1.0 - weights[edgei])*pf[edgeFaces[edgei][1]];
320  }
321 
322  for (label edgei = patch_.nInternalEdges(); edgei < edges.size(); edgei++)
323  {
324  result[edgei] = pf[edgeFaces[edgei][0]];
325  }
326 
327  return tresult;
328 }
329 
330 
331 template<class Patch>
332 template<class Type>
334 (
335  const tmp<Field<Type>>& tpf
336 ) const
337 {
338  tmp<Field<Type>> tint = faceToEdgeInterpolate(tpf());
339  tpf.clear();
340  return tint;
341 }
342 
343 
344 template<class Patch>
346 {
347  clearWeights();
348 
349  return true;
350 }
351 
352 
353 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
354 
355 } // End namespace Foam
356 
357 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::scalarListList
List< scalarList > scalarListList
A List of scalarList.
Definition: scalarList.H:66
Foam::scalarList
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:64
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::edgeList
List< edge > edgeList
A List of edges.
Definition: edgeList.H:63
PrimitivePatchInterpolation.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
faceList.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::Vector::centre
const Vector< Cmpt > & centre(const Foam::List< Vector< Cmpt >> &) const
Return *this (used for point which is a typedef to Vector<scalar>.
Definition: VectorI.H:119
Foam::FatalError
error FatalError
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
Foam::fv::ff
const FieldField< fvPatchField, Type > & ff(const FieldField< fvPatchField, Type > &bf)
Definition: CrankNicolsonDdtScheme.C:272
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::PrimitivePatchInterpolation
Interpolation class within a primitive patch. Allows interpolation from points to faces and vice vers...
Definition: PrimitivePatchInterpolation.H:53
Foam::List< labelList >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
N
const Vector< label > N(dict.get< Vector< label >>("N"))
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146