1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright held by original author
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
11 OpenFOAM is free software; you can redistribute it and/or modify it
12 under the terms of the GNU General Public License as published by the
13 Free Software Foundation; either version 2 of the License, or (at your
14 option) any later version.
16 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM; if not, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "PrimitivePatchInterpolation.H"
29 #include "demandDrivenData.H"
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
40 PrimitivePatchInterpolation<Patch>::faceToPointWeights() const
42 if (!faceToPointWeightsPtr_)
44 makeFaceToPointWeights();
47 return *faceToPointWeightsPtr_;
52 void PrimitivePatchInterpolation<Patch>::makeFaceToPointWeights() const
54 if (faceToPointWeightsPtr_)
58 "PrimitivePatchInterpolation<Patch>::makeFaceToPointWeights() const"
59 ) << "Face-to-edge weights already calculated"
63 const pointField& points = patch_.localPoints();
64 const faceList& faces = patch_.localFaces();
66 faceToPointWeightsPtr_ = new scalarListList(points.size());
67 scalarListList& weights = *faceToPointWeightsPtr_;
69 // get reference to addressing
70 const labelListList& pointFaces = patch_.pointFaces();
72 forAll(pointFaces, pointi)
74 const labelList& curFaces = pointFaces[pointi];
76 scalarList& pw = weights[pointi];
77 pw.setSize(curFaces.size());
81 forAll(curFaces, facei)
84 1.0/mag(faces[curFaces[facei]].centre(points) - points[pointi]);
88 forAll(curFaces, facei)
98 PrimitivePatchInterpolation<Patch>::faceToEdgeWeights() const
100 if (!faceToEdgeWeightsPtr_)
102 makeFaceToEdgeWeights();
105 return *faceToEdgeWeightsPtr_;
109 template<class Patch>
110 void PrimitivePatchInterpolation<Patch>::makeFaceToEdgeWeights() const
112 if (faceToEdgeWeightsPtr_)
116 "PrimitivePatchInterpolation<Patch>::makeFaceToEdgeWeights() const"
117 ) << "Face-to-edge weights already calculated"
118 << abort(FatalError);
121 const pointField& points = patch_.localPoints();
122 const faceList& faces = patch_.localFaces();
123 const edgeList& edges = patch_.edges();
124 const labelListList& edgeFaces = patch_.edgeFaces();
126 faceToEdgeWeightsPtr_ = new scalarList(patch_.nInternalEdges());
127 scalarList& weights = *faceToEdgeWeightsPtr_;
129 for (label edgei = 0; edgei < weights.size(); edgei++)
131 vector P = faces[edgeFaces[edgei][0]].centre(points);
132 vector N = faces[edgeFaces[edgei][1]].centre(points);
133 vector S = points[edges[edgei].start()];
134 vector e = edges[edgei].vec(points);
137 -(((N - P)^(S - P))&((N - P)^e))/(((N - P)^e )&((N - P)^e));
139 vector E = S + alpha*e;
141 weights[edgei] = mag(N - E)/(mag(N - E) + mag(E - P));
146 template<class Patch>
147 void PrimitivePatchInterpolation<Patch>::clearWeights()
149 deleteDemandDrivenData(faceToPointWeightsPtr_);
150 deleteDemandDrivenData(faceToEdgeWeightsPtr_);
154 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
156 template<class Patch>
157 PrimitivePatchInterpolation<Patch>::PrimitivePatchInterpolation(const Patch& p)
160 faceToPointWeightsPtr_(NULL),
161 faceToEdgeWeightsPtr_(NULL)
165 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
167 template<class Patch>
168 PrimitivePatchInterpolation<Patch>::~PrimitivePatchInterpolation()
174 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
176 template<class Patch>
178 tmp<Field<Type> > PrimitivePatchInterpolation<Patch>::faceToPointInterpolate
180 const Field<Type>& ff
183 // Check size of the given field
184 if (ff.size() != patch_.size())
188 "tmp<Field<Type> > PrimitivePatchInterpolation::"
189 "faceToPointInterpolate(const Field<Type> ff)"
190 ) << "given field does not correspond to patch. Patch size: "
191 << patch_.size() << " field size: " << ff.size()
192 << abort(FatalError);
195 tmp<Field<Type> > tresult
199 patch_.nPoints(), pTraits<Type>::zero
203 Field<Type>& result = tresult();
205 const labelListList& pointFaces = patch_.pointFaces();
206 const scalarListList& weights = faceToPointWeights();
208 forAll(pointFaces, pointi)
210 const labelList& curFaces = pointFaces[pointi];
211 const scalarList& w = weights[pointi];
213 forAll(curFaces, facei)
215 result[pointi] += w[facei]*ff[curFaces[facei]];
223 template<class Patch>
225 tmp<Field<Type> > PrimitivePatchInterpolation<Patch>::faceToPointInterpolate
227 const tmp<Field<Type> >& tff
230 tmp<Field<Type> > tint = faceToPointInterpolate(tff());
236 template<class Patch>
238 tmp<Field<Type> > PrimitivePatchInterpolation<Patch>::pointToFaceInterpolate
240 const Field<Type>& pf
243 if (pf.size() != patch_.nPoints())
247 "tmp<Field<Type> > PrimitivePatchInterpolation::"
248 "pointToFaceInterpolate(const Field<Type> pf)"
249 ) << "given field does not correspond to patch. Patch size: "
250 << patch_.nPoints() << " field size: " << pf.size()
251 << abort(FatalError);
254 tmp<Field<Type> > tresult
263 Field<Type>& result = tresult();
265 const faceList& localFaces = patch_.localFaces();
267 forAll(result, facei)
269 const labelList& curPoints = localFaces[facei];
271 forAll(curPoints, pointi)
273 result[facei] += pf[curPoints[pointi]];
276 result[facei] /= curPoints.size();
283 template<class Patch>
285 tmp<Field<Type> > PrimitivePatchInterpolation<Patch>::pointToFaceInterpolate
287 const tmp<Field<Type> >& tpf
290 tmp<Field<Type> > tint = pointToFaceInterpolate(tpf());
296 template<class Patch>
298 tmp<Field<Type> > PrimitivePatchInterpolation<Patch>::faceToEdgeInterpolate
300 const Field<Type>& pf
303 // Check size of the given field
304 if (pf.size() != patch_.size())
308 "tmp<Field<Type> > PrimitivePatchInterpolation::"
309 "faceToEdgeInterpolate(const Field<Type> ff)"
310 ) << "given field does not correspond to patch. Patch size: "
311 << patch_.size() << " field size: " << pf.size()
312 << abort(FatalError);
315 tmp<Field<Type> > tresult
317 new Field<Type>(patch_.nEdges(), pTraits<Type>::zero)
320 Field<Type>& result = tresult();
322 const pointField& points = patch_.localPoints();
323 const faceList& faces = patch_.localFaces();
324 const edgeList& edges = patch_.edges();
325 const labelListList& edgeFaces = patch_.edgeFaces();
327 const scalarList& weights = faceToEdgeWeights();
329 for (label edgei = 0; edgei < patch_.nInternalEdges(); edgei++)
332 weights[edgei]*pf[edgeFaces[edgei][0]]
333 + (1.0 - weights[edgei])*pf[edgeFaces[edgei][1]];
336 for (label edgei = patch_.nInternalEdges(); edgei < edges.size(); edgei++)
338 result[edgei] = pf[edgeFaces[edgei][0]];
345 template<class Patch>
347 tmp<Field<Type> > PrimitivePatchInterpolation<Patch>::faceToEdgeInterpolate
349 const tmp<Field<Type> >& tpf
352 tmp<Field<Type> > tint = faceToEdgeInterpolate(tpf());
358 template<class Patch>
359 bool PrimitivePatchInterpolation<Patch>::movePoints()
367 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
369 } // End namespace Foam
371 // ************************************************************************* //