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