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 -------------------------------------------------------------------------------
11 License
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 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(meshToMesh, 0);
42 }
43 
44 
45 const 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 
61 const Foam::Enum
62 <
64 >
66 {
67  { procMapMethod::pmAABB, "AABB" },
68  { procMapMethod::pmLOD, "LOD" },
69 };
70 
71 
72 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
73 
74 template<>
75 void Foam::meshToMesh::mapInternalSrcToTgt
76 (
78  const plusEqOp<sphericalTensor>& cop,
80  const bool secondOrder
81 ) const
82 {
83  mapSrcToTgt(field, cop, result.primitiveFieldRef());
84 }
85 
86 
87 template<>
88 void Foam::meshToMesh::mapInternalSrcToTgt
89 (
91  const minusEqOp<sphericalTensor>& cop,
93  const bool secondOrder
94 ) const
95 {
96  mapSrcToTgt(field, cop, result.primitiveFieldRef());
97 }
98 
99 
100 template<>
101 void Foam::meshToMesh::mapInternalSrcToTgt
102 (
104  const plusEqOp<symmTensor>& cop,
105  VolumeField<symmTensor>& result,
106  const bool secondOrder
107 ) const
108 {
109  mapSrcToTgt(field, cop, result.primitiveFieldRef());
110 }
111 
112 
113 template<>
114 void Foam::meshToMesh::mapInternalSrcToTgt
115 (
117  const minusEqOp<symmTensor>& cop,
118  VolumeField<symmTensor>& result,
119  const bool secondOrder
120 ) const
121 {
122  mapSrcToTgt(field, cop, result.primitiveFieldRef());
123 }
124 
125 
126 template<>
127 void Foam::meshToMesh::mapInternalSrcToTgt
128 (
129  const VolumeField<tensor>& field,
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 
139 template<>
140 void Foam::meshToMesh::mapInternalSrcToTgt
141 (
142  const VolumeField<tensor>& field,
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 
152 template<>
153 void 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 
165 template<>
166 void Foam::meshToMesh::mapInternalTgtToSrc
167 (
169  const minusEqOp<sphericalTensor>& cop,
171  const bool secondOrder
172 ) const
173 {
174  mapTgtToSrc(field, cop, result.primitiveFieldRef());
175 }
176 
177 
178 template<>
179 void Foam::meshToMesh::mapInternalTgtToSrc
180 (
182  const plusEqOp<symmTensor>& cop,
183  VolumeField<symmTensor>& result,
184  const bool secondOrder
185 ) const
186 {
187  mapTgtToSrc(field, cop, result.primitiveFieldRef());
188 }
189 
190 
191 template<>
192 void Foam::meshToMesh::mapInternalTgtToSrc
193 (
195  const minusEqOp<symmTensor>& cop,
196  VolumeField<symmTensor>& result,
197  const bool secondOrder
198 ) const
199 {
200  mapTgtToSrc(field, cop, result.primitiveFieldRef());
201 }
202 
203 
204 template<>
205 void Foam::meshToMesh::mapInternalTgtToSrc
206 (
207  const VolumeField<tensor>& field,
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 
217 template<>
218 void Foam::meshToMesh::mapInternalTgtToSrc
219 (
220  const VolumeField<tensor>& field,
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 
230 template<>
231 void Foam::meshToMesh::mapAndOpSrcToTgt
232 (
233  const AMIPatchToPatchInterpolation& AMI,
234  const Field<scalar>& srcField,
235  Field<scalar>& tgtField,
236  const plusEqOp<scalar>& cop
237 ) const
238 {}
239 
240 
241 template<>
242 void Foam::meshToMesh::mapAndOpSrcToTgt
243 (
244  const AMIPatchToPatchInterpolation& AMI,
245  const Field<vector>& srcField,
246  Field<vector>& tgtField,
247  const plusEqOp<vector>& cop
248 ) const
249 {}
250 
251 
252 template<>
253 void Foam::meshToMesh::mapAndOpSrcToTgt
254 (
255  const AMIPatchToPatchInterpolation& AMI,
256  const Field<sphericalTensor>& srcField,
257  Field<sphericalTensor>& tgtField,
259 ) const
260 {}
261 
262 
263 template<>
264 void Foam::meshToMesh::mapAndOpSrcToTgt
265 (
266  const AMIPatchToPatchInterpolation& AMI,
267  const Field<symmTensor>& srcField,
268  Field<symmTensor>& tgtField,
270 ) const
271 {}
272 
273 
274 template<>
275 void Foam::meshToMesh::mapAndOpSrcToTgt
276 (
277  const AMIPatchToPatchInterpolation& AMI,
278  const Field<tensor>& srcField,
279  Field<tensor>& tgtField,
280  const plusEqOp<tensor>& cop
281 ) const
282 {}
283 
284 
285 template<>
286 void Foam::meshToMesh::mapAndOpTgtToSrc
287 (
288  const AMIPatchToPatchInterpolation& AMI,
289  Field<scalar>& srcField,
290  const Field<scalar>& tgtField,
291  const plusEqOp<scalar>& cop
292 ) const
293 {}
294 
295 
296 template<>
297 void Foam::meshToMesh::mapAndOpTgtToSrc
298 (
299  const AMIPatchToPatchInterpolation& AMI,
300  Field<vector>& srcField,
301  const Field<vector>& tgtField,
302  const plusEqOp<vector>& cop
303 ) const
304 {}
305 
306 
307 template<>
308 void Foam::meshToMesh::mapAndOpTgtToSrc
309 (
310  const AMIPatchToPatchInterpolation& AMI,
311  Field<sphericalTensor>& srcField,
312  const Field<sphericalTensor>& tgtField,
314 ) const
315 {}
316 
317 
318 template<>
319 void Foam::meshToMesh::mapAndOpTgtToSrc
320 (
321  const AMIPatchToPatchInterpolation& AMI,
322  Field<symmTensor>& srcField,
323  const Field<symmTensor>& tgtField,
325 ) const
326 {}
327 
328 
329 template<>
330 void Foam::meshToMesh::mapAndOpTgtToSrc
331 (
332  const AMIPatchToPatchInterpolation& AMI,
333  Field<tensor>& srcField,
334  const Field<tensor>& tgtField,
335  const plusEqOp<tensor>& cop
336 ) const
337 {}
338 
339 
340 Foam::labelList Foam::meshToMesh::maskCells
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 
358  DynamicList<label> cells(src.size());
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 
377 void Foam::meshToMesh::normaliseWeights
378 (
379  const word& descriptor,
380  const labelListList& addr,
381  scalarListList& wght
382 ) const
383 {
384  if (returnReduceOr(wght.size()))
385  {
386  forAll(wght, celli)
387  {
388  scalarList& w = wght[celli];
389  scalar s = sum(w);
390 
391  forAll(w, i)
392  {
393  // note: normalise by s instead of cell volume since
394  // 1-to-1 methods duplicate contributions in parallel
395  w[i] /= s;
396  }
397  }
398  }
399 }
400 
401 
402 void Foam::meshToMesh::calcAddressing
403 (
404  const word& methodName,
405  const polyMesh& src,
406  const polyMesh& tgt
407 )
408 {
409  autoPtr<meshToMeshMethod> methodPtr
410  (
412  (
413  methodName,
414  src,
415  tgt
416  )
417  );
418 
419  methodPtr->calculate
420  (
421  srcToTgtCellAddr_,
422  srcToTgtCellWght_,
423  srcToTgtCellVec_,
424  tgtToSrcCellAddr_,
425  tgtToSrcCellWght_,
426  tgtToSrcCellVec_
427  );
428 
429  V_ = methodPtr->V();
430 
431  if (debug)
432  {
433  methodPtr->writeConnectivity(src, tgt, srcToTgtCellAddr_);
434  }
435 }
436 
437 
438 void Foam::meshToMesh::calculate(const word& methodName, const bool normalise)
439 {
440  Info<< "Creating mesh-to-mesh addressing for " << srcRegion_.name()
441  << " and " << tgtRegion_.name() << " regions using "
442  << methodName << endl;
443 
444  singleMeshProc_ = calcDistribution(srcRegion_, tgtRegion_);
445 
446  if (singleMeshProc_ == -1) // -> distributed()
447  {
448  // create global indexing for src and tgt meshes
449  globalIndex globalSrcCells(srcRegion_.nCells());
450  globalIndex globalTgtCells(tgtRegion_.nCells());
451 
452  // Create processor map of overlapping cells. This map gets
453  // (possibly remote) cells from the tgt mesh such that they (together)
454  // cover all of the src mesh
455  autoPtr<mapDistribute> mapPtr = calcProcMap(srcRegion_, tgtRegion_);
456  const mapDistribute& map = mapPtr();
457 
458  pointField newTgtPoints;
459  faceList newTgtFaces;
460  labelList newTgtFaceOwners;
461  labelList newTgtFaceNeighbours;
462  labelList newTgtCellIDs;
463 
464  distributeAndMergeCells
465  (
466  map,
467  tgtRegion_,
468  globalTgtCells,
469  newTgtPoints,
470  newTgtFaces,
471  newTgtFaceOwners,
472  newTgtFaceNeighbours,
473  newTgtCellIDs
474  );
475 
476 
477  // create a new target mesh
478  polyMesh newTgt
479  (
480  IOobject
481  (
482  "newTgt." + Foam::name(Pstream::myProcNo()),
483  tgtRegion_.time().timeName(),
484  tgtRegion_.time(),
487  ),
488  std::move(newTgtPoints),
489  std::move(newTgtFaces),
490  std::move(newTgtFaceOwners),
491  std::move(newTgtFaceNeighbours),
492  false // no parallel comms
493  );
494 
495  // create some dummy patch info
497  patches.set
498  (
499  0,
500  new polyPatch
501  (
502  "defaultFaces",
503  newTgt.nBoundaryFaces(),
504  newTgt.nInternalFaces(),
505  0,
506  newTgt.boundaryMesh(),
507  word::null
508  )
509  );
510 
511  newTgt.addPatches(patches);
512 
513  // force calculation of tet-base points used for point-in-cell
514  (void)newTgt.tetBasePtIs();
515 
516  // force construction of cell tree
517 // (void)newTgt.cellTree();
518 
519  if (debug)
520  {
521  Pout<< "Created newTgt mesh:" << nl
522  << " old cells = " << tgtRegion_.nCells()
523  << ", new cells = " << newTgt.nCells() << nl
524  << " old faces = " << tgtRegion_.nFaces()
525  << ", new faces = " << newTgt.nFaces() << endl;
526 
527  if (debug > 1)
528  {
529  Pout<< "Writing newTgt mesh: " << newTgt.name() << endl;
530  newTgt.write();
531  }
532  }
533 
534  calcAddressing(methodName, srcRegion_, newTgt);
535 
536  // Per source cell the target cell address in newTgt mesh
537  for (labelList& addressing : srcToTgtCellAddr_)
538  {
539  for (label& addr : addressing)
540  {
541  addr = newTgtCellIDs[addr];
542  }
543  }
544 
545  // Convert target addresses in newTgtMesh into global cell numbering
546  for (labelList& addressing : tgtToSrcCellAddr_)
547  {
548  globalSrcCells.inplaceToGlobal(addressing);
549  }
550 
551  // Set up as a reverse distribute
553  (
556  tgtRegion_.nCells(),
557  map.constructMap(),
558  false,
559  map.subMap(),
560  false,
561  tgtToSrcCellAddr_,
562  labelList(),
563  ListOps::appendEqOp<label>(),
564  flipOp(),
566  map.comm()
567  );
568 
569  // Set up as a reverse distribute
571  (
574  tgtRegion_.nCells(),
575  map.constructMap(),
576  false,
577  map.subMap(),
578  false,
579  tgtToSrcCellWght_,
580  scalarList(),
581  ListOps::appendEqOp<scalar>(),
582  flipOp(),
584  map.comm()
585  );
586 
587  // weights normalisation
588  if (normalise)
589  {
590  normaliseWeights
591  (
592  "source",
593  srcToTgtCellAddr_,
594  srcToTgtCellWght_
595  );
596 
597  normaliseWeights
598  (
599  "target",
600  tgtToSrcCellAddr_,
601  tgtToSrcCellWght_
602  );
603  }
604 
605  // cache maps and reset addresses
606  List<Map<label>> cMap;
607  srcMapPtr_.reset
608  (
609  new mapDistribute(globalSrcCells, tgtToSrcCellAddr_, cMap)
610  );
611  tgtMapPtr_.reset
612  (
613  new mapDistribute(globalTgtCells, srcToTgtCellAddr_, cMap)
614  );
615 
616  // collect volume intersection contributions
617  reduce(V_, sumOp<scalar>());
618  }
619  else
620  {
621  calcAddressing(methodName, srcRegion_, tgtRegion_);
622 
623  if (normalise)
624  {
625  normaliseWeights
626  (
627  "source",
628  srcToTgtCellAddr_,
629  srcToTgtCellWght_
630  );
631 
632  normaliseWeights
633  (
634  "target",
635  tgtToSrcCellAddr_,
636  tgtToSrcCellWght_
637  );
638  }
639  }
640 
641  Info<< " Overlap volume: " << V_ << endl;
642 }
643 
644 
646 (
647  const interpolationMethod method
648 )
649 {
650  switch (method)
651  {
652  case interpolationMethod::imDirect:
653  {
654  return nearestFaceAMI::typeName;
655  break;
656  }
657  case interpolationMethod::imMapNearest:
658  {
659  return nearestFaceAMI::typeName;
660  break;
661  }
662  case interpolationMethod::imCellVolumeWeight:
663  case interpolationMethod::imCorrectedCellVolumeWeight:
664  {
665  return faceAreaWeightAMI::typeName;
666  break;
667  }
668  default:
669  {
671  << "Unhandled enumeration " << interpolationMethodNames_[method]
672  << abort(FatalError);
673  }
674  }
675 
676  return nearestFaceAMI::typeName;
677 }
678 
679 
680 void Foam::meshToMesh::calculatePatchAMIs(const word& AMIMethodName)
681 {
682  if (!patchAMIs_.empty())
683  {
685  << "patch AMI already calculated"
686  << exit(FatalError);
687  }
688 
689  patchAMIs_.setSize(srcPatchID_.size());
690 
691  forAll(srcPatchID_, i)
692  {
693  label srcPatchi = srcPatchID_[i];
694  label tgtPatchi = tgtPatchID_[i];
695 
696  const polyPatch& srcPP = srcRegion_.boundaryMesh()[srcPatchi];
697  const polyPatch& tgtPP = tgtRegion_.boundaryMesh()[tgtPatchi];
698 
699  Info<< "Creating AMI between source patch " << srcPP.name()
700  << " and target patch " << tgtPP.name()
701  << " using " << AMIMethodName
702  << endl;
703 
704  Info<< incrIndent;
705 
706  patchAMIs_.set
707  (
708  i,
710  (
711  AMIMethodName,
712  false, // requireMatch
713  true, // flip target patch since patch normals are aligned
714  -1 // low weight correction
715  )
716  );
717 
718  patchAMIs_[i].calculate(srcPP, tgtPP);
719 
720  Info<< decrIndent;
721  }
722 }
723 
724 
725 void Foam::meshToMesh::constructNoCuttingPatches
726 (
727  const word& methodName,
728  const word& AMIMethodName,
729  const bool interpAllPatches
730 )
731 {
732  if (interpAllPatches)
733  {
734  const polyBoundaryMesh& srcBM = srcRegion_.boundaryMesh();
735  const polyBoundaryMesh& tgtBM = tgtRegion_.boundaryMesh();
736 
737  DynamicList<label> srcPatchID(srcBM.size());
738  DynamicList<label> tgtPatchID(tgtBM.size());
739  forAll(srcBM, patchi)
740  {
741  const polyPatch& pp = srcBM[patchi];
742 
743  // We want to map all the global patches, including constraint
744  // patches (since they might have mappable properties, e.g.
745  // jumpCyclic). We'll fix the value afterwards.
746  if (!isA<processorPolyPatch>(pp))
747  {
748  srcPatchID.append(pp.index());
749 
750  label tgtPatchi = tgtBM.findPatchID(pp.name());
751 
752  if (tgtPatchi != -1)
753  {
754  tgtPatchID.append(tgtPatchi);
755  }
756  else
757  {
759  << "Source patch " << pp.name()
760  << " not found in target mesh. "
761  << "Available target patches are " << tgtBM.names()
762  << exit(FatalError);
763  }
764  }
765  }
766 
767  srcPatchID_.transfer(srcPatchID);
768  tgtPatchID_.transfer(tgtPatchID);
769  }
770 
771  // calculate volume addressing and weights
772  calculate(methodName, true);
773 
774  // calculate patch addressing and weights
775  calculatePatchAMIs(AMIMethodName);
776 }
777 
778 
779 void Foam::meshToMesh::constructFromCuttingPatches
780 (
781  const word& methodName,
782  const word& AMIMethodName,
783  const HashTable<word>& patchMap,
784  const wordList& cuttingPatches,
785  const bool normalise
786 )
787 {
788  const polyBoundaryMesh& srcBm = srcRegion_.boundaryMesh();
789  const polyBoundaryMesh& tgtBm = tgtRegion_.boundaryMesh();
790 
791  // set IDs of cutting patches
792  cuttingPatches_.setSize(cuttingPatches.size());
793  forAll(cuttingPatches_, i)
794  {
795  const word& patchName = cuttingPatches[i];
796  label cuttingPatchi = srcBm.findPatchID(patchName);
797 
798  if (cuttingPatchi == -1)
799  {
801  << "Unable to find patch '" << patchName
802  << "' in mesh '" << srcRegion_.name() << "'. "
803  << " Available patches include:" << srcBm.names()
804  << exit(FatalError);
805  }
806 
807  cuttingPatches_[i] = cuttingPatchi;
808  }
809 
810  DynamicList<label> srcIDs(patchMap.size());
811  DynamicList<label> tgtIDs(patchMap.size());
812 
813  forAllConstIters(patchMap, iter)
814  {
815  const word& tgtPatchName = iter.key();
816  const word& srcPatchName = iter.val();
817 
818  const polyPatch& srcPatch = srcBm[srcPatchName];
819 
820  // We want to map all the global patches, including constraint
821  // patches (since they might have mappable properties, e.g.
822  // jumpCyclic). We'll fix the value afterwards.
823  if (!isA<processorPolyPatch>(srcPatch))
824  {
825  const polyPatch& tgtPatch = tgtBm[tgtPatchName];
826 
827  srcIDs.append(srcPatch.index());
828  tgtIDs.append(tgtPatch.index());
829  }
830  }
831 
832  srcPatchID_.transfer(srcIDs);
833  tgtPatchID_.transfer(tgtIDs);
834 
835  // calculate volume addressing and weights
836  calculate(methodName, normalise);
837 
838  // calculate patch addressing and weights
839  calculatePatchAMIs(AMIMethodName);
840 }
841 
842 
843 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
844 
845 Foam::meshToMesh::meshToMesh
846 (
847  const polyMesh& src,
848  const polyMesh& tgt,
849  const interpolationMethod method,
850  const procMapMethod mapMethod,
851  bool interpAllPatches
852 )
853 :
854  srcRegion_(src),
855  tgtRegion_(tgt),
856  procMapMethod_(mapMethod),
857  srcPatchID_(),
858  tgtPatchID_(),
859  patchAMIs_(),
860  cuttingPatches_(),
861  srcToTgtCellAddr_(),
862  tgtToSrcCellAddr_(),
863  srcToTgtCellWght_(),
864  tgtToSrcCellWght_(),
865  srcToTgtCellVec_(),
866  tgtToSrcCellVec_(),
867  V_(0.0),
868  singleMeshProc_(-1),
869  srcMapPtr_(nullptr),
870  tgtMapPtr_(nullptr)
871 {
872  constructNoCuttingPatches
873  (
876  interpAllPatches
877  );
878 }
879 
880 
881 Foam::meshToMesh::meshToMesh
882 (
883  const polyMesh& src,
884  const polyMesh& tgt,
885  const word& methodName,
886  const word& AMIMethodName,
887  const procMapMethod mapMethod,
888  bool interpAllPatches
889 )
890 :
891  srcRegion_(src),
892  tgtRegion_(tgt),
893  procMapMethod_(mapMethod),
894  srcPatchID_(),
895  tgtPatchID_(),
896  patchAMIs_(),
897  cuttingPatches_(),
898  srcToTgtCellAddr_(),
899  tgtToSrcCellAddr_(),
900  srcToTgtCellWght_(),
901  tgtToSrcCellWght_(),
902  srcToTgtCellVec_(),
903  tgtToSrcCellVec_(),
904  V_(0.0),
905  singleMeshProc_(-1),
906  srcMapPtr_(nullptr),
907  tgtMapPtr_(nullptr)
908 {
909  constructNoCuttingPatches(methodName, AMIMethodName, interpAllPatches);
910 }
911 
912 
913 Foam::meshToMesh::meshToMesh
914 (
915  const polyMesh& src,
916  const polyMesh& tgt,
917  const interpolationMethod method,
918  const HashTable<word>& patchMap,
919  const wordList& cuttingPatches,
920  const procMapMethod mapMethod,
921  const bool normalise
922 )
923 :
924  srcRegion_(src),
925  tgtRegion_(tgt),
926  procMapMethod_(mapMethod),
927  srcPatchID_(),
928  tgtPatchID_(),
929  patchAMIs_(),
930  cuttingPatches_(),
931  srcToTgtCellAddr_(),
932  tgtToSrcCellAddr_(),
933  srcToTgtCellWght_(),
934  tgtToSrcCellWght_(),
935  V_(0.0),
936  singleMeshProc_(-1),
937  srcMapPtr_(nullptr),
938  tgtMapPtr_(nullptr)
939 {
940  constructFromCuttingPatches
941  (
943  interpolationMethodAMI(method),
944  patchMap,
945  cuttingPatches,
946  normalise
947  );
948 }
949 
950 
951 Foam::meshToMesh::meshToMesh
952 (
953  const polyMesh& src,
954  const polyMesh& tgt,
955  const word& methodName, // internal mapping
956  const word& AMIMethodName, // boundary mapping
957  const HashTable<word>& patchMap,
958  const wordList& cuttingPatches,
959  const procMapMethod mapMethod,
960  const bool normalise
961 )
962 :
963  srcRegion_(src),
964  tgtRegion_(tgt),
965  procMapMethod_(mapMethod),
966  srcPatchID_(),
967  tgtPatchID_(),
968  patchAMIs_(),
969  cuttingPatches_(),
970  srcToTgtCellAddr_(),
971  tgtToSrcCellAddr_(),
972  srcToTgtCellWght_(),
973  tgtToSrcCellWght_(),
974  V_(0.0),
975  singleMeshProc_(-1),
976  srcMapPtr_(nullptr),
977  tgtMapPtr_(nullptr)
978 {
979  constructFromCuttingPatches
980  (
981  methodName,
982  AMIMethodName,
983  patchMap,
984  cuttingPatches,
985  normalise
986  );
987 }
988 
989 
990 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
991 
993 {}
994 
995 
996 // ************************************************************************* //
List< scalar > scalarList
List of scalar.
Definition: scalarList.H:32
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1)
rDeltaTY field()
List< cell > cellList
List of cell.
Definition: cellListFwd.H:39
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:517
static const Enum< interpolationMethod > interpolationMethodNames_
Definition: meshToMesh.H:77
static const Enum< procMapMethod > procMapMethodNames_
Definition: meshToMesh.H:88
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:1229
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1074
static const List< T > & null()
Return a null List.
Definition: ListI.H:130
virtual const fileName & name() const override
Get the name of the output serial stream. (eg, the name of the Fstream file name) ...
Definition: OSstream.H:128
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
procMapMethod
Enumeration specifying processor parallel map construction method.
Definition: meshToMesh.H:82
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
PtrList< polyPatch > polyPatchList
Store lists of polyPatch as a PtrList.
Definition: polyPatch.H:56
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const cellShapeList & cells
const pointField & points
List< scalarList > scalarListList
List of scalarList.
Definition: scalarList.H:35
A class for handling words, derived from Foam::string.
Definition: word.H:63
static const word null
An empty word.
Definition: word.H:84
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
A HashTable similar to std::unordered_map.
Definition: HashTable.H:108
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
errorManip< error > abort(error &err)
Definition: errorManip.H:139
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: PtrList.H:159
static word interpolationMethodAMI(const interpolationMethod method)
Conversion between mesh and patch interpolation methods.
Definition: meshToMesh.C:639
int debug
Static debugging option.
defineTypeNameAndDebug(combustionModel, 0)
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:511
List< word > wordList
List of word.
Definition: fileName.H:59
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
virtual ~meshToMesh()
Destructor.
Definition: meshToMesh.C:985
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
const polyBoundaryMesh & patches
Nothing to be read.
Automatically write from objectRegistry::writeObject()
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
"nonBlocking" : (MPI_Isend, MPI_Irecv)
messageStream Info
Information stream (stdout output on master, null elsewhere)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
static autoPtr< AMIInterpolation > New(const word &modelName, const dictionary &dict, const bool reverseTarget=false)
Selector for dictionary.
interpolationMethod
Enumeration specifying interpolation method.
Definition: meshToMesh.H:69
List< label > labelList
A List of labels.
Definition: List.H:62
void mapSrcToTgt(const UList< Type > &srcFld, const CombineOp &cop, List< Type > &result) const
Map field from src to tgt mesh with defined operation.
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))
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:502
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
static autoPtr< meshToMeshMethod > New(const word &methodName, const polyMesh &src, const polyMesh &tgt)
Selector.
static void distribute(const UPstream::commsTypes commsType, const List< labelPair > &schedule, const label constructSize, const labelListList &subMap, const bool subHasFlip, const labelListList &constructMap, const bool constructHasFlip, List< T > &field, const T &nullValue, const CombineOp &cop, const NegateOp &negOp, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Distribute combine data with specified combine operation and negate operator (for flips)...
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28