pointToPointPlanarInterpolation.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) 2012-2016 OpenFOAM Foundation
9  Copyright (C) 2019-2020 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 #include "boundBox.H"
31 #include "Random.H"
32 #include "vector2D.H"
33 #include "triSurface.H"
34 #include "triSurfaceTools.H"
35 #include "OBJstream.H"
36 #include "Time.H"
37 #include "matchPoints.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43  defineTypeNameAndDebug(pointToPointPlanarInterpolation, 0);
44 }
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
50 Foam::pointToPointPlanarInterpolation::calcCoordinateSystem
51 (
52  const pointField& points
53 ) const
54 {
55  if (points.size() < 3)
56  {
58  << "Only " << points.size() << " provided." << nl
59  << "Need at least three non-colinear points"
60  << " to be able to interpolate."
61  << exit(FatalError);
62  }
63 
64  const point& p0 = points[0];
65 
66  // Find furthest away point
67  vector e1;
68  label index1 = -1;
69  scalar maxDist = ROOTVSMALL;
70 
71  for (label i = 1; i < points.size(); i++)
72  {
73  const vector d = points[i] - p0;
74  scalar magD = mag(d);
75 
76  if (magD > maxDist)
77  {
78  e1 = d/magD;
79  index1 = i;
80  maxDist = magD;
81  }
82  }
83 
84  if (index1 == -1)
85  {
87  << "Cannot find any point that is different from first point"
88  << p0 << ". Are all your points coincident?"
89  << exit(FatalError);
90  }
91 
92 
93  // Find point that is furthest away from line p0-p1
94  const point& p1 = points[index1];
95 
96  label index2 = -1;
97  maxDist = ROOTVSMALL;
98  for (label i = 1; i < points.size(); i++)
99  {
100  if (i != index1)
101  {
102  const point& p2 = points[i];
103  vector e2(p2 - p0);
104  e2 -= (e2&e1)*e1;
105  scalar magE2 = mag(e2);
106 
107  if (magE2 > maxDist)
108  {
109  index2 = i;
110  maxDist = magE2;
111  }
112  }
113  }
114  if (index2 == -1)
115  {
117  << "Cannot find points that define a plane with a valid normal."
118  << nl << "Have so far points " << p0 << " and " << p1
119  << ". Are all your points on a single line instead of a plane?"
120  << exit(FatalError);
121  }
122 
123  const vector n = normalised(e1 ^ (points[index2]-p0));
124 
126  << " Used points " << p0 << ' ' << points[index1]
127  << ' ' << points[index2]
128  << " to define coordinate system with normal " << n << endl;
129 
130  return coordSystem::cartesian
131  (
132  p0, // origin
133  n, // normal
134  e1 // 0-axis
135  );
136 }
137 
138 
139 void Foam::pointToPointPlanarInterpolation::calcWeights
140 (
141  const pointField& sourcePoints,
142  const pointField& destPoints
143 )
144 {
145  if (nearestOnly_)
146  {
147  labelList destToSource;
148  bool fullMatch = matchPoints
149  (
150  destPoints,
151  sourcePoints,
152  scalarField(destPoints.size(), GREAT),
153  true, // verbose
154  destToSource
155  );
156 
157  if (!fullMatch)
158  {
160  << "Did not find a corresponding sourcePoint for every face"
161  << " centre" << exit(FatalError);
162  }
163 
164  nearestVertex_.setSize(destPoints.size());
165  nearestVertexWeight_.setSize(destPoints.size());
166  forAll(nearestVertex_, i)
167  {
168  nearestVertex_[i][0] = destToSource[i];
169  nearestVertex_[i][1] = -1;
170  nearestVertex_[i][2] = -1;
171 
172  nearestVertexWeight_[i][0] = 1.0;
173  nearestVertexWeight_[i][1] = 0.0;
174  nearestVertexWeight_[i][2] = 0.0;
175  }
176 
177  if (debug)
178  {
179  forAll(destPoints, i)
180  {
181  label v0 = nearestVertex_[i][0];
182 
183  Pout<< "For location " << destPoints[i]
184  << " sampling vertex " << v0
185  << " at:" << sourcePoints[v0]
186  << " distance:" << mag(sourcePoints[v0]-destPoints[i])
187  << endl;
188  }
189 
190  OBJstream str("destToSource.obj");
191  Pout<< "pointToPointPlanarInterpolation::calcWeights :"
192  << " Dumping lines from face centres to original points to "
193  << str.name() << endl;
194 
195  forAll(destPoints, i)
196  {
197  label v0 = nearestVertex_[i][0];
198  str.write(linePointRef(destPoints[i], sourcePoints[v0]));
199  }
200  }
201  }
202  else
203  {
204  auto tlocalVertices = referenceCS_.localPosition(sourcePoints);
205  auto& localVertices = tlocalVertices.ref();
206 
207  const boundBox bb(localVertices, true);
208  const point bbMid(bb.centre());
209 
211  << " Perturbing points with " << perturb_
212  << " fraction of a random position inside " << bb
213  << " to break any ties on regular meshes." << nl
214  << endl;
215 
216  Random rndGen(123456);
217  forAll(localVertices, i)
218  {
219  localVertices[i] +=
220  perturb_
221  *(rndGen.position(bb.min(), bb.max())-bbMid);
222  }
223 
224  // Determine triangulation
225  List<vector2D> localVertices2D(localVertices.size());
226  forAll(localVertices, i)
227  {
228  localVertices2D[i][0] = localVertices[i][0];
229  localVertices2D[i][1] = localVertices[i][1];
230  }
231 
232  triSurface s(triSurfaceTools::delaunay2D(localVertices2D));
233 
234  auto tlocalFaceCentres = referenceCS_.localPosition(destPoints);
235  const pointField& localFaceCentres = tlocalFaceCentres();
236 
237  if (debug)
238  {
239  Pout<< "pointToPointPlanarInterpolation::calcWeights :"
240  <<" Dumping triangulated surface to triangulation.stl" << endl;
241  s.write("triangulation.stl");
242 
243  OBJstream str("localFaceCentres.obj");
244  Pout<< "pointToPointPlanarInterpolation::calcWeights :"
245  << " Dumping face centres to " << str.name() << endl;
246 
247  forAll(localFaceCentres, i)
248  {
249  str.write(localFaceCentres[i]);
250  }
251  }
252 
253  // Determine interpolation onto face centres.
255  (
256  s,
257  localFaceCentres, // points to interpolate to
258  nearestVertex_,
259  nearestVertexWeight_
260  );
261 
262  if (debug)
263  {
264  forAll(sourcePoints, i)
265  {
266  Pout<< "source:" << i << " at:" << sourcePoints[i]
267  << " 2d:" << localVertices[i]
268  << endl;
269  }
270 
271  OBJstream str("stencil.obj");
272  Pout<< "pointToPointPlanarInterpolation::calcWeights :"
273  << " Dumping stencil to " << str.name() << endl;
274 
275  forAll(destPoints, i)
276  {
277  label v0 = nearestVertex_[i][0];
278  label v1 = nearestVertex_[i][1];
279  label v2 = nearestVertex_[i][2];
280 
281  Pout<< "For location " << destPoints[i]
282  << " 2d:" << localFaceCentres[i]
283  << " sampling vertices" << nl
284  << " " << v0
285  << " at:" << sourcePoints[v0]
286  << " weight:" << nearestVertexWeight_[i][0] << nl;
287 
288  str.write(linePointRef(destPoints[i], sourcePoints[v0]));
289 
290  if (v1 != -1)
291  {
292  Pout<< " " << v1
293  << " at:" << sourcePoints[v1]
294  << " weight:" << nearestVertexWeight_[i][1] << nl;
295  str.write(linePointRef(destPoints[i], sourcePoints[v1]));
296  }
297  if (v2 != -1)
298  {
299  Pout<< " " << v2
300  << " at:" << sourcePoints[v2]
301  << " weight:" << nearestVertexWeight_[i][2] << nl;
302  str.write(linePointRef(destPoints[i], sourcePoints[v2]));
303  }
304 
305  Pout<< endl;
306  }
307  }
308  }
309 }
310 
311 
312 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
313 
315 (
316  const pointField& sourcePoints,
317  const pointField& destPoints,
318  const scalar perturb,
319  const bool nearestOnly
320 )
321 :
322  perturb_(perturb),
323  nearestOnly_(nearestOnly),
324  referenceCS_(),
325  nPoints_(sourcePoints.size())
326 {
327  if (!nearestOnly)
328  {
329  referenceCS_ = calcCoordinateSystem(sourcePoints);
330  }
331  calcWeights(sourcePoints, destPoints);
332 }
333 
334 
336 (
337  const coordinateSystem& referenceCS,
338  const pointField& sourcePoints,
339  const pointField& destPoints,
340  const scalar perturb
341 )
342 :
343  perturb_(perturb),
344  nearestOnly_(false),
345  referenceCS_(referenceCS),
346  nPoints_(sourcePoints.size())
347 {
348  calcWeights(sourcePoints, destPoints);
349 }
350 
351 
353 (
354  const scalar perturb,
355  const bool nearestOnly,
356  const coordinateSystem& referenceCS,
357  const label sourceSize,
358  const List<FixedList<label, 3>>& nearestVertex,
359  const List<FixedList<scalar, 3>>& nearestVertexWeight
360 )
361 :
362  perturb_(perturb),
363  nearestOnly_(nearestOnly),
364  referenceCS_(referenceCS),
365  nPoints_(sourceSize),
366  nearestVertex_(nearestVertex),
367  nearestVertexWeight_(nearestVertexWeight)
368 {}
369 
370 
373 {
375  (
376  perturb_,
377  nearestOnly_,
378  referenceCS_,
379  nPoints_,
380  nearestVertex_,
381  nearestVertexWeight_
382  );
383 }
384 
385 
386 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
387 
389 (
390  const instantList& times
391 )
392 {
393  wordList names(times.size());
394 
395  forAll(times, i)
396  {
397  names[i] = times[i].name();
398  }
399  return names;
400 }
401 
402 
404 (
405  const instantList& times,
406  const label startSampleTime,
407  const scalar timeVal,
408  label& lo,
409  label& hi
410 )
411 {
412  lo = startSampleTime;
413  hi = -1;
414 
415  for (label i = startSampleTime+1; i < times.size(); i++)
416  {
417  if (times[i].value() > timeVal)
418  {
419  break;
420  }
421  else
422  {
423  lo = i;
424  }
425  }
426 
427  if (lo == -1)
428  {
429  //FatalErrorInFunction
430  // << "Cannot find starting sampling values for current time "
431  // << timeVal << nl
432  // << "Have sampling values for times "
433  // << timeNames(times) << nl
434  // << exit(FatalError);
435  return false;
436  }
437 
438  if (lo < times.size()-1)
439  {
440  hi = lo+1;
441  }
442 
443 
444  if (debug)
445  {
446  if (hi == -1)
447  {
448  Pout<< "findTime : Found time " << timeVal << " after"
449  << " index:" << lo << " time:" << times[lo].value()
450  << endl;
451  }
452  else
453  {
454  Pout<< "findTime : Found time " << timeVal << " inbetween"
455  << " index:" << lo << " time:" << times[lo].value()
456  << " and index:" << hi << " time:" << times[hi].value()
457  << endl;
458  }
459  }
460  return true;
461 }
462 
463 
464 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::autoPtr::New
static autoPtr< T > New(Args &&... args)
Construct autoPtr of T with forwarding arguments.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputSpray.H:25
Foam::triSurfaceTools::delaunay2D
static triSurface delaunay2D(const List< vector2D > &)
Do unconstrained Delaunay of points. Returns triSurface with 3D.
Definition: triSurfaceTools.C:2472
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
triSurface.H
Foam::pointToPointPlanarInterpolation::findTime
static bool findTime(const instantList &times, const label startSampleTime, const scalar timeVal, label &lo, label &hi)
Helper: find time. Return true if successful.
Definition: pointToPointPlanarInterpolation.C:404
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
n
label n
Definition: TABSMDCalcMethod2.H:31
matchPoints.H
Determine correspondence between points. See below.
Foam::Field< vector >
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:388
Foam::OSstream::name
virtual const fileName & name() const
Get the name of the stream.
Definition: OSstream.H:107
pointToPointPlanarInterpolation.H
vector2D.H
Foam::triSurfaceTools::calcInterpolationWeights
static void calcInterpolationWeights(const triPointRef &tri, const point &p, FixedList< scalar, 3 > &weights)
Calculate linear interpolation weights for point (guaranteed to be.
Definition: triSurfaceTools.C:2538
Foam::FatalError
error FatalError
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Random.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::normalised
VectorSpace< Form, Cmpt, Ncmpts > normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:487
boundBox.H
Time.H
Foam::autoPtr< Foam::pointToPointPlanarInterpolation >
Foam::linePointRef
line< point, const point & > linePointRef
A line using referred points.
Definition: linePointRef.H:47
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::coordSystem::cartesian
A Cartesian coordinate system.
Definition: cartesianCS.H:69
Foam::Random::position
Type position(const Type &start, const Type &end)
Return a sample on the interval [start,end].
Definition: RandomTemplates.C:62
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::FixedList< label, 3 >
points
const pointField & points
Definition: gmvOutputHeader.H:1
rndGen
Random rndGen
Definition: createFields.H:23
Foam::pointToPointPlanarInterpolation::clone
autoPtr< pointToPointPlanarInterpolation > clone() const
Construct and return a clone.
Definition: pointToPointPlanarInterpolation.C:372
p0
const volScalarField & p0
Definition: EEqn.H:36
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::PtrListOps::names
List< word > names(const UPtrList< T > &list, const UnaryMatchPredicate &matcher)
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
triSurfaceTools.H
Foam::matchPoints
bool matchPoints(const UList< point > &pts0, const UList< point > &pts1, const UList< scalar > &matchDistance, const bool verbose, labelList &from0To1, const point &origin=point::zero)
Determine correspondence between pointFields. Gets passed.
Definition: matchPoints.C:36
OBJstream.H
Foam::pointToPointPlanarInterpolation::timeNames
static wordList timeNames(const instantList &)
Helper: extract words of times.
Definition: pointToPointPlanarInterpolation.C:389
Foam::coordinateSystem
Base class for coordinate system specification, the default coordinate system type is cartesian .
Definition: coordinateSystem.H:132
Foam::pointToPointPlanarInterpolation::pointToPointPlanarInterpolation
pointToPointPlanarInterpolation(const pointField &sourcePoints, const pointField &destPoints, const scalar perturb, const bool nearestOnly=false)
Construct from 3D locations. Determines local coordinate system.
Definition: pointToPointPlanarInterpolation.C:315