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