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,2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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
36namespace Foam
37{
41 (
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/too far away/wrong orientation
137 if (!isCandidate(srcFacei, tgtFacei))
138 {
139 return;
140 }
141
142 const auto& srcPatch = this->srcPatch();
143 const auto& tgtPatch = this->tgtPatch();
144
145 const pointField& srcPoints = srcPatch.points();
146 const pointField& tgtPoints = tgtPatch.points();
147
148 const auto& srcTris = srcTris_[srcFacei];
149 const auto& tgtTris = tgtTris_[tgtFacei];
150
151 const auto& srcFaceNormals = srcPatch.faceNormals();
152
153 //triangle2D::debug = 1;
154
155 scalar area = 0;
156 vector centroid = Zero;
157
158 for (const auto& tris : srcTris)
159 {
160 const vector& origin = srcPoints[tris[0]];
161 const vector p10(srcPoints[tris[1]] - origin);
162 const vector p20(srcPoints[tris[2]] - origin);
163 const vector axis1(p10/(mag(p10) + ROOTVSMALL));
164 const vector axis2(srcFaceNormals[srcFacei]^axis1);
165
167 (
168 vector2D(0, 0),
169 vector2D(axis1&p10, axis2&p10),
170 vector2D(axis1&p20, axis2&p20)
171 );
172
173 for (const auto& trit : tgtTris)
174 {
175 // Triangle t has opposite orientation wrt triangle s
176 triangle2D t
177 (
178 tgtPoints[trit[0]] - origin,
179 tgtPoints[trit[2]] - origin,
180 tgtPoints[trit[1]] - origin,
181 axis1,
182 axis2
183 );
184
185 scalar da = 0;
186 vector2D c(Zero);
187
188 if (t.snapClosePoints(s) == 3)
189 {
190 c = s.centre();
191 da = mag(s.area());
192 }
193 else
194 {
195 s.interArea(t, c, da);
196 }
197
198 area += da;
199 centroid += da*(origin + c.x()*axis1 + c.y()*axis2);
200 }
201 }
202
203 //triangle2D::debug = 0;
204
205 if (area/srcMagSf_[srcFacei] > triangle2D::relTol)
206 {
207 centroid /= area + ROOTVSMALL;
208
209 srcAddr.append(tgtFacei);
210 srcWght.append(area);
211 srcCtr.append(centroid);
212
213 tgtAddr.append(srcFacei);
214 tgtWght.append(area);
215 }
216}
217
218
220(
221 const AABBTree<face>& tree,
222 const List<boundBox>& tgtFaceBbs,
223 const boundBox& srcFaceBb
224) const
225{
226 labelHashSet faceIds(16);
227
228 const auto& treeBb = tree.boundBoxes();
229 const auto& treeAddr = tree.addressing();
230
231 forAll(treeBb, boxi)
232 {
233 const auto& tbb = treeBb[boxi];
234
235 if (srcFaceBb.overlaps(tbb))
236 {
237 const auto& boxAddr = treeAddr[boxi];
238
239 for (const auto& tgtFacei : boxAddr)
240 {
241 if (srcFaceBb.overlaps(tgtFaceBbs[tgtFacei]))
242 {
243 faceIds.insert(tgtFacei);
244 }
245 }
246 }
247 }
248
249 return faceIds.toc();
250}
251
252
253// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
254
256(
257 const dictionary& dict,
258 const bool reverseTarget
259)
260:
261 advancingFrontAMI(dict, reverseTarget),
262 Cbb_(dict.getCheckOrDefault<scalar>("Cbb", 0.1, scalarMinMax::ge(SMALL)))
263{}
264
265
267(
268 const bool requireMatch,
269 const bool reverseTarget,
270 const scalar lowWeightCorrection,
272 const bool restartUncoveredSourceFace
273)
274:
276 (
277 requireMatch,
278 reverseTarget,
279 lowWeightCorrection,
280 triMode
281 ),
282 Cbb_(0.1)
283{}
284
285
287(
288 const faceAreaWeightAMI2D& ami
289)
290:
292 Cbb_(ami.Cbb_)
293{}
294
295
296// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
297
299(
300 const primitivePatch& srcPatch,
301 const primitivePatch& tgtPatch,
302 const autoPtr<searchableSurface>& surfPtr
303)
304{
305 if (upToDate_)
306 {
307 return false;
308 }
309
310 addProfiling(ami, "faceAreaWeightAMI2D::calculate");
311
312 advancingFrontAMI::calculate(srcPatch, tgtPatch, surfPtr);
313
314 const auto& src = this->srcPatch();
315 const auto& tgt = this->tgtPatch(); // might be the extended patch!
316
317 bool validSize = true;
318
319 // Check that patch sizes are valid
320 if (!src.size())
321 {
322 validSize = false;
323 }
324 else if (!tgt.size())
325 {
327 << src.size() << " source faces but no target faces" << endl;
328
329 validSize = false;
330 }
331
332 srcCentroids_.setSize(srcAddress_.size());
333
334 // Temporary storage for addressing and weights
335 List<DynamicList<label>> tgtAddr(tgt.size());
336 List<DynamicList<scalar>> tgtWght(tgt.size());
337
338 // Find an approximate face length scale
339 const scalar lf(Cbb_*Foam::sqrt(gAverage(srcMagSf_)));
340
341 // Expansion to apply to source face bounding box
342 const vector d(lf*vector::one);
343
344 if (validSize)
345 {
346 // Create the tgt tree
347 const bool equalBinSize = true;
348 const label maxLevel = 10;
349 const label minBinSize = 4;
350 AABBTree<face> tree
351 (
352 tgt,
353 tgt.points(),
354 equalBinSize,
355 maxLevel,
356 minBinSize
357 );
358
359 const auto& tgtPoints = tgt.points();
360 List<boundBox> tgtFaceBbs(tgt.size());
361 forAll(tgt, facei)
362 {
363 tgtFaceBbs[facei] = boundBox(tgtPoints, tgt[facei], false);
364 }
365
366 DynamicList<label> nonOverlapFaces;
367
368 const auto& srcPoints = src.points();
369
370 forAll(src, srcFacei)
371 {
372 const face& srcFace = src[srcFacei];
373
374 treeBoundBox srcFaceBb(srcPoints, srcFace);
375 srcFaceBb.min() -= d;
376 srcFaceBb.max() += d;
377
378 const labelList tgtFaces
379 (
380 overlappingTgtFaces(tree, tgtFaceBbs, srcFaceBb)
381 );
382
383 DynamicList<label> srcAddr(tgtFaces.size());
384 DynamicList<scalar> srcWght(tgtFaces.size());
385 DynamicList<point> srcCtr(tgtFaces.size());
386
387 for (const label tgtFacei : tgtFaces)
388 {
389 storeInterArea
390 (
391 srcFacei,
392 tgtFacei,
393 srcAddr,
394 srcWght,
395 srcCtr,
396 tgtAddr[tgtFacei],
397 tgtWght[tgtFacei]
398 );
399 }
400
401 if (mustMatchFaces() && srcAddr.empty())
402 {
403 if (debug) writeNoMatch(srcFacei, tgtFaces, srcFaceBb);
404
405 // FatalErrorInFunction
406 // << "Unable to find match for source face " << srcFace
407 // << exit(FatalError);
408 }
409
410 srcAddress_[srcFacei].transfer(srcAddr);
411 srcWeights_[srcFacei].transfer(srcWght);
412 srcCentroids_[srcFacei].transfer(srcCtr);
413 }
414
415 srcNonOverlap_.transfer(nonOverlapFaces);
416
417 if (debug && !srcNonOverlap_.empty())
418 {
419 Pout<< " AMI: " << srcNonOverlap_.size()
420 << " non-overlap faces identified"
421 << endl;
422 }
423 }
424
425 // Transfer data to persistent storage
426
427 forAll(tgtAddr, i)
428 {
429 tgtAddress_[i].transfer(tgtAddr[i]);
430 tgtWeights_[i].transfer(tgtWght[i]);
431 }
432
433 if (distributed())
434 {
435 const primitivePatch& srcPatch0 = this->srcPatch0();
436 const primitivePatch& tgtPatch0 = this->tgtPatch0();
437
438 // Create global indexing for each original patch
439 globalIndex globalSrcFaces(srcPatch0.size());
440 globalIndex globalTgtFaces(tgtPatch0.size());
441
442 for (labelList& addressing : srcAddress_)
443 {
444 for (label& addr : addressing)
445 {
446 addr = extendedTgtFaceIDs_[addr];
447 }
448 }
449
450 for (labelList& addressing : tgtAddress_)
451 {
452 globalSrcFaces.inplaceToGlobal(addressing);
453 }
454
455 // Send data back to originating procs. Note that contributions
456 // from different processors get added (ListOps::appendEqOp)
457
458 mapDistributeBase::distribute
459 (
462 tgtPatch0.size(),
463 extendedTgtMapPtr_->constructMap(),
464 false, // has flip
465 extendedTgtMapPtr_->subMap(),
466 false, // has flip
467 tgtAddress_,
468 labelList(),
470 flipOp()
471 );
472
473 mapDistributeBase::distribute
474 (
477 tgtPatch0.size(),
478 extendedTgtMapPtr_->constructMap(),
479 false,
480 extendedTgtMapPtr_->subMap(),
481 false,
482 tgtWeights_,
483 scalarList(),
485 flipOp()
486 );
487
488 // Note: using patch face areas calculated by the AMI method
489 extendedTgtMapPtr_->reverseDistribute(tgtPatch0.size(), tgtMagSf_);
490
491 // Cache maps and reset addresses
492 List<Map<label>> cMapSrc;
493 srcMapPtr_.reset
494 (
495 new mapDistribute(globalSrcFaces, tgtAddress_, cMapSrc)
496 );
497
498 List<Map<label>> cMapTgt;
499 tgtMapPtr_.reset
500 (
501 new mapDistribute(globalTgtFaces, srcAddress_, cMapTgt)
502 );
503 }
504
505 // Convert the weights from areas to normalised values
506 normaliseWeights(requireMatch_, true);
507
508 nonConformalCorrection();
509
510 upToDate_ = true;
511
512 return upToDate_;
513}
514
515
517{
519 os.writeEntryIfDifferent<scalar>("Cbb", 0.1, Cbb_);
520}
521
522
523// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Templated tree of axis-aligned bounding boxes (AABB)
Definition: AABBTree.H:73
const List< labelList > & addressing() const
Return the contents addressing.
Definition: AABBTree.C:449
const List< treeBoundBox > & boundBoxes() const
Return the bounding boxes making up the tree.
Definition: AABBTree.C:442
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:503
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
List< Key > toc() const
The table of contents (the keys) in unsorted order.
Definition: HashTable.C:122
Output to file stream, using an OSstream.
Definition: OFstream.H:57
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
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:251
A list of faces which address into the list of points.
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:427
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
@ nonBlocking
"nonBlocking"
Base class for Arbitrary Mesh Interface (AMI) methods.
const primitivePatch & tgtPatch() const
Return const access to the target patch.
const primitivePatch & srcPatch() const
Return const access to the source patch.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:64
bool overlaps(const boundBox &bb) const
Overlaps/touches boundingBox?
Definition: boundBoxI.H:221
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:97
tmp< pointField > points() const
Corner points in an order corresponding to a 'hex' cell.
Definition: boundBox.C:118
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Face area weighted Arbitrary Mesh Interface (AMI) method that performs the intersection of src and tg...
void storeInterArea(const label srcFacei, const label tgtFacei, DynamicList< label > &srcAddr, DynamicList< scalar > &srcWght, DynamicList< vector > &srcCtr, DynamicList< label > &tgtAddr, DynamicList< scalar > &tgtWght) const
virtual bool calculate(const primitivePatch &srcPatch, const primitivePatch &tgtPatch, const autoPtr< searchableSurface > &surfPtr=nullptr)
Update addressing, weights and (optional) centroids.
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.
void writeNoMatch(const label srcFacei, const labelList &tgtFaceCandidates, const boundBox &srcFaceBb) const
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
virtual bool write()
Write the output fields.
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
void inplaceToGlobal(labelUList &labels) const
From local to global index (inplace)
Definition: globalIndexI.H:303
Class containing processor-to-processor mapping information.
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition: one.H:62
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:89
2-D triangle and queries
Definition: triangle2D.H:55
label snapClosePoints(const triangle2D &triB)
Definition: triangle2D.C:93
static scalar relTol
Relative tolerance.
Definition: triangle2D.H:70
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
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))
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
Vector2D< scalar > vector2D
A 2D vector of scalars obtained from the generic Vector2D.
Definition: vector2D.H:51
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Type gAverage(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
#define addProfiling(name, descr)
Define profiling trigger with specified name and description string.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
List helper to append y elements onto the end of x.
Definition: ListOps.H:563
Functor to negate primitives. Dummy for most other types.
Definition: flipOp.H:69