Loading...
Searching...
No Matches
cyclicAMIPolyPatchTopologyChange.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) 2019-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
29#include "SubField.H"
30#include "vectorList.H"
31#include "polyTopoChange.H"
32
33// * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
34
36{
38
39 // Note: only used for topology update (createAMIFaces_ flag)
40 if (!createAMIFaces_)
41 {
43 << "Attempted to perform topology update when createAMIFaces_ "
44 << "flag is set to false"
45 << abort(FatalError);
46 }
47
48 if (boundaryMesh().mesh().hasCellVolumes())
49 {
51 << "Mesh already has volumes set!"
52 << endl;
53 }
54
57
59 << "Patch:" << name() << " before: sum(mag(faceAreas)):"
60 << gSum(mag(faceAreas)) << nl
61 << "Patch:" << name() << " before: sum(mag(faceAreas0)):"
62 << gSum(mag(faceAreas0_)) << endl;
63
66 {
67 DebugInfo << "Moving face centres" << endl;
69 }
70
71 faceAreas0_.clear();
72 faceCentres0_.clear();
73
75 << "Patch:" << name() << " after: sum(mag(faceAreas)):"
76 << gSum(mag(faceAreas)) << nl
77 << "Patch:" << name() << " after: sum(mag(faceAreas0)):"
78 << gSum(mag(faceAreas0_)) << endl;
79}
80
81
83{
85
86 // Note: only used for topology update (createAMIFaces_ flag)
87 if (!createAMIFaces_)
88 {
90 << "Attempted to perform topology update when createAMIFaces_ "
91 << "flag is set to false"
92 << abort(FatalError);
93 }
94
95 if (!owner())
96 {
97 return false;
98 }
99
100 bool changeRequired = false;
101
102 // Remove any faces that we inserted to create the 1-to-1 match...
103
104 const cyclicAMIPolyPatch& nbr = neighbPatch();
105
106 const label newSrcFaceStart = srcFaceIDs_.size();
107
108 if (newSrcFaceStart != 0)
109 {
110 for (label facei = newSrcFaceStart; facei < size(); ++facei)
111 {
112 changeRequired = true;
113 label meshFacei = start() + facei;
114 topoChange.removeFace(meshFacei, -1);
115 }
116 }
117
118 const label newTgtFaceStart = tgtFaceIDs_.size();
119
120 if (newTgtFaceStart != 0)
121 {
122 for (label facei = newTgtFaceStart; facei < nbr.size(); ++facei)
123 {
124 changeRequired = true;
125 label meshFacei = nbr.start() + facei;
126 topoChange.removeFace(meshFacei, -1);
127 }
128 }
129
130 srcFaceIDs_.clear();
131 tgtFaceIDs_.clear();
132
133 return changeRequired;
134}
135
136
138{
140
141 // Note: only used for topology update (createAMIFaces_ flag = true)
142 if (!createAMIFaces_)
143 {
145 << "Attempted to perform topology update when createAMIFaces_ "
146 << "flag is set to false"
147 << abort(FatalError);
148 }
149
150 bool changedFaces = false;
151 const cyclicAMIPolyPatch& nbr = neighbPatch();
152
153 polyMesh& mesh = const_cast<polyMesh&>(boundaryMesh().mesh());
154 const faceZoneMesh& faceZones = mesh.faceZones();
155
156 // First face address and weight are used to manipulate the
157 // original face - all other addresses and weights are used to
158 // create additional faces
159 const labelListList& srcToTgtAddr = AMI().srcAddress();
160 const labelListList& tgtToSrcAddr = AMI().tgtAddress();
161
162 const label nSrcFace = srcToTgtAddr.size();
163 const label nTgtFace = tgtToSrcAddr.size();
164
165 srcFaceIDs_.setSize(nSrcFace);
166 tgtFaceIDs_.setSize(nTgtFace);
167
168 label nNewSrcFaces = 0;
169 forAll(srcToTgtAddr, srcFacei)
170 {
171 const labelList& tgtAddr = srcToTgtAddr[srcFacei];
172
173 // No tgt faces linked to srcFacei (ACMI)
174 if (tgtAddr.empty()) continue;
175
176 srcFaceIDs_[srcFacei].setSize(tgtAddr.size());
177 srcFaceIDs_[srcFacei][0] = srcFacei;
178
179 const label meshFacei = start() + srcFacei;
180 for (label addri = 1; addri < tgtAddr.size(); ++addri)
181 {
182 changedFaces = true;
183
184 // Note: new faces reuse originating face points
185 // - but areas are scaled by the weights (later)
186
187 // New source face for each target face address
188 srcFaceIDs_[srcFacei][addri] = nNewSrcFaces + nSrcFace;
189 ++nNewSrcFaces;
190 (void)topoChange.addFace
191 (
192 mesh.faces()[meshFacei], // modified face
193 mesh.faceOwner()[meshFacei], // owner
194 -1, // neighbour
195 -1, // master point
196 -1, // master edge
197 meshFacei, // master face
198 false, // face flip
199 index(), // patch for face
200 faceZones.whichZone(meshFacei), // zone for original face
201 false // face flip in zone
202 );
203 }
204 }
205
206 label nNewTgtFaces = 0;
207 forAll(tgtToSrcAddr, tgtFacei)
208 {
209 const labelList& srcAddr = tgtToSrcAddr[tgtFacei];
210
211 // No src faces linked to tgtFacei (ACMI)
212 if (srcAddr.empty()) continue;
213
214 tgtFaceIDs_[tgtFacei].setSize(srcAddr.size());
215 tgtFaceIDs_[tgtFacei][0] = tgtFacei;
216
217 const label meshFacei = nbr.start() + tgtFacei;
218 for (label addri = 1; addri < srcAddr.size(); ++addri)
219 {
220 changedFaces = true;
221
222 // Note: new faces reuse originating face points
223 // - but areas are scaled by the weights (later)
224
225 // New target face for each source face address
226 tgtFaceIDs_[tgtFacei][addri] = nNewTgtFaces + nTgtFace;
227 ++nNewTgtFaces;
228
229 (void)topoChange.addFace
230 (
231 mesh.faces()[meshFacei], // modified face
232 mesh.faceOwner()[meshFacei], // owner
233 -1, // neighbour
234 -1, // master point
235 -1, // master edge
236 meshFacei, // master face
237 false, // face flip
238 nbr.index(), // patch for face
239 faceZones.whichZone(meshFacei), // zone for original face
240 false // face flip in zone
241 );
242 }
243 }
244
245 Info<< "AMI: Patch " << name() << " additional faces: "
246 << returnReduce(nNewSrcFaces, sumOp<label>()) << nl
247 << "AMI: Patch " << nbr.name() << " additional faces: "
248 << returnReduce(nNewTgtFaces, sumOp<label>())
249 << endl;
250
251 if (debug)
252 {
253 Pout<< "New faces - " << name() << ": " << nNewSrcFaces
254 << " " << nbr.name() << ": " << nNewTgtFaces << endl;
255 }
256
257 return returnReduceOr(changedFaces);
258}
259
260
262{
263 // Create new mesh faces so that there is a 1-to-1 correspondence
264 // between faces on each side of the AMI
265
266 // Note: only used for topology update (createAMIFaces_ flag)
267 if (!createAMIFaces_)
268 {
270 << "Attempted to perform topology update when createAMIFaces_ "
271 << "flag is set to false"
272 << abort(FatalError);
273 }
274
275
277
278 if (!owner())
279 {
280 return;
281 }
282
283 const cyclicAMIPolyPatch& nbr = neighbPatch();
284
285 vectorField& nbrFaceAreas0 = nbr.faceAreas0();
286 vectorField& nbrFaceCentres0 = nbr.faceCentres0();
287
288 // Scale the new face areas and set the centroids
289 // Note:
290 // - storing local copies so that they can be re-applied after the call to
291 // movePoints that will reset any changes to the areas and centroids
292 //
293 // - For AMI, src and tgt patches should be the same
294 // - For ACMI they are likely to be different!
295 faceAreas0_ = faceAreas();
296 faceCentres0_ = faceCentres();
297 nbrFaceAreas0 = nbr.faceAreas();
298 nbrFaceCentres0 = nbr.faceCentres();
299
300 // Original AMI info (based on the mesh state when the AMI was evaluated)
301 const labelListList& srcToTgtAddr0 = AMIPtr_->srcAddress();
302 const labelListList& tgtToSrcAddr0 = AMIPtr_->tgtAddress();
303 const pointListList& srcCtr0 = AMIPtr_->srcCentroids();
304 const scalarListList& srcToTgtWght0 = AMIPtr_->srcWeights();
305 const label singlePatchProc = AMIPtr_->singlePatchProc();
306
307 // New addressing on new mesh (extended by polyTopoChange)
308 labelListList srcToTgtAddr1(size(), labelList());
309 labelListList tgtToSrcAddr1(nbr.size(), labelList());
310
311 // Need to calc new parallel maps (mesh has changed since AMI was computed)
312 autoPtr<mapDistribute> srcToTgtMap1;
313 autoPtr<mapDistribute> tgtToSrcMap1;
314
315 if (AMIPtr_->distributed() && AMIPtr_().comm() != -1)
316 {
317 // Parallel running
318
319 // Global index based on old patch sizes (when AMI was computed)
320 globalIndex globalSrcFaces0(srcToTgtAddr0.size(), AMIPtr_().comm());
321 globalIndex globalTgtFaces0(tgtToSrcAddr0.size(), AMIPtr_().comm());
322
323 // Global index based on new patch sizes
324 globalIndex globalSrcFaces1(size(), AMIPtr_().comm());
325 globalIndex globalTgtFaces1(nbr.size(), AMIPtr_().comm());
326
327
328 // Gather source side info
329 // =======================
330
331 // Note: using new global index for addressing, and distributed using
332 // the old AMI map
333 labelListList newTgtGlobalFaces(tgtFaceIDs_);
334 forAll(newTgtGlobalFaces, tgtFacei)
335 {
336 globalTgtFaces1.inplaceToGlobal(newTgtGlobalFaces[tgtFacei]);
337 }
338 AMIPtr_->tgtMap().distribute(newTgtGlobalFaces);
339
340 // Now have new tgt face indices for each src face
341
342 labelList globalSrcFaceIDs(identity(srcToTgtAddr0.size()));
343 globalSrcFaces0.inplaceToGlobal(globalSrcFaceIDs);
344 AMIPtr_->srcMap().distribute(globalSrcFaceIDs);
345 // globalSrcFaceIDs now has remote data for each srcFacei0 known to the
346 // tgt patch
347
348 List<List<point>> globalSrcCtrs0(srcCtr0);
349 AMIPtr_->srcMap().distribute(globalSrcCtrs0);
350
351 labelList globalTgtFaceIDs(identity(tgtToSrcAddr0.size()));
352 globalTgtFaces0.inplaceToGlobal(globalTgtFaceIDs);
353 AMIPtr_->tgtMap().distribute(globalTgtFaceIDs);
354 // globalTgtFaceIDs now has remote data for each tgtFacei0 known to the
355 // src patch
356
357 // For debug - send tgt face centres and compare against mapped src
358 // face centres
359 //List<List<point>> globalTgtCtrs0(tgtCtr0);
360 //AMIPtr_->tgtMap().distribute(globalTgtCtrs0);
361
362 labelListList globalTgtToSrcAddr(tgtToSrcAddr0);
363 forAll(tgtToSrcAddr0, tgtFacei0)
364 {
365 forAll(tgtToSrcAddr0[tgtFacei0], addri)
366 {
367 const label globalSrcFacei =
368 globalSrcFaceIDs[tgtToSrcAddr0[tgtFacei0][addri]];
369 globalTgtToSrcAddr[tgtFacei0][addri] = globalSrcFacei;
370 }
371 }
372 AMIPtr_->tgtMap().distribute(globalTgtToSrcAddr);
373
374 labelListList globalSrcToTgtAddr(srcToTgtAddr0);
375 forAll(srcToTgtAddr0, srcFacei0)
376 {
377 forAll(srcToTgtAddr0[srcFacei0], addri)
378 {
379 const label globalTgtFacei =
380 globalTgtFaceIDs[srcToTgtAddr0[srcFacei0][addri]];
381 globalSrcToTgtAddr[srcFacei0][addri] = globalTgtFacei;
382 }
383 }
384 AMIPtr_->srcMap().distribute(globalSrcToTgtAddr);
385
386 label nError = 0;
387 forAll(srcToTgtAddr0, srcFacei0)
388 {
389 const labelList& newSrcFaces = srcFaceIDs_[srcFacei0];
390
391 forAll(newSrcFaces, i)
392 {
393 const label srcFacei1 = newSrcFaces[i];
394
395 // What index did srcFacei0 appear in tgtToSrc0 list?
396 // - if first index, all ok
397 // - else tgt face has been moved to according to tgtFaceIDs_
398 const label tgtFacei0 = srcToTgtAddr0[srcFacei0][i];
399 const label addri =
400 globalTgtToSrcAddr[tgtFacei0].find
401 (
402 globalSrcFaceIDs[srcFacei0]
403 );
404
405 if (addri == -1)
406 {
407 ++nError;
408 continue;
409
410 if (debug)
411 {
412 Pout<< "Unable to find global source face "
413 << globalSrcFaceIDs[srcFacei0]
414 << " in globalTgtToSrcAddr[" << tgtFacei0 << "]: "
415 << globalTgtToSrcAddr[tgtFacei0]
416 << endl;
417 }
418 }
419
420 const label tgtFacei1 = newTgtGlobalFaces[tgtFacei0][addri];
421
422 // Sanity check to see that we've picked the correct face
423 // point tgtCtr0(globalTgtCtrs0[tgtFacei0][addri]);
424 // Pout<< "srcCtr:" << srcCtr0[srcFacei0][i]
425 // << " tgtCtr:" << tgtCtr0 << endl;
426
427 srcToTgtAddr1[srcFacei1] = labelList(1, tgtFacei1);
428 faceAreas0_[srcFacei1] *= srcToTgtWght0[srcFacei0][i];
429 faceCentres0_[srcFacei1] = srcCtr0[srcFacei0][i];
430 }
431 }
432
433 if (nError)
434 {
436 << "Unable to find " << nError << " global source faces"
437 << abort(FatalError);
438 }
439
440
441 // Gather Target side info
442 // =======================
443
444 labelListList newSrcGlobalFaces(srcFaceIDs_);
445 forAll(newSrcGlobalFaces, srcFacei)
446 {
447 globalSrcFaces1.inplaceToGlobal(newSrcGlobalFaces[srcFacei]);
448 }
449
450 AMIPtr_->srcMap().distribute(newSrcGlobalFaces);
451
452 // Now have new src face indices for each tgt face
453 forAll(tgtToSrcAddr0, tgtFacei0)
454 {
455 const labelList& newTgtFaces = tgtFaceIDs_[tgtFacei0];
456 forAll(newTgtFaces, i)
457 {
458 const label srcFacei0 = tgtToSrcAddr0[tgtFacei0][i];
459
460 const label addri =
461 globalSrcToTgtAddr[srcFacei0].find
462 (
463 globalTgtFaceIDs[tgtFacei0]
464 );
465
466 if (addri == -1)
467 {
468 ++nError;
469 continue;
470
471 if (debug)
472 {
473 Pout<< "Unable to find global target face "
474 << globalTgtFaceIDs[tgtFacei0]
475 << " in globalSrcToTgtAddr[" << srcFacei0 << "]: "
476 << globalSrcToTgtAddr[srcFacei0]
477 << endl;
478 }
479 }
480
481 const label srcFacei1 = newSrcGlobalFaces[srcFacei0][addri];
482
483 // Sanity check to see that we've picked the correct face
484 point srcCtr0(globalSrcCtrs0[srcFacei0][addri]);
485 reverseTransformPosition(srcCtr0, srcFacei0);
486
487 const label tgtFacei1 = newTgtFaces[i];
488 tgtToSrcAddr1[tgtFacei1] = labelList(1, srcFacei1);
489 nbrFaceCentres0[tgtFacei1] = srcCtr0;
490 }
491 }
492
493 if (nError)
494 {
496 << "Unable to find " << nError << " global target faces"
497 << abort(FatalError);
498 }
499
500 // Update the maps
501 {
502 List<Map<label>> cMap;
503 srcToTgtMap1.reset
504 (
505 new mapDistribute(globalSrcFaces1, tgtToSrcAddr1, cMap)
506 );
507 }
508 {
509 List<Map<label>> cMap;
510 tgtToSrcMap1.reset
511 (
512 new mapDistribute(globalTgtFaces1, srcToTgtAddr1, cMap)
513 );
514 }
515
516 // Reset tgt patch areas using the new map
517 vectorList newSrcGlobalFaceAreas(faceAreas0_);
518
519 srcToTgtMap1->distribute(newSrcGlobalFaceAreas);
520 forAll(nbrFaceAreas0, tgtFacei)
521 {
522 if (!tgtToSrcAddr1[tgtFacei].empty())
523 {
524 const label srcFacei = tgtToSrcAddr1[tgtFacei][0];
525 nbrFaceAreas0[tgtFacei] = -newSrcGlobalFaceAreas[srcFacei];
526 }
527 }
528 }
529 else
530 {
531 label nError = 0;
532 forAll(srcToTgtAddr0, srcFacei0)
533 {
534 const labelList& srcFaceTgtAddr0 = srcToTgtAddr0[srcFacei0];
535 const scalarList& srcFaceTgtWght0 = srcToTgtWght0[srcFacei0];
536 const pointList& srcFaceTgtCtr0 = srcCtr0[srcFacei0];
537 forAll(srcFaceTgtAddr0, addri)
538 {
539 const label srcFacei1 = srcFaceIDs_[srcFacei0][addri];
540
541 // Find which slot srcFacei0 appears in tgt->src addressing
542 const label tgtFacei0 = srcFaceTgtAddr0[addri];
543 const label tgtAddri0 =
544 tgtToSrcAddr0[tgtFacei0].find(srcFacei0);
545
546 if (tgtAddri0 == -1)
547 {
548 ++nError;
549 continue;
550
551 if (debug)
552 {
553 Pout<< "Unable to find source face " << srcFacei0
554 << " in tgtToSrcAddr0[" << tgtFacei0 << "]: "
555 << tgtToSrcAddr0[tgtFacei0]
556 << endl;
557 }
558 }
559
560 const label tgtFacei1 = tgtFaceIDs_[tgtFacei0][tgtAddri0];
561
562 faceAreas0_[srcFacei1] *= srcFaceTgtWght0[addri];
563 nbrFaceAreas0[tgtFacei1] = -faceAreas0_[srcFacei1];
564
565 point pt(srcFaceTgtCtr0[addri]);
566 faceCentres0_[srcFacei1] = pt;
567 reverseTransformPosition(pt, srcFacei0);
568 nbrFaceCentres0[tgtFacei1] = pt;
569
570 // SANITY CHECK
571 // Info<< "srcPt:" << srcFaceCentres[srcFacei1]
572 // << " tgtPt:" << tgtFaceCentres[tgtFacei1] << endl;
573
574 srcToTgtAddr1[srcFacei1] = labelList(1, tgtFacei1);
575 tgtToSrcAddr1[tgtFacei1] = labelList(1, srcFacei1);
576 }
577 }
578
579 if (nError)
580 {
582 << "Unable to find " << nError
583 << " source faces in tgtToSrcAddr0"
584 << abort(FatalError);
585 }
586 }
587
588 scalarListList newSrcToTgtWeights(srcToTgtAddr1.size());
589 forAll(srcToTgtAddr1, facei)
590 {
591 if (srcToTgtAddr1[facei].size())
592 {
593 newSrcToTgtWeights[facei] = scalarList(1, scalar(1));
594 }
595 else
596 {
597 // No connection - effect of face removed by setting area to a
598 // a small value
599 faceAreas0_[facei] *= tolerance_;
600 }
601 }
602
603 scalarListList newTgtToSrcWeights(tgtToSrcAddr1.size());
604 forAll(tgtToSrcAddr1, facei)
605 {
606 if (tgtToSrcAddr1[facei].size())
607 {
608 newTgtToSrcWeights[facei] = scalarList(1, scalar(1));
609 }
610 else
611 {
612 // No connection - effect of face removed by setting area to a
613 // a small value
614 nbrFaceAreas0[facei] *= tolerance_;
615 }
616 }
617
618 // Reset the AMI addressing and weights to reflect the new 1-to-1
619 // correspondence
620 AMIPtr_->reset
621 (
622 std::move(srcToTgtMap1),
623 std::move(tgtToSrcMap1),
624 std::move(srcToTgtAddr1),
625 std::move(newSrcToTgtWeights),
626 std::move(tgtToSrcAddr1),
627 std::move(newTgtToSrcWeights),
628 singlePatchProc
629 );
630
631 // Need to set areas, e.g. for agglomeration to (re-)normalisation weights
632 AMIPtr_->srcMagSf() = mag(faceAreas0_);
633 AMIPtr_->tgtMagSf() = mag(nbrFaceAreas0);
634
635 if (debug)
636 {
637 Pout<< "cyclicAMIPolyPatch : " << name()
638 << " constructed AMI with " << nl
639 << " " << "srcAddress:" << AMIPtr_().srcAddress().size()
640 << nl
641 << " " << "tgAddress :" << AMIPtr_().tgtAddress().size()
642 << nl << endl;
643 }
644}
645
646
647// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
648
650{
653 createAMIFaces_ = true;
654
655 return true;
656}
657
658
660{
662
663 if (createAMIFaces_ && owner())
664 {
665 // Calculate the AMI using the new points
666 // Note: mesh still has old points
667 resetAMI(topoChange.points());
668
669 removeAMIFaces(topoChange);
670
671 addAMIFaces(topoChange);
672
673 return true;
674 }
675
676 return false;
677}
678
679
680// ************************************************************************* //
SubField< vector > subField
Definition Field.H:183
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void setSize(label n)
Alias for resize().
Definition List.H:536
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
void size(const label n)
Definition UList.H:118
label find(const T &val) const
Find index of the first occurrence of the value.
Definition UList.C:160
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition ZoneMesh.C:410
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Definition autoPtrI.H:37
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
const bMesh & mesh() const
Cyclic patch for Arbitrary Mesh Interface (AMI).
vectorField & faceAreas0() const
Return access to the initial face areas.
bool createAMIFaces_
Flag to indicate that new AMI faces will created.
virtual bool owner() const
Does this side own the patch?
virtual bool addAMIFaces(polyTopoChange &topoChange)
Collect faces to add in the topoChange container.
vectorField faceCentres0_
Temporary storage for AMI face centres.
virtual void restoreScaledGeometry()
Helper to re-apply the geometric scaling lost during mesh updates.
virtual bool removeAMIFaces(polyTopoChange &topoChange)
Collect faces to remove in the topoChange container.
virtual bool setTopology(polyTopoChange &topoChange)
Set topology changes in the polyTopoChange object.
vectorField faceAreas0_
Temporary storage for AMI face areas.
const AMIPatchToPatchInterpolation & AMI() const
Return a reference to the AMI interpolator.
virtual bool changeTopology() const
Return true if this patch changes the mesh topology.
virtual const cyclicAMIPolyPatch & neighbPatch() const
Return a reference to the neighbour patch.
virtual void setAMIFaces()
Set properties of newly inserted faces after topological changes.
static const scalar tolerance_
Tolerance used e.g. for area calculations/limits.
vectorField & faceCentres0() const
Return access to the initial face centres.
cyclicAMIPolyPatch(const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm, const word &patchType, const transformType transform=UNKNOWN, const word &defaultAMIMethod=faceAreaWeightAMI::typeName)
Construct from (base coupled patch) components.
virtual void resetAMI(const UList< point > &points) const
Reset the AMI interpolator, supply patch points.
virtual void reverseTransformPosition(point &l, const label facei) const
Transform a patch-based position from this side to nbr side.
autoPtr< AMIPatchToPatchInterpolation > AMIPtr_
AMI interpolation class.
bool moveFaceCentres_
Move face centres (default = no).
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
void inplaceToGlobal(const label proci, labelUList &labels) const
From local to global index on proci (inplace).
Class containing processor-to-processor mapping information.
label index() const noexcept
The index of this patch in the boundaryMesh.
const word & name() const noexcept
The patch name.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
virtual const faceList & faces() const
Return raw faces.
Definition polyMesh.C:1088
virtual const labelList & faceOwner() const
Return face owner.
Definition polyMesh.C:1101
const vectorField::subField faceAreas() const
Return face normals.
Definition polyPatch.C:326
const vectorField::subField faceCentres() const
Return face centres.
Definition polyPatch.C:320
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundaryMesh reference.
Definition polyPatch.C:295
label start() const noexcept
Return start label of this patch in the polyMesh face list.
Definition polyPatch.H:446
Direct mesh changes based on v1.3 polyTopoChange syntax.
void removeFace(const label facei, const label mergeFacei)
Remove/merge face.
label addFace(const face &f, const label own, const label nei, const label masterPointID, const label masterEdgeID, const label masterFaceID, const bool flipFaceFlux, const label patchID, const label zoneID, const bool zoneFlip)
Add face to cells. Return new face label.
const DynamicList< point > & points() const
Points. Shrunk after constructing mesh (or calling of compact()).
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
auto & name
#define DebugInfo
Report an information message using Foam::Info.
#define WarningInFunction
Report a warning using Foam::Warning.
#define DebugInFunction
Report an information message using Foam::Info.
Namespace for handling debugging switches.
Definition debug.C:45
List< scalarList > scalarListList
List of scalarList.
Definition scalarList.H:35
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
Type gSum(const FieldField< Field, Type > &f)
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
messageStream Info
Information stream (stdout output on master, null elsewhere).
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with faceZone content on a polyMesh.
List< vector > vectorList
List of vector.
Definition vectorList.H:32
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
List< point > pointList
List of point.
Definition pointList.H:32
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
List< pointList > pointListList
List of pointList.
Definition pointList.H:35
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299