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-2022 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 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
48 
50 Foam::pointToPointPlanarInterpolation::calcCoordinateSystem
51 (
52  const pointField& points
53 )
54 {
55  if (points.size() < 3)
56  {
58  << "Need at least 3 non-collinear points for planar interpolation,"
59  << " but only had " << points.size() << " points" << nl
60  << exit(FatalError);
61  }
62 
63  const point& p0 = points[0];
64 
65  // Find furthest away point
66  label index1 = -1;
67  scalar maxDistSqr = ROOTVSMALL;
68 
69  for (label i = 1; i < points.size(); ++i)
70  {
71  const scalar mag2 = magSqr(points[i] - p0);
72 
73  if (maxDistSqr < mag2)
74  {
75  maxDistSqr = mag2;
76  index1 = i;
77  }
78  }
79  if (index1 == -1)
80  {
82  << "Cannot find any point that is different from first point"
83  << p0 << ". Are all your points coincident?"
84  << exit(FatalError);
85  }
86 
87  const vector e1(normalised(points[index1] - p0));
88 
89  // Find point that is furthest perpendicular distance from the p0-p1 line
90  label index2 = -1;
91  maxDistSqr = ROOTVSMALL;
92  for (label i = 1; i < points.size(); i++)
93  {
94  if (i != index1)
95  {
96  vector e2(points[i] - p0);
97  e2.removeCollinear(e1);
98 
99  const scalar mag2 = magSqr(e2);
100 
101  if (maxDistSqr < mag2)
102  {
103  maxDistSqr = mag2;
104  index2 = i;
105  }
106  }
107  }
108  if (index2 == -1)
109  {
111  << "Cannot find points that define a plane with a valid normal."
112  << nl << "Have so far points " << p0 << " and " << points[index1]
113  << ". Are all your points on a single line instead of a plane?"
114  << exit(FatalError);
115  }
116 
117  const vector n = normalised(e1 ^ (points[index2]-p0));
118 
120  << " Used points "
121  << p0 << ' ' << points[index1] << ' ' << points[index2]
122  << " to define coordinate system with normal " << n << endl;
123 
124  return coordSystem::cartesian
125  (
126  p0, // origin
127  n, // normal
128  e1 // 0-axis
129  );
130 }
131 
132 
133 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
134 
135 void Foam::pointToPointPlanarInterpolation::calcWeights
136 (
137  const pointField& sourcePoints,
138  const pointField& destPoints
139 )
140 {
141  if (nearestOnly_)
142  {
143  labelList destToSource;
144  bool fullMatch = matchPoints
145  (
146  destPoints,
147  sourcePoints,
148  scalarField(destPoints.size(), GREAT),
149  true, // verbose
150  destToSource
151  );
152 
153  if (!fullMatch)
154  {
156  << "Did not find a corresponding sourcePoint for every face"
157  << " centre" << exit(FatalError);
158  }
159 
160  nearestVertex_.resize(destPoints.size());
161  nearestVertexWeight_.resize(destPoints.size());
162  forAll(nearestVertex_, i)
163  {
164  nearestVertex_[i][0] = destToSource[i];
165  nearestVertex_[i][1] = -1;
166  nearestVertex_[i][2] = -1;
167 
168  nearestVertexWeight_[i][0] = 1.0;
169  nearestVertexWeight_[i][1] = 0.0;
170  nearestVertexWeight_[i][2] = 0.0;
171  }
172 
173  if (debug)
174  {
175  forAll(destPoints, i)
176  {
177  label v0 = nearestVertex_[i][0];
178 
179  Pout<< "For location " << destPoints[i]
180  << " sampling vertex " << v0
181  << " at:" << sourcePoints[v0]
182  << " distance:" << mag(sourcePoints[v0]-destPoints[i])
183  << endl;
184  }
185 
186  OBJstream str
187  (
188  "destToSource_" + Foam::name(UPstream::myProcNo()) + ".obj"
189  );
190  Pout<< "pointToPointPlanarInterpolation::calcWeights :"
191  << " Dumping lines from face centres to original points to "
192  << str.name() << endl;
193 
194  forAll(destPoints, i)
195  {
196  const label v0 = nearestVertex_[i][0];
197  str.writeLine(destPoints[i], sourcePoints[v0]);
198  }
199  }
200  }
201  else
202  {
203  auto tlocalVertices = referenceCS_.localPosition(sourcePoints);
204  auto& localVertices = tlocalVertices.ref();
205 
206  const boundBox bb(localVertices, true);
207  const point bbMid(bb.centre());
208 
210  << " Perturbing points with " << perturb_
211  << " fraction of a random position inside " << bb
212  << " to break any ties on regular meshes." << nl
213  << endl;
214 
215  Random rndGen(123456);
216  forAll(localVertices, i)
217  {
218  localVertices[i] +=
219  perturb_
220  *(rndGen.position(bb.min(), bb.max())-bbMid);
221  }
222 
223  // Determine triangulation
224  List<vector2D> localVertices2D(localVertices.size());
225  forAll(localVertices, i)
226  {
227  localVertices2D[i][0] = localVertices[i][0];
228  localVertices2D[i][1] = localVertices[i][1];
229  }
230 
231  triSurface s(triSurfaceTools::delaunay2D(localVertices2D));
232 
233  auto tlocalFaceCentres = referenceCS_.localPosition(destPoints);
234  const pointField& localFaceCentres = tlocalFaceCentres();
235 
236  if (debug)
237  {
238  fileName outName
239  (
240  "triangulation_" + Foam::name(UPstream::myProcNo()) + ".obj"
241  );
242 
243  Pout<< "pointToPointPlanarInterpolation::calcWeights :"
244  <<" Dumping triangulated surface to " << outName << endl;
245 
246  s.write(outName);
247 
248  OBJstream os
249  (
250  "localFaceCentres_" + Foam::name(UPstream::myProcNo()) + ".obj"
251  );
252  Pout<< "pointToPointPlanarInterpolation::calcWeights :"
253  << " Dumping face centres to " << os.name() << endl;
254 
255  os.write(localFaceCentres);
256  }
257 
258  // Determine interpolation onto face centres.
260  (
261  s,
262  localFaceCentres, // points to interpolate to
263  nearestVertex_,
264  nearestVertexWeight_
265  );
266 
267  if (debug)
268  {
269  forAll(sourcePoints, i)
270  {
271  Pout<< "source:" << i << " at:" << sourcePoints[i]
272  << " 2d:" << localVertices[i]
273  << endl;
274  }
275 
276  OBJstream str
277  (
278  "stencil_" + Foam::name(UPstream::myProcNo()) + ".obj"
279  );
280  Pout<< "pointToPointPlanarInterpolation::calcWeights :"
281  << " Dumping stencil to " << str.name() << endl;
282 
283  forAll(destPoints, i)
284  {
285  const label v0 = nearestVertex_[i][0];
286  const label v1 = nearestVertex_[i][1];
287  const label v2 = nearestVertex_[i][2];
288 
289  Pout<< "For location " << destPoints[i]
290  << " 2d:" << localFaceCentres[i]
291  << " sampling vertices" << nl
292  << " " << v0
293  << " at:" << sourcePoints[v0]
294  << " weight:" << nearestVertexWeight_[i][0] << nl;
295 
296  str.writeLine(destPoints[i], sourcePoints[v0]);
297 
298  if (v1 != -1)
299  {
300  Pout<< " " << v1
301  << " at:" << sourcePoints[v1]
302  << " weight:" << nearestVertexWeight_[i][1] << nl;
303  str.writeLine(destPoints[i], sourcePoints[v1]);
304  }
305  if (v2 != -1)
306  {
307  Pout<< " " << v2
308  << " at:" << sourcePoints[v2]
309  << " weight:" << nearestVertexWeight_[i][2] << nl;
310  str.writeLine(destPoints[i], sourcePoints[v2]);
311  }
312 
313  Pout<< endl;
314  }
315  }
316  }
317 }
318 
319 
320 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
321 
323 (
324  const pointField& sourcePoints,
325  const pointField& destPoints,
326  const scalar perturb,
327  const bool nearestOnly
328 )
329 :
330  perturb_(perturb),
331  nearestOnly_(nearestOnly),
332  referenceCS_(),
333  nPoints_(sourcePoints.size())
334 {
335  if (!nearestOnly)
336  {
337  referenceCS_ = calcCoordinateSystem(sourcePoints);
338  }
339  calcWeights(sourcePoints, destPoints);
340 }
341 
342 
344 (
345  const coordinateSystem& referenceCS,
346  const pointField& sourcePoints,
347  const pointField& destPoints,
348  const scalar perturb
349 )
350 :
351  perturb_(perturb),
352  nearestOnly_(false),
353  referenceCS_(referenceCS),
354  nPoints_(sourcePoints.size())
355 {
356  calcWeights(sourcePoints, destPoints);
357 }
358 
359 
361 (
362  const scalar perturb,
363  const bool nearestOnly,
364  const coordinateSystem& referenceCS,
365  const label sourceSize,
366  List<FixedList<label, 3>>&& nearestVertex,
367  List<FixedList<scalar, 3>>&& nearestVertexWeight
368 )
369 :
370  perturb_(perturb),
371  nearestOnly_(nearestOnly),
372  referenceCS_(referenceCS),
373  nPoints_(sourceSize),
374  nearestVertex_(std::move(nearestVertex)),
375  nearestVertexWeight_(std::move(nearestVertexWeight))
376 {}
377 
378 
379 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
380 
382 (
383  const instantList& times
384 )
385 {
386  wordList names(times.size());
387 
388  forAll(times, i)
389  {
390  names[i] = times[i].name();
391  }
392  return names;
393 }
394 
395 
396 // ************************************************************************* //
Base class for coordinate system specification, the default coordinate system type is cartesian ...
static wordList timeNames(const instantList &times)
Helper: extract words of times.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
List< word > names(const UPtrList< T > &list, const UnaryMatchPredicate &matcher)
List of names generated by calling name() for each list item and filtered for matches.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
virtual Ostream & write(const char c) override
Write character.
Definition: OBJstream.C:69
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
virtual const fileName & name() const override
Read/write access to the name of the stream.
Definition: OSstream.H:128
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Random rndGen
Definition: createFields.H:23
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:674
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1074
virtual const fileName & name() const override
Get the name of the output serial stream. (eg, the name of the Fstream file name) ...
Definition: OSstream.H:128
static triSurface delaunay2D(const List< vector2D > &)
Do unconstrained Delaunay of points. Returns triSurface with 3D.
Determine correspondence between points. See below.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
static void calcInterpolationWeights(const triPointRef &tri, const point &p, FixedList< scalar, 3 > &weights)
Calculate linear interpolation weights for point (guaranteed to be.
A Cartesian coordinate system.
Definition: cartesianCS.H:65
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:29
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const pointField & points
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
#define DebugInFunction
Report an information message using Foam::Info.
Vector< scalar > vector
Definition: vector.H:57
int debug
Static debugging option.
OBJstream os(runTime.globalPath()/outputName)
defineTypeNameAndDebug(combustionModel, 0)
bool nearestOnly() const noexcept
Whether to use nearest point only (avoids triangulation, projection)
vector point
Point is a vector.
Definition: point.H:37
label n
List< label > labelList
A List of labels.
Definition: List.H:62
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))
pointToPointPlanarInterpolation(const pointToPointPlanarInterpolation &)=default
Copy construct.
Type position(const Type &start, const Type &end)
Return a sample on the interval [start,end].
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
const volScalarField & p0
Definition: EEqn.H:36
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.