faceAreaWeightAMI2D.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) 2020 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "faceAreaWeightAMI2D.H"
29 #include "profiling.H"
30 #include "OBJstream.H"
32 #include "triangle2D.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(faceAreaWeightAMI2D, 0);
39  addToRunTimeSelectionTable(AMIInterpolation, faceAreaWeightAMI2D, dict);
41  (
42  AMIInterpolation,
43  faceAreaWeightAMI2D,
44  component
45  );
46 }
47 
48 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
49 
51 (
52  const label srcFacei,
53  const labelList& tgtFaceCandidates,
54  const boundBox& srcFaceBb
55 ) const
56 {
57  Info<< "NO MATCH for source face " << srcFacei << endl;
58  {
59  const auto& src = this->srcPatch();
60  const auto& tgt = this->tgtPatch(); // might be the extended patch!
61 
62  OFstream os("no_match_" + Foam::name(srcFacei) + ".obj");
63 
64  const pointField& srcPoints = src.points();
65  const pointField& tgtPoints = tgt.points();
66 
67  label np = 0;
68 
69  // Write source face
70  const face& srcF = src[srcFacei];
71  string faceStr = "f";
72  for (const label pointi : srcF)
73  {
74  const point& p = srcPoints[pointi];
75  os << "v " << p.x() << " " << p.y() << " " << p.z() << nl;
76  ++np;
77  faceStr += " " + Foam::name(np);
78  }
79  os << faceStr.c_str() << " " << (np - srcF.size() + 1) << nl;
80 
81  // Write target faces as lines
82  for (const label tgtFacei : tgtFaceCandidates)
83  {
84  const face& tgtF = tgt[tgtFacei];
85  forAll(tgtF, pointi)
86  {
87  const point& p = tgtPoints[tgtF[pointi]];
88  os << "v " << p.x() << " " << p.y() << " " << p.z() << nl;
89  ++np;
90  if (pointi)
91  {
92  os << "l " << np-1 << " " << np << nl;
93  }
94  }
95  os << "l " << (np - tgtF.size() + 1) << " " << np << nl;
96  }
97  }
98 
99  {
100  OFstream os("no_match_" + Foam::name(srcFacei) + "_bb.obj");
101 
102  const pointField points(srcFaceBb.points());
103  for (const point& p : points)
104  {
105  os << "v " << p.x() << " " << p.y() << " " << p.z() << endl;
106  }
107  os << "l 1 2" << nl;
108  os << "l 2 4" << nl;
109  os << "l 4 3" << nl;
110  os << "l 3 1" << nl;
111  os << "l 5 6" << nl;
112  os << "l 6 8" << nl;
113  os << "l 8 7" << nl;
114  os << "l 7 5" << nl;
115  os << "l 5 1" << nl;
116  os << "l 6 2" << nl;
117  os << "l 8 4" << nl;
118  os << "l 7 3" << nl;
119  }
120 }
121 
122 
124 (
125  const label srcFacei,
126  const label tgtFacei,
127  DynamicList<label>& srcAddr,
128  DynamicList<scalar>& srcWght,
129  DynamicList<vector>& srcCtr,
130  DynamicList<label>& tgtAddr,
131  DynamicList<scalar>& tgtWght
132 ) const
133 {
134  addProfiling(ami, "faceAreaWeightAMI2D::calcInterArea");
135 
136  // Quick reject if either face has zero area
137  if
138  (
139  (srcMagSf_[srcFacei] < ROOTVSMALL)
140  || (tgtMagSf_[tgtFacei] < ROOTVSMALL)
141  )
142  {
143  return;
144  }
145 
146  const auto& srcPatch = this->srcPatch();
147  const auto& tgtPatch = this->tgtPatch();
148 
149  const pointField& srcPoints = srcPatch.points();
150  const pointField& tgtPoints = tgtPatch.points();
151 
152  const auto& srcTris = srcTris_[srcFacei];
153  const auto& tgtTris = tgtTris_[tgtFacei];
154 
155  const auto& srcFaceNormals = srcPatch.faceNormals();
156 
157  //triangle2D::debug = 1;
158 
159  scalar area = 0;
160  vector centroid = Zero;
161 
162  for (const auto& tris : srcTris)
163  {
164  const vector& origin = srcPoints[tris[0]];
165  const vector p10(srcPoints[tris[1]] - origin);
166  const vector p20(srcPoints[tris[2]] - origin);
167  const vector axis1(p10/(mag(p10) + ROOTVSMALL));
168  const vector axis2(srcFaceNormals[srcFacei]^axis1);
169 
170  triangle2D s
171  (
172  vector2D(0, 0),
173  vector2D(axis1&p10, axis2&p10),
174  vector2D(axis1&p20, axis2&p20)
175  );
176 
177  for (const auto& trit : tgtTris)
178  {
179  // Triangle t has opposite orientation wrt triangle s
180  triangle2D t
181  (
182  tgtPoints[trit[0]] - origin,
183  tgtPoints[trit[2]] - origin,
184  tgtPoints[trit[1]] - origin,
185  axis1,
186  axis2
187  );
188 
189  scalar da = 0;
190  vector2D c(Zero);
191 
192  if (t.snapClosePoints(s) == 3)
193  {
194  c = s.centre();
195  da = mag(s.area());
196  }
197  else
198  {
199  s.interArea(t, c, da);
200  }
201 
202  area += da;
203  centroid += da*(origin + c.x()*axis1 + c.y()*axis2);
204  }
205  }
206 
207  //triangle2D::debug = 0;
208 
209  if (area/srcMagSf_[srcFacei] > triangle2D::relTol)
210  {
211  centroid /= area + ROOTVSMALL;
212 
213  srcAddr.append(tgtFacei);
214  srcWght.append(area);
215  srcCtr.append(centroid);
216 
217  tgtAddr.append(srcFacei);
218  tgtWght.append(area);
219  }
220 }
221 
222 
224 (
225  const AABBTree<face>& tree,
226  const List<boundBox>& tgtFaceBbs,
227  const boundBox& srcFaceBb
228 ) const
229 {
230  labelHashSet faceIds(16);
231 
232  const auto& treeBb = tree.boundBoxes();
233  const auto& treeAddr = tree.addressing();
234 
235  forAll(treeBb, boxi)
236  {
237  const auto& tbb = treeBb[boxi];
238 
239  if (srcFaceBb.overlaps(tbb))
240  {
241  const auto& boxAddr = treeAddr[boxi];
242 
243  for (const auto& tgtFacei : boxAddr)
244  {
245  if (srcFaceBb.overlaps(tgtFaceBbs[tgtFacei]))
246  {
247  faceIds.insert(tgtFacei);
248  }
249  }
250  }
251  }
252 
253  return faceIds.toc();
254 }
255 
256 
257 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
258 
260 (
261  const dictionary& dict,
262  const bool reverseTarget
263 )
264 :
265  advancingFrontAMI(dict, reverseTarget),
266  Cbb_(dict.getCheckOrDefault<scalar>("Cbb", 0.1, scalarMinMax::ge(SMALL)))
267 {}
268 
269 
271 (
272  const bool requireMatch,
273  const bool reverseTarget,
274  const scalar lowWeightCorrection,
276  const bool restartUncoveredSourceFace
277 )
278 :
280  (
281  requireMatch,
282  reverseTarget,
283  lowWeightCorrection,
284  triMode
285  ),
286  Cbb_(0.1)
287 {}
288 
289 
291 (
292  const faceAreaWeightAMI2D& ami
293 )
294 :
295  advancingFrontAMI(ami),
296  Cbb_(ami.Cbb_)
297 {}
298 
299 
300 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
301 
303 (
304  const primitivePatch& srcPatch,
305  const primitivePatch& tgtPatch,
306  const autoPtr<searchableSurface>& surfPtr
307 )
308 {
309  if (upToDate_)
310  {
311  return false;
312  }
313 
314  addProfiling(ami, "faceAreaWeightAMI2D::calculate");
315 
316  advancingFrontAMI::calculate(srcPatch, tgtPatch, surfPtr);
317 
318  const auto& src = this->srcPatch();
319  const auto& tgt = this->tgtPatch(); // might be the extended patch!
320 
321  bool validSize = true;
322 
323  // Check that patch sizes are valid
324  if (!src.size())
325  {
326  validSize = false;
327  }
328  else if (!tgt.size())
329  {
331  << src.size() << " source faces but no target faces" << endl;
332 
333  validSize = false;
334  }
335 
336  srcCentroids_.setSize(srcAddress_.size());
337 
338  // Temporary storage for addressing and weights
339  List<DynamicList<label>> tgtAddr(tgt.size());
340  List<DynamicList<scalar>> tgtWght(tgt.size());
341 
342  // Find an approximate face length scale
343  const scalar lf(Cbb_*Foam::sqrt(gAverage(srcMagSf_)));
344 
345  // Expansion to apply to source face bounding box
346  const vector d(lf*vector::one);
347 
348  if (validSize)
349  {
350  // Create the tgt tree
351  const bool equalBinSize = true;
352  const label maxLevel = 10;
353  const label minBinSize = 4;
354  AABBTree<face> tree
355  (
356  tgt,
357  tgt.points(),
358  equalBinSize,
359  maxLevel,
360  minBinSize
361  );
362 
363  const auto& tgtPoints = tgt.points();
364  List<boundBox> tgtFaceBbs(tgt.size());
365  forAll(tgt, facei)
366  {
367  tgtFaceBbs[facei] = boundBox(tgtPoints, tgt[facei], false);
368  }
369 
370  DynamicList<label> nonOverlapFaces;
371 
372  const auto& srcPoints = src.points();
373 
374  forAll(src, srcFacei)
375  {
376  const face& srcFace = src[srcFacei];
377 
378  treeBoundBox srcFaceBb(srcPoints, srcFace);
379  srcFaceBb.min() -= d;
380  srcFaceBb.max() += d;
381 
382  const labelList tgtFaces
383  (
384  overlappingTgtFaces(tree, tgtFaceBbs, srcFaceBb)
385  );
386 
387  DynamicList<label> srcAddr(tgtFaces.size());
388  DynamicList<scalar> srcWght(tgtFaces.size());
389  DynamicList<point> srcCtr(tgtFaces.size());
390 
391  for (const label tgtFacei : tgtFaces)
392  {
393  storeInterArea
394  (
395  srcFacei,
396  tgtFacei,
397  srcAddr,
398  srcWght,
399  srcCtr,
400  tgtAddr[tgtFacei],
401  tgtWght[tgtFacei]
402  );
403  }
404 
405  if (mustMatchFaces() && srcAddr.empty())
406  {
407  if (debug) writeNoMatch(srcFacei, tgtFaces, srcFaceBb);
408 
409  // FatalErrorInFunction
410  // << "Unable to find match for source face " << srcFace
411  // << exit(FatalError);
412  }
413 
414  srcAddress_[srcFacei].transfer(srcAddr);
415  srcWeights_[srcFacei].transfer(srcWght);
416  srcCentroids_[srcFacei].transfer(srcCtr);
417  }
418 
419  srcNonOverlap_.transfer(nonOverlapFaces);
420 
421  if (debug && !srcNonOverlap_.empty())
422  {
423  Pout<< " AMI: " << srcNonOverlap_.size()
424  << " non-overlap faces identified"
425  << endl;
426  }
427  }
428 
429  // Transfer data to persistent storage
430 
431  forAll(tgtAddr, i)
432  {
433  tgtAddress_[i].transfer(tgtAddr[i]);
434  tgtWeights_[i].transfer(tgtWght[i]);
435  }
436 
437  if (distributed())
438  {
439  const primitivePatch& srcPatch0 = this->srcPatch0();
440  const primitivePatch& tgtPatch0 = this->tgtPatch0();
441 
442  // Create global indexing for each original patch
443  globalIndex globalSrcFaces(srcPatch0.size());
444  globalIndex globalTgtFaces(tgtPatch0.size());
445 
446  for (labelList& addressing : srcAddress_)
447  {
448  for (label& addr : addressing)
449  {
450  addr = extendedTgtFaceIDs_[addr];
451  }
452  }
453 
454  for (labelList& addressing : tgtAddress_)
455  {
456  globalSrcFaces.inplaceToGlobal(addressing);
457  }
458 
459  // Send data back to originating procs. Note that contributions
460  // from different processors get added (ListOps::appendEqOp)
461 
462  mapDistributeBase::distribute
463  (
464  Pstream::commsTypes::nonBlocking,
465  List<labelPair>(),
466  tgtPatch0.size(),
467  extendedTgtMapPtr_->constructMap(),
468  false, // has flip
469  extendedTgtMapPtr_->subMap(),
470  false, // has flip
471  tgtAddress_,
472  labelList(),
474  flipOp()
475  );
476 
477  mapDistributeBase::distribute
478  (
479  Pstream::commsTypes::nonBlocking,
480  List<labelPair>(),
481  tgtPatch0.size(),
482  extendedTgtMapPtr_->constructMap(),
483  false,
484  extendedTgtMapPtr_->subMap(),
485  false,
486  tgtWeights_,
487  scalarList(),
489  flipOp()
490  );
491 
492  // Note: using patch face areas calculated by the AMI method
493  extendedTgtMapPtr_->reverseDistribute(tgtPatch0.size(), tgtMagSf_);
494 
495  // Cache maps and reset addresses
496  List<Map<label>> cMapSrc;
497  srcMapPtr_.reset
498  (
499  new mapDistribute(globalSrcFaces, tgtAddress_, cMapSrc)
500  );
501 
502  List<Map<label>> cMapTgt;
503  tgtMapPtr_.reset
504  (
505  new mapDistribute(globalTgtFaces, srcAddress_, cMapTgt)
506  );
507  }
508 
509  // Convert the weights from areas to normalised values
510  normaliseWeights(requireMatch_, true);
511 
512  nonConformalCorrection();
513 
514  upToDate_ = true;
515 
516  return upToDate_;
517 }
518 
519 
521 {
523 }
524 
525 
526 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::PrimitivePatch::points
const Field< point_type > & points() const noexcept
Return reference to global points.
Definition: PrimitivePatch.H:299
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
profiling.H
Foam::scalarList
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:64
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::component
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
Definition: FieldFieldFunctions.C:44
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::faceAreaIntersect::triangulationMode
triangulationMode
Definition: faceAreaIntersect.H:63
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::DynamicList< label >
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:604
Foam::treeBoundBox
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:86
Foam::AABBTree::boundBoxes
const List< treeBoundBox > & boundBoxes() const
Return the bounding boxes making up the tree.
Definition: AABBTree.C:442
Foam::faceAreaWeightAMI2D::storeInterArea
void storeInterArea(const label srcFacei, const label tgtFacei, DynamicList< label > &srcAddr, DynamicList< scalar > &srcWght, DynamicList< vector > &srcCtr, DynamicList< label > &tgtAddr, DynamicList< scalar > &tgtWght) const
Definition: faceAreaWeightAMI2D.C:124
Foam::faceAreaWeightAMI2D
Face area weighted Arbitrary Mesh Interface (AMI) method that performs the intersection of src and tg...
Definition: faceAreaWeightAMI2D.H:52
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::globalIndex::inplaceToGlobal
void inplaceToGlobal(labelList &labels) const
From local to global index (inplace)
Definition: globalIndexI.H:283
Foam::boundBox::max
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:97
Foam::Vector2D< scalar >
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::HashSet< label, Hash< label > >
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::boundBox::overlaps
bool overlaps(const boundBox &bb) const
Overlaps/touches boundingBox?
Definition: boundBoxI.H:221
Foam::Field< vector >
Foam::faceAreaWeightAMI2D::faceAreaWeightAMI2D
faceAreaWeightAMI2D(const dictionary &dict, const bool reverseTarget=false)
Construct from dictionary.
Definition: faceAreaWeightAMI2D.C:260
Foam::faceAreaWeightAMI2D::writeNoMatch
void writeNoMatch(const label srcFacei, const labelList &tgtFaceCandidates, const boundBox &srcFaceBb) const
Definition: faceAreaWeightAMI2D.C:51
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:511
Foam::mapDistribute
Class containing processor-to-processor mapping information.
Definition: mapDistribute.H:163
Foam::AMIInterpolation::write
virtual void write(Ostream &os) const
Write.
Definition: AMIInterpolation.C:1243
Foam::faceAreaWeightAMI2D::calculate
virtual bool calculate(const primitivePatch &srcPatch, const primitivePatch &tgtPatch, const autoPtr< searchableSurface > &surfPtr=nullptr)
Update addressing, weights and (optional) centroids.
Definition: faceAreaWeightAMI2D.C:303
Foam::triangle2D
2-D triangle and queries
Definition: triangle2D.H:52
dict
dictionary dict
Definition: searchingEngine.H:14
addProfiling
#define addProfiling(name, descr)
Define profiling trigger with specified name and description string.
Definition: profilingTrigger.H:113
faceAreaWeightAMI2D.H
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:53
Foam::faceAreaWeightAMI2D::overlappingTgtFaces
labelList overlappingTgtFaces(const AABBTree< face > &tree, const List< boundBox > &tgtFaceBbs, const boundBox &srcFaceBb) const
Return the set of tgt face IDs that overlap the src face bb.
Definition: faceAreaWeightAMI2D.C:224
Foam::ListOps::appendEqOp
List helper to append y elements onto the end of x.
Definition: ListOps.H:570
Foam::faceAreaWeightAMI2D::Cbb_
scalar Cbb_
Face bounding box factor.
Definition: faceAreaWeightAMI2D.H:61
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::faceAreaWeightAMI2D::write
virtual void write(Ostream &os) const
Write.
Definition: faceAreaWeightAMI2D.C:520
Foam::AABBTree::addressing
const List< labelList > & addressing() const
Return the contents addressing.
Definition: AABBTree.C:449
Foam::vector2D
Vector2D< scalar > vector2D
A 2D vector of scalars obtained from the generic Vector2D.
Definition: vector2D.H:51
Foam::Vector< scalar >
Foam::AABBTree
Templated tree of axis-aligned bounding boxes (AABB)
Definition: AABBTree.H:59
Foam::List< label >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
Foam::advancingFrontAMI
Base class for Arbitrary Mesh Interface (AMI) methods.
Definition: advancingFrontAMI.H:55
Foam::PrimitivePatch::faceNormals
const Field< point_type > & faceNormals() const
Return face unit normals for patch.
Definition: PrimitivePatch.C:445
Foam::triangle2D::snapClosePoints
label snapClosePoints(const triangle2D &triB)
Definition: triangle2D.C:93
Foam::fieldTypes::area
const wordList area
Standard area field types (scalar, vector, tensor, etc)
Foam::boundBox
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
triangle2D.H
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::boundBox::points
tmp< pointField > points() const
Corner points in an order corresponding to a 'hex' cell.
Definition: boundBox.C:118
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
OBJstream.H
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:79
Foam::flipOp
Functor to negate primitives. Dummy for most other types.
Definition: flipOp.H:53