cyclicPeriodicAMIPolyPatch.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) 2015-2021 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
30
31// For debugging
32#include "OBJstream.H"
33#include "PatchTools.H"
34#include "Time.H"
35// Note: cannot use vtkSurfaceWriter here - circular linkage
36// but foamVtkSurfaceWriter (vtk::surfaceWriter) would be okay.
37//
38// #include "foamVtkSurfaceWriter.H"
39
40// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41
42namespace Foam
43{
45
48 (
52 );
53}
54
55
56// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
57
58void Foam::cyclicPeriodicAMIPolyPatch::syncTransforms() const
59{
60 if (owner())
61 {
62 // At this point we guarantee that the transformations have been
63 // updated. There is one particular case, where if the periodic patch
64 // locally has zero faces then its transformation will not be set. We
65 // need to synchronise the transforms over the zero-sized patches as
66 // well.
67 //
68 // We can't put the logic into cyclicPolyPatch as
69 // processorCyclicPolyPatch uses cyclicPolyPatch functionality.
70 // Synchronisation will fail because processor-type patches do not exist
71 // on all processors.
72 //
73 // The use in cyclicPeriodicAMI is special; we use the patch even
74 // though we have no faces in it. Ideally, in future, the transformation
75 // logic will be abstracted out, and we won't need a periodic patch
76 // here. Until then, we can just fix the behaviour of the zero-sized
77 // coupled patches here
78
79 // Get the periodic patch
80 const coupledPolyPatch& periodicPatch
81 (
82 refCast<const coupledPolyPatch>
83 (
84 boundaryMesh()[periodicPatchID()]
85 )
86 );
87
88 // If there are any zero-sized periodic patches
89 if (returnReduce((size() && !periodicPatch.size()), orOp<bool>()))
90 {
91 if (periodicPatch.separation().size() > 1)
92 {
94 << "Periodic patch " << periodicPatchName_
95 << " has non-uniform separation vector "
96 << periodicPatch.separation()
97 << "This is not allowed inside " << type()
98 << " patch " << name()
99 << exit(FatalError);
100 }
101
102 if (periodicPatch.forwardT().size() > 1)
103 {
105 << "Periodic patch " << periodicPatchName_
106 << " has non-uniform transformation tensor "
107 << periodicPatch.forwardT()
108 << "This is not allowed inside " << type()
109 << " patch " << name()
110 << exit(FatalError);
111 }
112
113 // Note that zero-sized patches will have zero-sized fields for the
114 // separation vector, forward and reverse transforms. These need
115 // replacing with the transformations from other processors.
116
117 // Parallel in this context refers to a parallel transformation,
118 // rather than a rotational transformation.
119
120 // Note that a cyclic with zero faces is considered parallel so
121 // explicitly check for that.
122 bool isParallel =
123 (
124 periodicPatch.size()
125 && periodicPatch.parallel()
126 );
127 reduce(isParallel, orOp<bool>());
128
129 if (isParallel)
130 {
131 // Sync a list of separation vectors
132 List<vectorField> sep(Pstream::nProcs());
133 sep[Pstream::myProcNo()] = periodicPatch.separation();
135
136 List<boolList> coll(Pstream::nProcs());
137 coll[Pstream::myProcNo()] = periodicPatch.collocated();
139
140 // If locally we have zero faces pick the first one that has a
141 // separation vector
142 if (!periodicPatch.size())
143 {
144 forAll(sep, procI)
145 {
146 if (sep[procI].size())
147 {
148 const_cast<vectorField&>
149 (
150 periodicPatch.separation()
151 ) = sep[procI];
152 const_cast<boolList&>
153 (
154 periodicPatch.collocated()
155 ) = coll[procI];
156
157 break;
158 }
159 }
160 }
161 }
162 else
163 {
164 // Sync a list of forward and reverse transforms
165 List<tensorField> forwardT(Pstream::nProcs());
166 forwardT[Pstream::myProcNo()] = periodicPatch.forwardT();
168
169 List<tensorField> reverseT(Pstream::nProcs());
170 reverseT[Pstream::myProcNo()] = periodicPatch.reverseT();
172
173 // If locally we have zero faces pick the first one that has a
174 // transformation vector
175 if (!periodicPatch.size())
176 {
177 forAll(forwardT, procI)
178 {
179 if (forwardT[procI].size())
180 {
181 const_cast<tensorField&>
182 (
183 periodicPatch.forwardT()
184 ) = forwardT[procI];
185 const_cast<tensorField&>
186 (
187 periodicPatch.reverseT()
188 ) = reverseT[procI];
189
190 break;
191 }
192 }
193 }
194 }
195 }
196 }
197}
198
199
200void Foam::cyclicPeriodicAMIPolyPatch::writeOBJ
201(
202 const primitivePatch& p,
203 OBJstream& str
204) const
205{
206 // Collect faces and points
208 faceList allFaces;
209 labelList pointMergeMap;
211 (
212 -1.0, // do not merge points
213 p,
214 allPoints,
215 allFaces,
216 pointMergeMap
217 );
218
219 if (Pstream::master())
220 {
221 // Write base geometry
222 str.write(allFaces, allPoints);
223 }
224}
225
226
227void Foam::cyclicPeriodicAMIPolyPatch::resetAMI() const
228{
229 if (owner())
230 {
231 // Get the periodic patch
232 const coupledPolyPatch& periodicPatch
233 (
234 refCast<const coupledPolyPatch>
235 (
236 boundaryMesh()[periodicPatchID()]
237 )
238 );
239
240 // Synchronise the transforms
241 syncTransforms();
242
243 // Create copies of both patches' points, transformed to the owner
244 pointField thisPoints0(localPoints());
245 pointField nbrPoints0(neighbPatch().localPoints());
246 transformPosition(nbrPoints0);
247
248 // Reset the stored number of periodic transformations to a lower
249 // absolute value if possible
250 if (nSectors_ > 0)
251 {
252 if (nTransforms_ > nSectors_/2)
253 {
254 nTransforms_ -= nSectors_;
255 }
256 else if (nTransforms_ < - nSectors_/2)
257 {
258 nTransforms_ += nSectors_;
259 }
260 }
261
262 // Apply the stored number of periodic transforms
263 for (label i = 0; i < nTransforms_; ++ i)
264 {
265 periodicPatch.transformPosition(thisPoints0);
266 }
267 for (label i = 0; i > nTransforms_; -- i)
268 {
269 periodicPatch.transformPosition(nbrPoints0);
270 }
271
272 autoPtr<OBJstream> ownStr;
273 autoPtr<OBJstream> neiStr;
274 if (debug)
275 {
276 const Time& runTime = boundaryMesh().mesh().time();
277
278 fileName dir(runTime.globalPath());
279 fileName postfix("_" + runTime.timeName()+"_expanded.obj");
280
281 ownStr.reset(new OBJstream(dir/name() + postfix));
282 neiStr.reset(new OBJstream(dir/neighbPatch().name() + postfix));
283
285 << "patch:" << name()
286 << " writing accumulated AMI to " << ownStr().name()
287 << " and " << neiStr().name() << endl;
288 }
289
290 // Create another copy
291 pointField thisPoints(thisPoints0);
292 pointField nbrPoints(nbrPoints0);
293
294 // Create patches for all the points
295
296 // Source patch at initial location
297 const primitivePatch thisPatch0
298 (
299 SubList<face>(localFaces(), size()),
300 thisPoints0
301 );
302 // Source patch that gets moved
303 primitivePatch thisPatch
304 (
305 SubList<face>(localFaces(), size()),
306 thisPoints
307 );
308 // Target patch at initial location
309 const primitivePatch nbrPatch0
310 (
311 SubList<face>(neighbPatch().localFaces(), neighbPatch().size()),
312 nbrPoints0
313 );
314 // Target patch that gets moved
315 primitivePatch nbrPatch
316 (
317 SubList<face>(neighbPatch().localFaces(), neighbPatch().size()),
318 nbrPoints
319 );
320
321 // Construct a new AMI interpolation between the initial patch locations
322 AMIPtr_->setRequireMatch(false);
323 AMIPtr_->calculate(thisPatch0, nbrPatch0, surfPtr());
324
325 // Number of geometry replications
326 label iter(0);
327 label nTransformsOld(nTransforms_);
328
329 if (ownStr)
330 {
331 writeOBJ(thisPatch0, *ownStr);
332 }
333 if (neiStr)
334 {
335 writeOBJ(nbrPatch0, *neiStr);
336 }
337
338
339 // Weight sum averages
340 scalar srcSum(gAverage(AMIPtr_->srcWeightsSum()));
341 scalar tgtSum(gAverage(AMIPtr_->tgtWeightsSum()));
342 // Direction (or rather side of AMI : this or nbr patch) of
343 // geometry replication
344 bool direction = nTransforms_ >= 0;
345
346 // Increase in the source weight sum for the last iteration in the
347 // opposite direction. If the current increase is less than this, the
348 // direction (= side of AMI to transform) is reversed.
349 // We switch the side to replicate instead of reversing the transform
350 // since at the coupledPolyPatch level there is no
351 // 'reverseTransformPosition' functionality.
352 scalar srcSumDiff = 0;
353
355 << "patch:" << name()
356 << " srcSum:" << srcSum
357 << " tgtSum:" << tgtSum
358 << " direction:" << direction
359 << endl;
360
361 // Loop, replicating the geometry
362 while
363 (
364 (iter < maxIter_)
365 && (
366 (1 - srcSum > matchTolerance())
367 || (1 - tgtSum > matchTolerance())
368 )
369 )
370 {
371 if (direction)
372 {
373 periodicPatch.transformPosition(thisPoints);
374
376 << "patch:" << name()
377 << " moving this side from:"
378 << gAverage(thisPatch.points())
379 << " to:" << gAverage(thisPoints) << endl;
380
381 thisPatch.movePoints(thisPoints);
382
384 << "patch:" << name()
385 << " appending weights with untransformed slave side"
386 << endl;
387
388 AMIPtr_->append(thisPatch, nbrPatch0);
389
390 if (ownStr)
391 {
392 writeOBJ(thisPatch, *ownStr);
393 }
394 }
395 else
396 {
397 periodicPatch.transformPosition(nbrPoints);
398
400 << "patch:" << name()
401 << " moving neighbour side from:"
402 << gAverage(nbrPatch.points())
403 << " to:" << gAverage(nbrPoints) << endl;
404
405 nbrPatch.movePoints(nbrPoints);
406
407 AMIPtr_->append(thisPatch0, nbrPatch);
408
409 if (neiStr)
410 {
411 writeOBJ(nbrPatch, *neiStr);
412 }
413 }
414
415 const scalar srcSumNew = gAverage(AMIPtr_->srcWeightsSum());
416 const scalar srcSumDiffNew = srcSumNew - srcSum;
417
418 if (srcSumDiffNew < srcSumDiff || srcSumDiffNew < SMALL)
419 {
421
422 srcSumDiff = srcSumDiffNew;
423 }
424
425 srcSum = srcSumNew;
426 tgtSum = gAverage(AMIPtr_->tgtWeightsSum());
427
428 nTransforms_ += direction ? +1 : -1;
429
430 ++iter;
431
433 << "patch:" << name()
434 << " iteration:" << iter
435 << " srcSum:" << srcSum
436 << " tgtSum:" << tgtSum
437 << " direction:" << direction
438 << endl;
439 }
440
441
442 // Close debug streams
443 ownStr.reset(nullptr);
444 neiStr.reset(nullptr);
445
446
447 // Average the number of transformations
448 nTransforms_ = (nTransforms_ + nTransformsOld)/2;
449
450 // Check that the match is complete
451 if (iter == maxIter_)
452 {
453 // The matching algorithm has exited without getting the
454 // srcSum and tgtSum above 1. This can happen because
455 // - of an incorrect setup
456 // - or because of non-exact meshes and truncation errors
457 // (transformation, accumulation of cutting errors)
458 // so for now this situation is flagged as a SeriousError instead of
459 // a FatalError since the default matchTolerance is quite strict
460 // (0.001) and can get triggered far into the simulation.
462 << "Patches " << name() << " and " << neighbPatch().name()
463 << " do not couple to within a tolerance of "
464 << matchTolerance()
465 << " when transformed according to the periodic patch "
466 << periodicPatch.name() << "." << nl
467 << "The current sum of weights are for owner " << name()
468 << " : " << srcSum << " and for neighbour "
469 << neighbPatch().name() << " : " << tgtSum << nl
470 << "This is only acceptable during post-processing"
471 << "; not during running. Improve your mesh or increase"
472 << " the 'matchTolerance' setting in the patch specification."
473 << endl;
474 }
475
476 // Check that both patches have replicated an integer number of times
477 if
478 (
479 mag(srcSum - floor(srcSum + 0.5)) > srcSum*matchTolerance()
480 || mag(tgtSum - floor(tgtSum + 0.5)) > tgtSum*matchTolerance()
481 )
482 {
483 // This condition is currently enforced until there is more
484 // experience with the matching algorithm and numerics.
485 // This check means that e.g. different numbers of stator and
486 // rotor partitions are not allowed.
487 // Again see the section above about tolerances.
489 << "Patches " << name() << " and " << neighbPatch().name()
490 << " do not overlap an integer number of times when transformed"
491 << " according to the periodic patch "
492 << periodicPatch.name() << "." << nl
493 << "The current matchTolerance : " << matchTolerance()
494 << ", sum of owner weights : " << srcSum
495 << ", sum of neighbour weights : " << tgtSum
496 << "." << nl
497 << "This is only acceptable during post-processing"
498 << "; not during running. Improve your mesh or increase"
499 << " the 'matchTolerance' setting in the patch specification."
500 << endl;
501 }
502
503 // Print some statistics
504 const label nFace = returnReduce(size(), sumOp<label>());
505
506 if (nFace)
507 {
508 scalarField srcWghtSum(size(), Zero);
509 forAll(srcWghtSum, faceI)
510 {
511 srcWghtSum[faceI] = sum(AMIPtr_->srcWeights()[faceI]);
512 }
513 scalarField tgtWghtSum(neighbPatch().size(), Zero);
514 forAll(tgtWghtSum, faceI)
515 {
516 tgtWghtSum[faceI] = sum(AMIPtr_->tgtWeights()[faceI]);
517 }
518
519 Info<< indent
520 << "AMI: Patch " << name()
521 << " sum(weights)"
522 << " min:" << gMin(srcWghtSum)
523 << " max:" << gMax(srcWghtSum)
524 << " average:" << gAverage(srcWghtSum) << nl;
525 Info<< indent
526 << "AMI: Patch " << neighbPatch().name()
527 << " sum(weights)"
528 << " min:" << gMin(tgtWghtSum)
529 << " max:" << gMax(tgtWghtSum)
530 << " average:" << gAverage(tgtWghtSum) << nl;
531 }
532 }
533}
534
535
536// * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
537
539(
540 const word& name,
541 const label size,
542 const label start,
543 const label index,
544 const polyBoundaryMesh& bm,
545 const word& patchType,
547)
548:
550 (
551 name,
552 size,
553 start,
554 index,
555 bm,
556 patchType,
557 transform,
558 faceAreaWeightAMI::typeName
559 ),
560 nTransforms_(0),
561 nSectors_(0),
562 maxIter_(36)
563{
564 AMIPtr_->setRequireMatch(false);
565}
566
567
569(
570 const word& name,
571 const dictionary& dict,
572 const label index,
573 const polyBoundaryMesh& bm,
574 const word& patchType
575)
576:
578 (
579 name,
580 dict,
581 index,
582 bm,
583 patchType,
584 faceAreaWeightAMI::typeName
585 ),
586 nTransforms_(dict.getOrDefault<label>("nTransforms", 0)),
587 nSectors_(dict.getOrDefault<label>("nSectors", 0)),
588 maxIter_(dict.getOrDefault<label>("maxIter", 36))
589{
590 AMIPtr_->setRequireMatch(false);
591}
592
593
595(
597 const polyBoundaryMesh& bm
598)
599:
600 cyclicAMIPolyPatch(pp, bm),
601 nTransforms_(pp.nTransforms_),
602 nSectors_(pp.nSectors_),
603 maxIter_(pp.maxIter_)
604{
605 AMIPtr_->setRequireMatch(false);
606}
607
608
610(
612 const polyBoundaryMesh& bm,
613 const label index,
614 const label newSize,
615 const label newStart,
616 const word& nbrPatchName
617)
618:
619 cyclicAMIPolyPatch(pp, bm, index, newSize, newStart, nbrPatchName),
620 nTransforms_(pp.nTransforms_),
621 nSectors_(pp.nSectors_),
622 maxIter_(pp.maxIter_)
623{
624 AMIPtr_->setRequireMatch(false);
625}
626
627
629(
631 const polyBoundaryMesh& bm,
632 const label index,
633 const labelUList& mapAddressing,
634 const label newStart
635)
636:
637 cyclicAMIPolyPatch(pp, bm, index, mapAddressing, newStart),
638 nTransforms_(pp.nTransforms_),
639 nSectors_(pp.nSectors_),
640 maxIter_(pp.maxIter_)
641{
642 AMIPtr_->setRequireMatch(false);
643}
644
645
646// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
647
649{}
650
651
652// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
653
655{
657
658 os.writeEntryIfDifferent<label>("nTransforms", 0, nTransforms_);
659 os.writeEntryIfDifferent<label>("nSectors", 0, nSectors_);
660 os.writeEntryIfDifferent<label>("maxIter", 36, maxIter_);
661}
662
663
664// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
virtual const fileName & name() const
Return the name of the stream.
Definition: IOstream.C:40
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
static void gatherAndMerge(const scalar mergeDist, const PrimitivePatch< FaceList, PointField > &pp, Field< typename PrimitivePatch< FaceList, PointField >::point_type > &mergedPoints, List< typename PrimitivePatch< FaceList, PointField >::face_type > &mergedFaces, labelList &pointMergeMap, const bool useLocal=false)
Gather points and faces onto master and merge into single patch.
label nProcs() const noexcept
Number of ranks associated with PstreamBuffers.
static void allGatherList(const List< commsStruct > &comms, List< T > &values, const int tag, const label comm)
fileName globalPath() const
Return global path for the case.
Definition: TimePathsI.H:80
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
virtual const tensorField & reverseT() const
Return neighbour-cell transformation tensor.
virtual const tensorField & forwardT() const
Return face transformation tensor.
Cyclic patch for Arbitrary Mesh Interface (AMI)
virtual bool owner() const
Does this side own the patch?
word periodicPatchName_
Periodic patch name.
autoPtr< AMIPatchToPatchInterpolation > AMIPtr_
AMI interpolation class.
label periodicPatchID() const
Periodic patch ID (or -1)
Cyclic patch for periodic Arbitrary Mesh Interface (AMI)
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.
virtual bool write()
Write the output fields.
const Time & time() const noexcept
Return time registry.
const word & name() const noexcept
The patch name.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
int myProcNo() const noexcept
Return processor number.
splitCell * master() const
Definition: splitCell.H:113
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
engineTime & runTime
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
#define DebugInFunction
Report an information message using Foam::Info.
#define InfoInFunction
Report an information message using Foam::Info.
#define SeriousErrorInFunction
Report an error message using Foam::SeriousError.
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:203
Namespace for OpenFOAM.
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:536
List< label > labelList
A List of labels.
Definition: List.H:66
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
messageStream Info
Information stream (stdout output on master, null elsewhere)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:342
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Field< vector > vectorField
Specialisation of Field<T> for vector.
uint8_t direction
Definition: direction.H:56
List< bool > boolList
A List of bools.
Definition: List.H:64
Type gAverage(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
error FatalError
Type gMin(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Type gMax(const FieldField< Field, Type > &f)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333