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 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 
125  if (debug)
126  {
128  << " Used points " << p0 << ' ' << points[index1]
129  << ' ' << points[index2]
130  << " to define coordinate system with normal " << n << endl;
131  }
132 
133  return coordSystem::cartesian
134  (
135  p0, // origin
136  n, // normal
137  e1 // 0-axis
138  );
139 }
140 
141 
142 void Foam::pointToPointPlanarInterpolation::calcWeights
143 (
144  const pointField& sourcePoints,
145  const pointField& destPoints
146 )
147 {
148  if (nearestOnly_)
149  {
150  labelList destToSource;
151  bool fullMatch = matchPoints
152  (
153  destPoints,
154  sourcePoints,
155  scalarField(destPoints.size(), GREAT),
156  true, // verbose
157  destToSource
158  );
159 
160  if (!fullMatch)
161  {
163  << "Did not find a corresponding sourcePoint for every face"
164  << " centre" << exit(FatalError);
165  }
166 
167  nearestVertex_.setSize(destPoints.size());
168  nearestVertexWeight_.setSize(destPoints.size());
169  forAll(nearestVertex_, i)
170  {
171  nearestVertex_[i][0] = destToSource[i];
172  nearestVertex_[i][1] = -1;
173  nearestVertex_[i][2] = -1;
174 
175  nearestVertexWeight_[i][0] = 1.0;
176  nearestVertexWeight_[i][1] = 0.0;
177  nearestVertexWeight_[i][2] = 0.0;
178  }
179 
180  if (debug)
181  {
182  forAll(destPoints, i)
183  {
184  label v0 = nearestVertex_[i][0];
185 
186  Pout<< "For location " << destPoints[i]
187  << " sampling vertex " << v0
188  << " at:" << sourcePoints[v0]
189  << " distance:" << mag(sourcePoints[v0]-destPoints[i])
190  << endl;
191  }
192 
193  OBJstream str("destToSource.obj");
194  Pout<< "pointToPointPlanarInterpolation::calcWeights :"
195  << " Dumping lines from face centres to original points to "
196  << str.name() << endl;
197 
198  forAll(destPoints, i)
199  {
200  label v0 = nearestVertex_[i][0];
201  str.write(linePointRef(destPoints[i], sourcePoints[v0]));
202  }
203  }
204  }
205  else
206  {
207  auto tlocalVertices = referenceCS_.localPosition(sourcePoints);
208  auto& localVertices = tlocalVertices.ref();
209 
210  const boundBox bb(localVertices, true);
211  const point bbMid(bb.centre());
212 
213  if (debug)
214  {
216  << " Perturbing points with " << perturb_
217  << " fraction of a random position inside " << bb
218  << " to break any ties on regular meshes."
219  << nl << endl;
220  }
221 
222  Random rndGen(123456);
223  forAll(localVertices, i)
224  {
225  localVertices[i] +=
226  perturb_
227  *(rndGen.position(bb.min(), bb.max())-bbMid);
228  }
229 
230  // Determine triangulation
231  List<vector2D> localVertices2D(localVertices.size());
232  forAll(localVertices, i)
233  {
234  localVertices2D[i][0] = localVertices[i][0];
235  localVertices2D[i][1] = localVertices[i][1];
236  }
237 
238  triSurface s(triSurfaceTools::delaunay2D(localVertices2D));
239 
240  auto tlocalFaceCentres = referenceCS_.localPosition(destPoints);
241  const pointField& localFaceCentres = tlocalFaceCentres();
242 
243  if (debug)
244  {
245  Pout<< "pointToPointPlanarInterpolation::calcWeights :"
246  <<" Dumping triangulated surface to triangulation.stl" << endl;
247  s.write("triangulation.stl");
248 
249  OBJstream str("localFaceCentres.obj");
250  Pout<< "pointToPointPlanarInterpolation::calcWeights :"
251  << " Dumping face centres to " << str.name() << endl;
252 
253  forAll(localFaceCentres, i)
254  {
255  str.write(localFaceCentres[i]);
256  }
257  }
258 
259  // Determine interpolation onto face centres.
261  (
262  s,
263  localFaceCentres, // points to interpolate to
264  nearestVertex_,
265  nearestVertexWeight_
266  );
267 
268  if (debug)
269  {
270  forAll(sourcePoints, i)
271  {
272  Pout<< "source:" << i << " at:" << sourcePoints[i]
273  << " 2d:" << localVertices[i]
274  << endl;
275  }
276 
277  OBJstream str("stencil.obj");
278  Pout<< "pointToPointPlanarInterpolation::calcWeights :"
279  << " Dumping stencil to " << str.name() << endl;
280 
281  forAll(destPoints, i)
282  {
283  label v0 = nearestVertex_[i][0];
284  label v1 = nearestVertex_[i][1];
285  label v2 = nearestVertex_[i][2];
286 
287  Pout<< "For location " << destPoints[i]
288  << " 2d:" << localFaceCentres[i]
289  << " sampling vertices" << nl
290  << " " << v0
291  << " at:" << sourcePoints[v0]
292  << " weight:" << nearestVertexWeight_[i][0] << nl;
293 
294  str.write(linePointRef(destPoints[i], sourcePoints[v0]));
295 
296  if (v1 != -1)
297  {
298  Pout<< " " << v1
299  << " at:" << sourcePoints[v1]
300  << " weight:" << nearestVertexWeight_[i][1] << nl;
301  str.write(linePointRef(destPoints[i], sourcePoints[v1]));
302  }
303  if (v2 != -1)
304  {
305  Pout<< " " << v2
306  << " at:" << sourcePoints[v2]
307  << " weight:" << nearestVertexWeight_[i][2] << nl;
308  str.write(linePointRef(destPoints[i], sourcePoints[v2]));
309  }
310 
311  Pout<< endl;
312  }
313  }
314  }
315 }
316 
317 
318 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
319 
321 (
322  const pointField& sourcePoints,
323  const pointField& destPoints,
324  const scalar perturb,
325  const bool nearestOnly
326 )
327 :
328  perturb_(perturb),
329  nearestOnly_(nearestOnly),
330  referenceCS_(),
331  nPoints_(sourcePoints.size())
332 {
333  if (!nearestOnly)
334  {
335  referenceCS_ = calcCoordinateSystem(sourcePoints);
336  }
337  calcWeights(sourcePoints, destPoints);
338 }
339 
340 
342 (
343  const coordinateSystem& referenceCS,
344  const pointField& sourcePoints,
345  const pointField& destPoints,
346  const scalar perturb
347 )
348 :
349  perturb_(perturb),
350  nearestOnly_(false),
351  referenceCS_(referenceCS),
352  nPoints_(sourcePoints.size())
353 {
354  calcWeights(sourcePoints, destPoints);
355 }
356 
357 
359 (
360  const scalar perturb,
361  const bool nearestOnly,
362  const coordinateSystem& referenceCS,
363  const label sourceSize,
364  const List<FixedList<label, 3>>& nearestVertex,
365  const List<FixedList<scalar, 3>>& nearestVertexWeight
366 )
367 :
368  perturb_(perturb),
369  nearestOnly_(nearestOnly),
370  referenceCS_(referenceCS),
371  nPoints_(sourceSize),
372  nearestVertex_(nearestVertex),
373  nearestVertexWeight_(nearestVertexWeight)
374 {}
375 
376 
379 {
381  (
382  perturb_,
383  nearestOnly_,
384  referenceCS_,
385  nPoints_,
386  nearestVertex_,
387  nearestVertexWeight_
388  );
389 }
390 
391 
392 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
393 
395 (
396  const instantList& times
397 )
398 {
399  wordList names(times.size());
400 
401  forAll(times, i)
402  {
403  names[i] = times[i].name();
404  }
405  return names;
406 }
407 
408 
410 (
411  const instantList& times,
412  const label startSampleTime,
413  const scalar timeVal,
414  label& lo,
415  label& hi
416 )
417 {
418  lo = startSampleTime;
419  hi = -1;
420 
421  for (label i = startSampleTime+1; i < times.size(); i++)
422  {
423  if (times[i].value() > timeVal)
424  {
425  break;
426  }
427  else
428  {
429  lo = i;
430  }
431  }
432 
433  if (lo == -1)
434  {
435  //FatalErrorInFunction
436  // << "Cannot find starting sampling values for current time "
437  // << timeVal << nl
438  // << "Have sampling values for times "
439  // << timeNames(times) << nl
440  // << exit(FatalError);
441  return false;
442  }
443 
444  if (lo < times.size()-1)
445  {
446  hi = lo+1;
447  }
448 
449 
450  if (debug)
451  {
452  if (hi == -1)
453  {
454  Pout<< "findTime : Found time " << timeVal << " after"
455  << " index:" << lo << " time:" << times[lo].value()
456  << endl;
457  }
458  else
459  {
460  Pout<< "findTime : Found time " << timeVal << " inbetween"
461  << " index:" << lo << " time:" << times[lo].value()
462  << " and index:" << hi << " time:" << times[hi].value()
463  << endl;
464  }
465  }
466  return true;
467 }
468 
469 
470 // ************************************************************************* //
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:74
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
InfoInFunction
#define InfoInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:316
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:2475
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
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:410
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
n
label n
Definition: TABSMDCalcMethod2.H:31
matchPoints.H
Determine correspondence between points. See below.
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< vector >
Foam::OSstream::name
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:91
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:2541
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:496
boundBox.H
Time.H
Foam::autoPtr< Foam::pointToPointPlanarInterpolation >
Foam::linePointRef
line< point, const point & > linePointRef
Line using referred points.
Definition: linePointRef.H:47
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::nl
constexpr char nl
Definition: Ostream.H:372
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: HashTable.H:102
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:378
p0
const volScalarField & p0
Definition: EEqn.H:36
Foam::point
vector point
Point is a vector.
Definition: point.H:43
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:35
OBJstream.H
Foam::pointToPointPlanarInterpolation::timeNames
static wordList timeNames(const instantList &)
Helper: extract words of times.
Definition: pointToPointPlanarInterpolation.C:395
Foam::coordinateSystem
Base class for coordinate system specification, the default coordinate system type is cartesian .
Definition: coordinateSystem.H:122
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:321