advancingFrontAMI.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) 2013-2016 OpenFOAM Foundation
9  Copyright (C) 2015-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 
29 #include "advancingFrontAMI.H"
30 #include "meshTools.H"
31 #include "mapDistribute.H"
32 #include "unitConversion.H"
33 
34 #include "findNearestMaskedOp.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
41 }
42 
43 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
44 
46 {
47  const auto& src = srcPatch();
48  const auto& tgt = tgtPatch();
49 
50  if (debug && (!src.size() || !tgt.size()))
51  {
52  Pout<< "AMI: Patches not on processor: Source faces = "
53  << src.size() << ", target faces = " << tgt.size()
54  << endl;
55  }
56 
57 
58  if (requireMatch_)
59  {
60  const scalar maxBoundsError = 0.05;
61 
62  // Check bounds of source and target
63  boundBox bbSrc(src.points(), src.meshPoints(), true);
64  boundBox bbTgt(tgt.points(), tgt.meshPoints(), true);
65 
66  boundBox bbTgtInf(bbTgt);
67  bbTgtInf.inflate(maxBoundsError);
68 
69  if (!bbTgtInf.contains(bbSrc))
70  {
72  << "Source and target patch bounding boxes are not similar"
73  << nl
74  << " source box span : " << bbSrc.span() << nl
75  << " target box span : " << bbTgt.span() << nl
76  << " source box : " << bbSrc << nl
77  << " target box : " << bbTgt << nl
78  << " inflated target box : " << bbTgtInf << endl;
79  }
80  }
81 }
82 
83 
85 (
86  const label srcFacei,
87  const label tgtFacei
88 ) const
89 {
90  const auto& srcPatch = this->srcPatch();
91  const auto& tgtPatch = this->tgtPatch();
92 
93  if
94  (
95  (srcMagSf_[srcFacei] < ROOTVSMALL)
96  || (tgtMagSf_[tgtFacei] < ROOTVSMALL)
97  )
98  {
99  return false;
100  }
101 
102  if (maxDistance2_ > 0)
103  {
104  const point& srcFc = srcPatch.faceCentres()[srcFacei];
105  const point& tgtFc = tgtPatch.faceCentres()[tgtFacei];
106  const vector& srcN = srcPatch.faceNormals()[srcFacei];
107 
108  const scalar normalDist((tgtFc-srcFc)&srcN);
109  //if (magSqr(srcFc-tgtFc) >= maxDistance2_)
110  if (sqr(normalDist) >= maxDistance2_)
111  {
112  return false;
113  }
114  }
115 
116  if (minCosAngle_ > -1)
117  {
118  const vector& srcN = srcPatch.faceNormals()[srcFacei];
119  vector tgtN = tgtPatch.faceNormals()[tgtFacei];
120  if (!reverseTarget_)
121  {
122  tgtN = -tgtN;
123  }
124 
125  if ((srcN & tgtN) <= minCosAngle_)
126  {
127  return false;
128  }
129  }
130 
131  return true;
132 }
133 
134 
136 {
137  // Create processor map of extended cells. This map gets (possibly
138  // remote) cells from the src mesh such that they (together) cover
139  // all of tgt
140  extendedTgtMapPtr_.reset(calcProcMap(srcPatch0(), tgtPatch0()));
141  const mapDistribute& map = extendedTgtMapPtr_();
142 
143  // Original faces from tgtPatch
144  // Note: in globalIndexing since might be remote
145  globalIndex globalTgtFaces(tgtPatch0().size());
146  distributeAndMergePatches
147  (
148  map,
149  tgtPatch0(),
150  globalTgtFaces,
151  extendedTgtFaces_,
152  extendedTgtPoints_,
153  extendedTgtFaceIDs_
154  );
155 
156  // Create a representation of the tgt patch that is extended to overlap
157  // the src patch
158  extendedTgtPatchPtr_.reset
159  (
161  (
162  SubList<face>(extendedTgtFaces_),
163  extendedTgtPoints_
164  )
165  );
166 }
167 
168 
169 bool Foam::advancingFrontAMI::initialiseWalk(label& srcFacei, label& tgtFacei)
170 {
171  const auto& src = this->srcPatch();
172  const auto& tgt = this->tgtPatch();
173 
174  bool foundFace = false;
175 
176  // Check that patch sizes are valid
177  if (!src.size())
178  {
179  return foundFace;
180  }
181  else if (!tgt.size())
182  {
184  << src.size() << " source faces but no target faces" << endl;
185 
186  return foundFace;
187  }
188 
189  // Reset the octree
190  treePtr_.reset(createTree(tgt));
191 
192  // Find initial face match using brute force/octree search
193  if ((srcFacei == -1) || (tgtFacei == -1))
194  {
195  srcFacei = 0;
196  tgtFacei = 0;
197  forAll(src, facei)
198  {
199  tgtFacei = findTargetFace(facei);
200  if (tgtFacei >= 0)
201  {
202  srcFacei = facei;
203  foundFace = true;
204  break;
205  }
206  }
207 
208  if (!foundFace)
209  {
210  if (requireMatch_)
211  {
213  << "Unable to find initial target face"
214  << abort(FatalError);
215  }
216 
217  return foundFace;
218  }
219  }
220 
221  if (debug)
222  {
223  Pout<< "AMI: initial target face = " << tgtFacei << endl;
224  }
225 
226  return true;
227 }
228 
229 
231 (
232  const scalar area,
233  const face& f1,
234  const face& f2,
235  const pointField& f1Points,
236  const pointField& f2Points
237 ) const
238 {
239  static label count = 1;
240 
241  const pointField f1pts = f1.points(f1Points);
242  const pointField f2pts = f2.points(f2Points);
243 
244  Pout<< "Face intersection area (" << count << "):" << nl
245  << " f1 face = " << f1 << nl
246  << " f1 pts = " << f1pts << nl
247  << " f2 face = " << f2 << nl
248  << " f2 pts = " << f2pts << nl
249  << " area = " << area
250  << endl;
251 
252  OFstream os("areas" + name(count) + ".obj");
253 
254  for (const point& pt : f1pts)
255  {
256  meshTools::writeOBJ(os, pt);
257  }
258  os<< "l";
259  forAll(f1pts, i)
260  {
261  os<< " " << i + 1;
262  }
263  os<< " 1" << endl;
264 
265  for (const point& pt : f2pts)
266  {
267  meshTools::writeOBJ(os, pt);
268  }
269  os<< "l";
270  const label n = f1pts.size();
271  forAll(f2pts, i)
272  {
273  os<< " " << n + i + 1;
274  }
275  os<< " " << n + 1 << endl;
276 
277  ++count;
278 }
279 
280 
282 (
283  const label srcFacei,
284  const UList<label>& excludeFaces,
285  const label srcFacePti
286 ) const
287 {
288  const auto& src = srcPatch();
289 
290  label targetFacei = -1;
291 
292  const pointField& srcPts = src.points();
293  const face& srcFace = src[srcFacei];
294 
295  findNearestMaskedOp<primitivePatch> fnOp(*treePtr_, excludeFaces);
296 
297  const boundBox bb(srcPts, srcFace, false);
298 
299  const point srcPt =
300  srcFacePti == -1 ? bb.centre() : srcPts[srcFace[srcFacePti]];
301 
302  const pointIndexHit sample =
303  treePtr_->findNearest(srcPt, 0.25*bb.magSqr(), fnOp);
304 
305  if (sample.hit() && isCandidate(srcFacei, sample.index()))
306  {
307  targetFacei = sample.index();
308 
309  if (debug)
310  {
311  Pout<< "Source point = " << srcPt << ", Sample point = "
312  << sample.point() << ", Sample index = " << sample.index()
313  << endl;
314  }
315  }
316 
317  return targetFacei;
318 }
319 
320 
322 (
323  const label facei,
324  const primitivePatch& patch,
325  const DynamicList<label>& visitedFaces,
326  DynamicList<label>& faceIDs
327 ) const
328 {
329  static const scalar thetaCos = Foam::cos(degToRad(89.0));
330 
331  const labelList& nbrFaces = patch.faceFaces()[facei];
332 
333  // Filter out faces already visited from face neighbours
334  for (const label nbrFacei : nbrFaces)
335  {
336  // Prevent addition of face if it is not on the same plane-ish
337  if (!visitedFaces.found(nbrFacei) && !faceIDs.found(nbrFacei))
338  {
339  const vector& n1 = patch.faceNormals()[facei];
340  const vector& n2 = patch.faceNormals()[nbrFacei];
341 
342  const scalar cosI = n1 & n2;
343 
344  if (cosI > thetaCos)
345  {
346  faceIDs.append(nbrFacei);
347  }
348  }
349  }
350 }
351 
352 
354 (
355  const primitivePatch& patch,
356  List<DynamicList<face>>& tris,
357  List<scalar>& magSf
358 ) const
359 {
360  const pointField& points = patch.points();
361  tris.setSize(patch.size());
362  magSf.setSize(patch.size());
363 
364  const auto& faceNormals = patch.faceNormals();
365 
366  // Using methods that index into existing points
367  forAll(patch, facei)
368  {
369  tris[facei].clear();
370 
371  switch (triMode_)
372  {
374  {
375  faceAreaIntersect::triangleFan(patch[facei], tris[facei]);
376  break;
377  }
379  {
380  patch[facei].triangles(points, tris[facei]);
381  break;
382  }
383  }
384 
385  const DynamicList<face>& triFaces = tris[facei];
386  magSf[facei] = 0;
387  for (const face& f : triFaces)
388  {
389  magSf[facei] +=
391  (
392  points[f[0]],
393  points[f[1]],
394  points[f[2]]
395  ).areaNormal()
396  & faceNormals[facei];
397  }
398  }
399 }
400 
401 
403 {
404  if (!requireMatch_ && distributed())
405  {
406  scalarList newTgtMagSf(std::move(tgtMagSf_));
407 
408  // Assign default sizes. Override selected values with calculated
409  // values. This is to support ACMI where some of the target faces
410  // are never used (so never get sent over and hence never assigned
411  // to)
412  tgtMagSf_ = tgtPatch0().magFaceAreas();
413 
414  for (const labelList& smap : this->extendedTgtMapPtr_->subMap())
415  {
416  UIndirectList<scalar>(tgtMagSf_, smap) =
417  UIndirectList<scalar>(newTgtMagSf, smap);
418  }
419  }
420 }
421 
422 
423 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
424 
426 (
427  const dictionary& dict,
428  const bool reverseTarget
429 )
430 :
431  AMIInterpolation(dict, reverseTarget),
432  maxDistance2_(dict.getOrDefault<scalar>("maxDistance2", -1)),
433  minCosAngle_(dict.getOrDefault<scalar>("minCosAngle", -1)),
434  srcTris_(),
435  tgtTris_(),
436  extendedTgtPatchPtr_(nullptr),
437  extendedTgtFaces_(),
438  extendedTgtPoints_(),
439  extendedTgtFaceIDs_(),
440  extendedTgtMapPtr_(nullptr),
441  srcNonOverlap_(),
442  triMode_
443  (
444  faceAreaIntersect::triangulationModeNames_.getOrDefault
445  (
446  "triMode",
447  dict,
449  )
450  )
451 {}
452 
453 
455 (
456  const bool requireMatch,
457  const bool reverseTarget,
458  const scalar lowWeightCorrection,
460 )
461 :
462  AMIInterpolation(requireMatch, reverseTarget, lowWeightCorrection),
463  maxDistance2_(-1),
464  minCosAngle_(-1),
465  srcTris_(),
466  tgtTris_(),
467  extendedTgtPatchPtr_(nullptr),
468  extendedTgtFaces_(),
469  extendedTgtPoints_(),
470  extendedTgtFaceIDs_(),
471  extendedTgtMapPtr_(nullptr),
472  srcNonOverlap_(),
473  triMode_(triMode)
474 {}
475 
476 
478 :
479  AMIInterpolation(ami),
480  maxDistance2_(ami.maxDistance2_),
481  minCosAngle_(ami.minCosAngle_),
482  srcTris_(),
483  tgtTris_(),
484  extendedTgtPatchPtr_(nullptr),
485  extendedTgtFaces_(),
486  extendedTgtPoints_(),
487  extendedTgtFaceIDs_(),
488  extendedTgtMapPtr_(nullptr),
489  srcNonOverlap_(),
490  triMode_(ami.triMode_)
491 {}
492 
493 
494 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
495 
497 (
498  const primitivePatch& srcPatch,
499  const primitivePatch& tgtPatch,
500  const autoPtr<searchableSurface>& surfPtr
501 )
502 {
503  if (AMIInterpolation::calculate(srcPatch, tgtPatch, surfPtr))
504  {
505  // Create a representation of the target patch that covers the source
506  // patch
507  if (distributed())
508  {
509  createExtendedTgtPatch();
510  }
511 
512  const auto& src = this->srcPatch();
513  const auto& tgt = this->tgtPatch();
514 
515 
516  if (maxDistance2_ > 0)
517  {
518  // Early trigger face centre calculation
519  (void)src.faceCentres();
520  (void)tgt.faceCentres();
521  // Early trigger face normals calculation
522  (void)src.faceNormals();
523  (void)tgt.faceNormals();
524  }
525  if (minCosAngle_ > -1)
526  {
527  // Early trigger face normals calculation
528  (void)src.faceNormals();
529  (void)tgt.faceNormals();
530  }
531 
532 
533  // Initialise area magnitudes
534  srcMagSf_.setSize(src.size(), 1.0);
535  tgtMagSf_.setSize(tgt.size(), 1.0);
536 
537  // Source and target patch triangulations
538  triangulatePatch(src, srcTris_, srcMagSf_);
539  triangulatePatch(tgt, tgtTris_, tgtMagSf_);
540 
541  checkPatches();
542 
543  // Set initial sizes for weights and addressing - must be done even if
544  // returns false below
545  srcAddress_.setSize(src.size());
546  srcWeights_.setSize(src.size());
547  tgtAddress_.setSize(tgt.size());
548  tgtWeights_.setSize(tgt.size());
549 
550  return true;
551  }
552 
553  return false;
554 }
555 
556 
557 void Foam::advancingFrontAMI::write(Ostream& os) const
558 {
560  os.writeEntryIfDifferent<scalar>("maxDistance2", -1, maxDistance2_);
561  os.writeEntryIfDifferent<scalar>("minCosAngle", -1, minCosAngle_);
563  (
564  "triMode",
567  );
568 }
569 
570 
571 // ************************************************************************* //
List< scalar > scalarList
List of scalar.
Definition: scalarList.H:32
label findTargetFace(const label srcFacei, const UList< label > &excludeFaces=UList< label >::null(), const label srcFacePti=-1) const
dictionary dict
void writeIntersectionOBJ(const scalar area, const face &f1, const face &f2, const pointField &f1Points, const pointField &f2Points) const
Write triangle intersection to OBJ file.
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
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:578
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Unit conversion functions.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
bool initialiseWalk(label &srcFacei, label &tgtFacei)
Initialise walk and return true if all ok.
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
This class describes the interaction of an object (often a face) and a point. It carries the info of ...
Definition: pointIndexHit.H:44
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:196
void createExtendedTgtPatch()
Create a map that extends tgtPatch so that it covers srcPatch.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:414
virtual void write(Ostream &os) const
Write.
const primitivePatch & tgtPatch() const
Return const access to the target patch.
const Field< point_type > & faceNormals() const
Return face unit normals for patch.
scalar magSqr() const
The magnitude/length squared of bounding box diagonal.
Definition: boundBoxI.H:204
A list of faces which address into the list of points.
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of &#39;true&#39; entries.
Definition: BitOps.H:73
const point_type & point() const noexcept
Return point, no checks.
A List obtained as a section of another List.
Definition: SubList.H:50
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
virtual bool calculate(const primitivePatch &srcPatch, const primitivePatch &tgtPatch, const autoPtr< searchableSurface > &surfPtr=nullptr)
Update addressing, weights and (optional) centroids.
dimensionedScalar cos(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
const pointField & points
static autoPtr< indexedOctree< treeDataPoint > > createTree(const pointField &points)
Construct search tree for points.
const wordList area
Standard area field types (scalar, vector, tensor, etc)
const Field< point_type > & faceCentres() const
Return face centres for patch.
Vector< scalar > vector
Definition: vector.H:57
label index() const noexcept
Return the hit index.
errorManip< error > abort(error &err)
Definition: errorManip.H:139
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:327
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
bool isCandidate(const label srcFacei, const label tgtFacei) const
Is source/target a valid pair (i.e. not too far/different.
int debug
Static debugging option.
Face intersection class.
OBJstream os(runTime.globalPath()/outputName)
defineTypeNameAndDebug(combustionModel, 0)
static void triangleFan(const face &f, DynamicList< face > &faces)
Decompose face into triangle fan.
labelList f(nPoints)
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field...
void checkPatches() const
Check AMI patch coupling.
void appendNbrFaces(const label facei, const primitivePatch &patch, const DynamicList< label > &visitedFaces, DynamicList< label > &faceIDs) const
Add faces neighbouring facei to the ID list.
virtual void nonConformalCorrection()
Correction for non-conformal interpolations, e.g. for ACMI.
virtual bool calculate(const primitivePatch &srcPatch, const primitivePatch &tgtPatch, const autoPtr< searchableSurface > &surfPtr=nullptr)
Update addressing, weights and (optional) centroids.
bool hit() const noexcept
Is there a hit?
bool requireMatch_
Flag to indicate that the two patches must be matched/an overlap exists between them.
vector point
Point is a vector.
Definition: point.H:37
#define WarningInFunction
Report a warning using Foam::Warning.
void triangulatePatch(const primitivePatch &patch, List< DynamicList< face >> &tris, List< scalar > &magSf) const
Helper function to decompose a patch.
Base class for Arbitrary Mesh Interface (AMI) methods.
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
const std::string patch
OpenFOAM patch number as a std::string.
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition: triangleFwd.H:39
label n
static const Enum< triangulationMode > triangulationModeNames_
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
List< label > labelList
A List of labels.
Definition: List.H:62
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition: autoPtr.H:178
point centre() const
The centre (midpoint) of the bounding box.
Definition: boundBoxI.H:186
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
const primitivePatch & srcPatch() const
Return const access to the source patch.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Namespace for OpenFOAM.
virtual void write(Ostream &os) const
Write.
advancingFrontAMI(const dictionary &dict, const bool reverseTarget)
Construct from components.