meshToMesh.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) 2012-2017 OpenFOAM Foundation
9 Copyright (C) 2015-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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 "meshToMesh.H"
30#include "Time.H"
31#include "globalIndex.H"
32#include "meshToMeshMethod.H"
33#include "nearestFaceAMI.H"
34#include "processorPolyPatch.H"
35#include "faceAreaWeightAMI.H"
36
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
40{
42}
43
44
45const Foam::Enum
46<
48>
50({
51 { interpolationMethod::imDirect, "direct" },
52 { interpolationMethod::imMapNearest, "mapNearest" },
53 { interpolationMethod::imCellVolumeWeight, "cellVolumeWeight" },
54 {
55 interpolationMethod::imCorrectedCellVolumeWeight,
56 "correctedCellVolumeWeight"
57 },
58});
59
60
61const Foam::Enum
62<
64>
66{
67 { procMapMethod::pmAABB, "AABB" },
68 { procMapMethod::pmLOD, "LOD" },
69};
70
71
72// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
73
74template<>
75void Foam::meshToMesh::mapInternalSrcToTgt
76(
80 const bool secondOrder
81) const
82{
83 mapSrcToTgt(field, cop, result.primitiveFieldRef());
84}
85
86
87template<>
88void Foam::meshToMesh::mapInternalSrcToTgt
89(
93 const bool secondOrder
94) const
95{
96 mapSrcToTgt(field, cop, result.primitiveFieldRef());
97}
98
99
100template<>
101void Foam::meshToMesh::mapInternalSrcToTgt
102(
104 const plusEqOp<symmTensor>& cop,
106 const bool secondOrder
107) const
108{
109 mapSrcToTgt(field, cop, result.primitiveFieldRef());
110}
111
112
113template<>
114void Foam::meshToMesh::mapInternalSrcToTgt
115(
117 const minusEqOp<symmTensor>& cop,
119 const bool secondOrder
120) const
121{
122 mapSrcToTgt(field, cop, result.primitiveFieldRef());
123}
124
125
126template<>
127void Foam::meshToMesh::mapInternalSrcToTgt
128(
130 const plusEqOp<tensor>& cop,
131 VolumeField<tensor>& result,
132 const bool secondOrder
133) const
134{
135 mapSrcToTgt(field, cop, result.primitiveFieldRef());
136}
137
138
139template<>
140void Foam::meshToMesh::mapInternalSrcToTgt
141(
143 const minusEqOp<tensor>& cop,
144 VolumeField<tensor>& result,
145 const bool secondOrder
146) const
147{
148 mapSrcToTgt(field, cop, result.primitiveFieldRef());
149}
150
151
152template<>
153void Foam::meshToMesh::mapInternalTgtToSrc
154(
156 const plusEqOp<sphericalTensor>& cop,
158 const bool secondOrder
159) const
160{
161 mapTgtToSrc(field, cop, result.primitiveFieldRef());
162}
163
164
165template<>
166void Foam::meshToMesh::mapInternalTgtToSrc
167(
171 const bool secondOrder
172) const
173{
174 mapTgtToSrc(field, cop, result.primitiveFieldRef());
175}
176
177
178template<>
179void Foam::meshToMesh::mapInternalTgtToSrc
180(
182 const plusEqOp<symmTensor>& cop,
184 const bool secondOrder
185) const
186{
187 mapTgtToSrc(field, cop, result.primitiveFieldRef());
188}
189
190
191template<>
192void Foam::meshToMesh::mapInternalTgtToSrc
193(
195 const minusEqOp<symmTensor>& cop,
197 const bool secondOrder
198) const
199{
200 mapTgtToSrc(field, cop, result.primitiveFieldRef());
201}
202
203
204template<>
205void Foam::meshToMesh::mapInternalTgtToSrc
206(
208 const plusEqOp<tensor>& cop,
209 VolumeField<tensor>& result,
210 const bool secondOrder
211) const
212{
213 mapTgtToSrc(field, cop, result.primitiveFieldRef());
214}
215
216
217template<>
218void Foam::meshToMesh::mapInternalTgtToSrc
219(
221 const minusEqOp<tensor>& cop,
222 VolumeField<tensor>& result,
223 const bool secondOrder
224) const
225{
226 mapTgtToSrc(field, cop, result.primitiveFieldRef());
227}
228
229
230template<>
231void Foam::meshToMesh::mapAndOpSrcToTgt
232(
234 const Field<scalar>& srcField,
235 Field<scalar>& tgtField,
236 const plusEqOp<scalar>& cop
237) const
238{}
239
240
241template<>
242void Foam::meshToMesh::mapAndOpSrcToTgt
243(
245 const Field<vector>& srcField,
246 Field<vector>& tgtField,
247 const plusEqOp<vector>& cop
248) const
249{}
250
251
252template<>
253void Foam::meshToMesh::mapAndOpSrcToTgt
254(
256 const Field<sphericalTensor>& srcField,
257 Field<sphericalTensor>& tgtField,
259) const
260{}
261
262
263template<>
264void Foam::meshToMesh::mapAndOpSrcToTgt
265(
267 const Field<symmTensor>& srcField,
268 Field<symmTensor>& tgtField,
269 const plusEqOp<symmTensor>& cop
270) const
271{}
272
273
274template<>
275void Foam::meshToMesh::mapAndOpSrcToTgt
276(
278 const Field<tensor>& srcField,
279 Field<tensor>& tgtField,
280 const plusEqOp<tensor>& cop
281) const
282{}
283
284
285template<>
286void Foam::meshToMesh::mapAndOpTgtToSrc
287(
289 Field<scalar>& srcField,
290 const Field<scalar>& tgtField,
291 const plusEqOp<scalar>& cop
292) const
293{}
294
295
296template<>
297void Foam::meshToMesh::mapAndOpTgtToSrc
298(
300 Field<vector>& srcField,
301 const Field<vector>& tgtField,
302 const plusEqOp<vector>& cop
303) const
304{}
305
306
307template<>
308void Foam::meshToMesh::mapAndOpTgtToSrc
309(
311 Field<sphericalTensor>& srcField,
312 const Field<sphericalTensor>& tgtField,
314) const
315{}
316
317
318template<>
319void Foam::meshToMesh::mapAndOpTgtToSrc
320(
322 Field<symmTensor>& srcField,
323 const Field<symmTensor>& tgtField,
324 const plusEqOp<symmTensor>& cop
325) const
326{}
327
328
329template<>
330void Foam::meshToMesh::mapAndOpTgtToSrc
331(
333 Field<tensor>& srcField,
334 const Field<tensor>& tgtField,
335 const plusEqOp<tensor>& cop
336) const
337{}
338
339
341(
342 const polyMesh& src,
343 const polyMesh& tgt
344) const
345{
346 boundBox intersectBb
347 (
348 max(src.bounds().min(), tgt.bounds().min()),
349 min(src.bounds().max(), tgt.bounds().max())
350 );
351
352 intersectBb.inflate(0.01);
353
354 const cellList& srcCells = src.cells();
355 const faceList& srcFaces = src.faces();
356 const pointField& srcPts = src.points();
357
359 forAll(srcCells, srcI)
360 {
361 boundBox cellBb(srcCells[srcI].points(srcFaces, srcPts), false);
362 if (intersectBb.overlaps(cellBb))
363 {
364 cells.append(srcI);
365 }
366 }
367
368 if (debug)
369 {
370 Pout<< "participating source mesh cells: " << cells.size() << endl;
371 }
372
373 return cells;
374}
375
376
377void Foam::meshToMesh::normaliseWeights
378(
379 const word& descriptor,
380 const labelListList& addr,
381 scalarListList& wght
382) const
383{
384 const label nCell = returnReduce(wght.size(), sumOp<label>());
385
386 if (nCell > 0)
387 {
388 forAll(wght, celli)
389 {
390 scalarList& w = wght[celli];
391 scalar s = sum(w);
392
393 forAll(w, i)
394 {
395 // note: normalise by s instead of cell volume since
396 // 1-to-1 methods duplicate contributions in parallel
397 w[i] /= s;
398 }
399 }
400 }
401}
402
403
404void Foam::meshToMesh::calcAddressing
405(
406 const word& methodName,
407 const polyMesh& src,
408 const polyMesh& tgt
409)
410{
411 autoPtr<meshToMeshMethod> methodPtr
412 (
414 (
415 methodName,
416 src,
417 tgt
418 )
419 );
420
421 methodPtr->calculate
422 (
423 srcToTgtCellAddr_,
424 srcToTgtCellWght_,
425 srcToTgtCellVec_,
426 tgtToSrcCellAddr_,
427 tgtToSrcCellWght_,
428 tgtToSrcCellVec_
429 );
430
431 V_ = methodPtr->V();
432
433 if (debug)
434 {
435 methodPtr->writeConnectivity(src, tgt, srcToTgtCellAddr_);
436 }
437}
438
439
440void Foam::meshToMesh::calculate(const word& methodName, const bool normalise)
441{
442 Info<< "Creating mesh-to-mesh addressing for " << srcRegion_.name()
443 << " and " << tgtRegion_.name() << " regions using "
444 << methodName << endl;
445
446 singleMeshProc_ = calcDistribution(srcRegion_, tgtRegion_);
447
448 if (singleMeshProc_ == -1)
449 {
450 // create global indexing for src and tgt meshes
451 globalIndex globalSrcCells(srcRegion_.nCells());
452 globalIndex globalTgtCells(tgtRegion_.nCells());
453
454 // Create processor map of overlapping cells. This map gets
455 // (possibly remote) cells from the tgt mesh such that they (together)
456 // cover all of the src mesh
457 autoPtr<mapDistribute> mapPtr = calcProcMap(srcRegion_, tgtRegion_);
458 const mapDistribute& map = mapPtr();
459
460 pointField newTgtPoints;
461 faceList newTgtFaces;
462 labelList newTgtFaceOwners;
463 labelList newTgtFaceNeighbours;
464 labelList newTgtCellIDs;
465
466 distributeAndMergeCells
467 (
468 map,
469 tgtRegion_,
470 globalTgtCells,
471 newTgtPoints,
472 newTgtFaces,
473 newTgtFaceOwners,
474 newTgtFaceNeighbours,
475 newTgtCellIDs
476 );
477
478
479 // create a new target mesh
480 polyMesh newTgt
481 (
482 IOobject
483 (
484 "newTgt." + Foam::name(Pstream::myProcNo()),
485 tgtRegion_.time().timeName(),
486 tgtRegion_.time(),
489 ),
490 std::move(newTgtPoints),
491 std::move(newTgtFaces),
492 std::move(newTgtFaceOwners),
493 std::move(newTgtFaceNeighbours),
494 false // no parallel comms
495 );
496
497 // create some dummy patch info
498 List<polyPatch*> patches(1);
499 patches[0] = new polyPatch
500 (
501 "defaultFaces",
502 newTgt.nBoundaryFaces(),
503 newTgt.nInternalFaces(),
504 0,
505 newTgt.boundaryMesh(),
507 );
508
509 newTgt.addPatches(patches);
510
511 // force calculation of tet-base points used for point-in-cell
512 (void)newTgt.tetBasePtIs();
513
514 // force construction of cell tree
515// (void)newTgt.cellTree();
516
517 if (debug)
518 {
519 Pout<< "Created newTgt mesh:" << nl
520 << " old cells = " << tgtRegion_.nCells()
521 << ", new cells = " << newTgt.nCells() << nl
522 << " old faces = " << tgtRegion_.nFaces()
523 << ", new faces = " << newTgt.nFaces() << endl;
524
525 if (debug > 1)
526 {
527 Pout<< "Writing newTgt mesh: " << newTgt.name() << endl;
528 newTgt.write();
529 }
530 }
531
532 calcAddressing(methodName, srcRegion_, newTgt);
533
534 // Per source cell the target cell address in newTgt mesh
535 for (labelList& addressing : srcToTgtCellAddr_)
536 {
537 for (label& addr : addressing)
538 {
539 addr = newTgtCellIDs[addr];
540 }
541 }
542
543 // Convert target addresses in newTgtMesh into global cell numbering
544 for (labelList& addressing : tgtToSrcCellAddr_)
545 {
546 globalSrcCells.inplaceToGlobal(addressing);
547 }
548
549 // Set up as a reverse distribute
550 mapDistributeBase::distribute
551 (
554 tgtRegion_.nCells(),
555 map.constructMap(),
556 false,
557 map.subMap(),
558 false,
559 tgtToSrcCellAddr_,
560 labelList(),
561 ListOps::appendEqOp<label>(),
562 flipOp(),
564 map.comm()
565 );
566
567 // Set up as a reverse distribute
568 mapDistributeBase::distribute
569 (
572 tgtRegion_.nCells(),
573 map.constructMap(),
574 false,
575 map.subMap(),
576 false,
577 tgtToSrcCellWght_,
578 scalarList(),
579 ListOps::appendEqOp<scalar>(),
580 flipOp(),
582 map.comm()
583 );
584
585 // weights normalisation
586 if (normalise)
587 {
588 normaliseWeights
589 (
590 "source",
591 srcToTgtCellAddr_,
592 srcToTgtCellWght_
593 );
594
595 normaliseWeights
596 (
597 "target",
598 tgtToSrcCellAddr_,
599 tgtToSrcCellWght_
600 );
601 }
602
603 // cache maps and reset addresses
604 List<Map<label>> cMap;
605 srcMapPtr_.reset
606 (
607 new mapDistribute(globalSrcCells, tgtToSrcCellAddr_, cMap)
608 );
609 tgtMapPtr_.reset
610 (
611 new mapDistribute(globalTgtCells, srcToTgtCellAddr_, cMap)
612 );
613
614 // collect volume intersection contributions
615 reduce(V_, sumOp<scalar>());
616 }
617 else
618 {
619 calcAddressing(methodName, srcRegion_, tgtRegion_);
620
621 if (normalise)
622 {
623 normaliseWeights
624 (
625 "source",
626 srcToTgtCellAddr_,
627 srcToTgtCellWght_
628 );
629
630 normaliseWeights
631 (
632 "target",
633 tgtToSrcCellAddr_,
634 tgtToSrcCellWght_
635 );
636 }
637 }
638
639 Info<< " Overlap volume: " << V_ << endl;
640}
641
642
644(
645 const interpolationMethod method
646)
647{
648 switch (method)
649 {
650 case interpolationMethod::imDirect:
651 {
653 break;
654 }
655 case interpolationMethod::imMapNearest:
656 {
658 break;
659 }
660 case interpolationMethod::imCellVolumeWeight:
661 case interpolationMethod::imCorrectedCellVolumeWeight:
662 {
664 break;
665 }
666 default:
667 {
669 << "Unhandled enumeration " << interpolationMethodNames_[method]
670 << abort(FatalError);
671 }
672 }
673
675}
676
677
678void Foam::meshToMesh::calculatePatchAMIs(const word& AMIMethodName)
679{
680 if (!patchAMIs_.empty())
681 {
683 << "patch AMI already calculated"
684 << exit(FatalError);
685 }
686
687 patchAMIs_.setSize(srcPatchID_.size());
688
689 forAll(srcPatchID_, i)
690 {
691 label srcPatchi = srcPatchID_[i];
692 label tgtPatchi = tgtPatchID_[i];
693
694 const polyPatch& srcPP = srcRegion_.boundaryMesh()[srcPatchi];
695 const polyPatch& tgtPP = tgtRegion_.boundaryMesh()[tgtPatchi];
696
697 Info<< "Creating AMI between source patch " << srcPP.name()
698 << " and target patch " << tgtPP.name()
699 << " using " << AMIMethodName
700 << endl;
701
703
704 patchAMIs_.set
705 (
706 i,
708 (
709 AMIMethodName,
710 false, // requireMatch
711 true, // flip target patch since patch normals are aligned
712 -1 // low weight correction
713 )
714 );
715
716 patchAMIs_[i].calculate(srcPP, tgtPP);
717
719 }
720}
721
722
723void Foam::meshToMesh::constructNoCuttingPatches
724(
725 const word& methodName,
726 const word& AMIMethodName,
727 const bool interpAllPatches
728)
729{
730 if (interpAllPatches)
731 {
732 const polyBoundaryMesh& srcBM = srcRegion_.boundaryMesh();
733 const polyBoundaryMesh& tgtBM = tgtRegion_.boundaryMesh();
734
735 DynamicList<label> srcPatchID(srcBM.size());
736 DynamicList<label> tgtPatchID(tgtBM.size());
737 forAll(srcBM, patchi)
738 {
739 const polyPatch& pp = srcBM[patchi];
740
741 // We want to map all the global patches, including constraint
742 // patches (since they might have mappable properties, e.g.
743 // jumpCyclic). We'll fix the value afterwards.
744 if (!isA<processorPolyPatch>(pp))
745 {
746 srcPatchID.append(pp.index());
747
748 label tgtPatchi = tgtBM.findPatchID(pp.name());
749
750 if (tgtPatchi != -1)
751 {
752 tgtPatchID.append(tgtPatchi);
753 }
754 else
755 {
757 << "Source patch " << pp.name()
758 << " not found in target mesh. "
759 << "Available target patches are " << tgtBM.names()
760 << exit(FatalError);
761 }
762 }
763 }
764
765 srcPatchID_.transfer(srcPatchID);
766 tgtPatchID_.transfer(tgtPatchID);
767 }
768
769 // calculate volume addressing and weights
770 calculate(methodName, true);
771
772 // calculate patch addressing and weights
773 calculatePatchAMIs(AMIMethodName);
774}
775
776
777void Foam::meshToMesh::constructFromCuttingPatches
778(
779 const word& methodName,
780 const word& AMIMethodName,
781 const HashTable<word>& patchMap,
782 const wordList& cuttingPatches,
783 const bool normalise
784)
785{
786 const polyBoundaryMesh& srcBm = srcRegion_.boundaryMesh();
787 const polyBoundaryMesh& tgtBm = tgtRegion_.boundaryMesh();
788
789 // set IDs of cutting patches
790 cuttingPatches_.setSize(cuttingPatches.size());
791 forAll(cuttingPatches_, i)
792 {
793 const word& patchName = cuttingPatches[i];
794 label cuttingPatchi = srcBm.findPatchID(patchName);
795
796 if (cuttingPatchi == -1)
797 {
799 << "Unable to find patch '" << patchName
800 << "' in mesh '" << srcRegion_.name() << "'. "
801 << " Available patches include:" << srcBm.names()
802 << exit(FatalError);
803 }
804
805 cuttingPatches_[i] = cuttingPatchi;
806 }
807
808 DynamicList<label> srcIDs(patchMap.size());
809 DynamicList<label> tgtIDs(patchMap.size());
810
811 forAllConstIters(patchMap, iter)
812 {
813 const word& tgtPatchName = iter.key();
814 const word& srcPatchName = iter.val();
815
816 const polyPatch& srcPatch = srcBm[srcPatchName];
817
818 // We want to map all the global patches, including constraint
819 // patches (since they might have mappable properties, e.g.
820 // jumpCyclic). We'll fix the value afterwards.
821 if (!isA<processorPolyPatch>(srcPatch))
822 {
823 const polyPatch& tgtPatch = tgtBm[tgtPatchName];
824
825 srcIDs.append(srcPatch.index());
826 tgtIDs.append(tgtPatch.index());
827 }
828 }
829
830 srcPatchID_.transfer(srcIDs);
831 tgtPatchID_.transfer(tgtIDs);
832
833 // calculate volume addressing and weights
834 calculate(methodName, normalise);
835
836 // calculate patch addressing and weights
837 calculatePatchAMIs(AMIMethodName);
838}
839
840
841// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
842
844(
845 const polyMesh& src,
846 const polyMesh& tgt,
847 const interpolationMethod& method,
848 const procMapMethod& mapMethod,
849 bool interpAllPatches
850)
851:
852 srcRegion_(src),
853 tgtRegion_(tgt),
854 procMapMethod_(mapMethod),
855 srcPatchID_(),
856 tgtPatchID_(),
857 patchAMIs_(),
858 cuttingPatches_(),
859 srcToTgtCellAddr_(),
860 tgtToSrcCellAddr_(),
861 srcToTgtCellWght_(),
862 tgtToSrcCellWght_(),
863 srcToTgtCellVec_(),
864 tgtToSrcCellVec_(),
865 V_(0.0),
866 singleMeshProc_(-1),
867 srcMapPtr_(nullptr),
868 tgtMapPtr_(nullptr)
869{
870 constructNoCuttingPatches
871 (
874 interpAllPatches
875 );
876}
877
878
880(
881 const polyMesh& src,
882 const polyMesh& tgt,
883 const word& methodName,
884 const word& AMIMethodName,
885 const procMapMethod& mapMethod,
886 bool interpAllPatches
887)
888:
889 srcRegion_(src),
890 tgtRegion_(tgt),
891 procMapMethod_(mapMethod),
892 srcPatchID_(),
893 tgtPatchID_(),
894 patchAMIs_(),
895 cuttingPatches_(),
896 srcToTgtCellAddr_(),
897 tgtToSrcCellAddr_(),
898 srcToTgtCellWght_(),
899 tgtToSrcCellWght_(),
900 srcToTgtCellVec_(),
901 tgtToSrcCellVec_(),
902 V_(0.0),
903 singleMeshProc_(-1),
904 srcMapPtr_(nullptr),
905 tgtMapPtr_(nullptr)
906{
907 constructNoCuttingPatches(methodName, AMIMethodName, interpAllPatches);
908}
909
910
912(
913 const polyMesh& src,
914 const polyMesh& tgt,
915 const interpolationMethod& method,
916 const HashTable<word>& patchMap,
917 const wordList& cuttingPatches,
918 const procMapMethod& mapMethod,
919 const bool normalise
920)
921:
922 srcRegion_(src),
923 tgtRegion_(tgt),
924 procMapMethod_(mapMethod),
925 srcPatchID_(),
926 tgtPatchID_(),
927 patchAMIs_(),
928 cuttingPatches_(),
929 srcToTgtCellAddr_(),
930 tgtToSrcCellAddr_(),
931 srcToTgtCellWght_(),
932 tgtToSrcCellWght_(),
933 V_(0.0),
934 singleMeshProc_(-1),
935 srcMapPtr_(nullptr),
936 tgtMapPtr_(nullptr)
937{
938 constructFromCuttingPatches
939 (
942 patchMap,
943 cuttingPatches,
944 normalise
945 );
946}
947
948
950(
951 const polyMesh& src,
952 const polyMesh& tgt,
953 const word& methodName, // internal mapping
954 const word& AMIMethodName, // boundary mapping
955 const HashTable<word>& patchMap,
956 const wordList& cuttingPatches,
957 const procMapMethod& mapMethod,
958 const bool normalise
959)
960:
961 srcRegion_(src),
962 tgtRegion_(tgt),
963 procMapMethod_(mapMethod),
964 srcPatchID_(),
965 tgtPatchID_(),
966 patchAMIs_(),
967 cuttingPatches_(),
968 srcToTgtCellAddr_(),
969 tgtToSrcCellAddr_(),
970 srcToTgtCellWght_(),
971 tgtToSrcCellWght_(),
972 V_(0.0),
973 singleMeshProc_(-1),
974 srcMapPtr_(nullptr),
975 tgtMapPtr_(nullptr)
976{
977 constructFromCuttingPatches
978 (
979 methodName,
980 AMIMethodName,
981 patchMap,
982 cuttingPatches,
983 normalise
984 );
985}
986
987
988// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
989
991{}
992
993
994// ************************************************************************* //
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
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: Enum.H:61
Generic templated field type.
Definition: Field.H:82
Generic GeometricField class.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
A HashTable similar to std::unordered_map.
Definition: HashTable.H:123
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:52
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
static const List< T > & null()
Return a null List.
Definition: ListI.H:109
virtual const fileName & name() const
Get the name of the stream.
Definition: OSstream.H:107
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
@ nonBlocking
"nonBlocking"
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:556
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:64
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
labelList maskCells() const
Return src cell IDs for the overlap region.
Class to calculate the cell-addressing between two overlapping meshes.
Definition: meshToMesh.H:65
virtual ~meshToMesh()
Destructor.
Definition: meshToMesh.C:990
interpolationMethod
Enumeration specifying interpolation method.
Definition: meshToMesh.H:72
procMapMethod
Enumeration specifying processor parallel map construction method.
Definition: meshToMesh.H:83
static word interpolationMethodAMI(const interpolationMethod method)
Conversion between mesh and patch interpolation methods.
Definition: meshToMesh.C:644
void mapSrcToTgt(const UList< Type > &srcFld, const CombineOp &cop, List< Type > &result) const
Map field from src to tgt mesh with defined operation.
static const Enum< procMapMethod > procMapMethodNames_
Definition: meshToMesh.H:88
static const Enum< interpolationMethod > interpolationMethodNames_
Definition: meshToMesh.H:79
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1108
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:462
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1083
const cellList & cells() const
int myProcNo() const noexcept
Return processor number.
A class for handling words, derived from Foam::string.
Definition: word.H:68
static const word null
An empty word.
Definition: word.H:80
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
rDeltaTY field()
const polyBoundaryMesh & patches
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const pointField & points
const cellShapeList & cells
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))
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:63
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
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
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:64
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:349
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
errorManip< error > abort(error &err)
Definition: errorManip.H:144
error FatalError
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.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
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
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:356
List< scalarList > scalarListList
A List of scalarList.
Definition: scalarList.H:66
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278
static const char *const typeName
The type name used in ensight case files.