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  mapTgtToSrc(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  mapTgtToSrc(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  mapTgtToSrc(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  mapTgtToSrc(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  mapTgtToSrc(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  mapTgtToSrc(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  labelList(),
560  ListOps::appendEqOp<label>(),
561  flipOp(),
563  map.comm()
564  );
565 
566  // Set up as a reverse distribute
568  (
570  List<labelPair>(),
571  tgtRegion_.nCells(),
572  map.constructMap(),
573  false,
574  map.subMap(),
575  false,
576  tgtToSrcCellWght_,
577  scalarList(),
578  ListOps::appendEqOp<scalar>(),
579  flipOp(),
581  map.comm()
582  );
583 
584  // weights normalisation
585  if (normalise)
586  {
587  normaliseWeights
588  (
589  "source",
590  srcToTgtCellAddr_,
591  srcToTgtCellWght_
592  );
593 
594  normaliseWeights
595  (
596  "target",
597  tgtToSrcCellAddr_,
598  tgtToSrcCellWght_
599  );
600  }
601 
602  // cache maps and reset addresses
603  List<Map<label>> cMap;
604  srcMapPtr_.reset
605  (
606  new mapDistribute(globalSrcCells, tgtToSrcCellAddr_, cMap)
607  );
608  tgtMapPtr_.reset
609  (
610  new mapDistribute(globalTgtCells, srcToTgtCellAddr_, cMap)
611  );
612 
613  // collect volume intersection contributions
614  reduce(V_, sumOp<scalar>());
615  }
616  else
617  {
618  calcAddressing(methodName, srcRegion_, tgtRegion_);
619 
620  if (normalise)
621  {
622  normaliseWeights
623  (
624  "source",
625  srcToTgtCellAddr_,
626  srcToTgtCellWght_
627  );
628 
629  normaliseWeights
630  (
631  "target",
632  tgtToSrcCellAddr_,
633  tgtToSrcCellWght_
634  );
635  }
636  }
637 
638  Info<< " Overlap volume: " << V_ << endl;
639 }
640 
641 
643 (
644  const interpolationMethod method
645 )
646 {
647  switch (method)
648  {
649  case interpolationMethod::imDirect:
650  {
651  return nearestFaceAMI::typeName;
652  break;
653  }
654  case interpolationMethod::imMapNearest:
655  {
656  return nearestFaceAMI::typeName;
657  break;
658  }
659  case interpolationMethod::imCellVolumeWeight:
660  case interpolationMethod::imCorrectedCellVolumeWeight:
661  {
662  return faceAreaWeightAMI::typeName;
663  break;
664  }
665  default:
666  {
668  << "Unhandled enumeration " << interpolationMethodNames_[method]
669  << abort(FatalError);
670  }
671  }
672 
673  return nearestFaceAMI::typeName;
674 }
675 
676 
677 void Foam::meshToMesh::calculatePatchAMIs(const word& AMIMethodName)
678 {
679  if (!patchAMIs_.empty())
680  {
682  << "patch AMI already calculated"
683  << exit(FatalError);
684  }
685 
686  patchAMIs_.setSize(srcPatchID_.size());
687 
688  forAll(srcPatchID_, i)
689  {
690  label srcPatchi = srcPatchID_[i];
691  label tgtPatchi = tgtPatchID_[i];
692 
693  const polyPatch& srcPP = srcRegion_.boundaryMesh()[srcPatchi];
694  const polyPatch& tgtPP = tgtRegion_.boundaryMesh()[tgtPatchi];
695 
696  Info<< "Creating AMI between source patch " << srcPP.name()
697  << " and target patch " << tgtPP.name()
698  << " using " << AMIMethodName
699  << endl;
700 
701  Info<< incrIndent;
702 
703  patchAMIs_.set
704  (
705  i,
707  (
708  AMIMethodName,
709  false, // requireMatch
710  true, // flip target patch since patch normals are aligned
711  -1 // low weight correction
712  )
713  );
714 
715  patchAMIs_[i].calculate(srcPP, tgtPP);
716 
717  Info<< decrIndent;
718  }
719 }
720 
721 
722 void Foam::meshToMesh::constructNoCuttingPatches
723 (
724  const word& methodName,
725  const word& AMIMethodName,
726  const bool interpAllPatches
727 )
728 {
729  if (interpAllPatches)
730  {
731  const polyBoundaryMesh& srcBM = srcRegion_.boundaryMesh();
732  const polyBoundaryMesh& tgtBM = tgtRegion_.boundaryMesh();
733 
734  DynamicList<label> srcPatchID(srcBM.size());
735  DynamicList<label> tgtPatchID(tgtBM.size());
736  forAll(srcBM, patchi)
737  {
738  const polyPatch& pp = srcBM[patchi];
739 
740  // We want to map all the global patches, including constraint
741  // patches (since they might have mappable properties, e.g.
742  // jumpCyclic). We'll fix the value afterwards.
743  if (!isA<processorPolyPatch>(pp))
744  {
745  srcPatchID.append(pp.index());
746 
747  label tgtPatchi = tgtBM.findPatchID(pp.name());
748 
749  if (tgtPatchi != -1)
750  {
751  tgtPatchID.append(tgtPatchi);
752  }
753  else
754  {
756  << "Source patch " << pp.name()
757  << " not found in target mesh. "
758  << "Available target patches are " << tgtBM.names()
759  << exit(FatalError);
760  }
761  }
762  }
763 
764  srcPatchID_.transfer(srcPatchID);
765  tgtPatchID_.transfer(tgtPatchID);
766  }
767 
768  // calculate volume addressing and weights
769  calculate(methodName, true);
770 
771  // calculate patch addressing and weights
772  calculatePatchAMIs(AMIMethodName);
773 }
774 
775 
776 void Foam::meshToMesh::constructFromCuttingPatches
777 (
778  const word& methodName,
779  const word& AMIMethodName,
780  const HashTable<word>& patchMap,
781  const wordList& cuttingPatches,
782  const bool normalise
783 )
784 {
785  const polyBoundaryMesh& srcBm = srcRegion_.boundaryMesh();
786  const polyBoundaryMesh& tgtBm = tgtRegion_.boundaryMesh();
787 
788  // set IDs of cutting patches
789  cuttingPatches_.setSize(cuttingPatches.size());
790  forAll(cuttingPatches_, i)
791  {
792  const word& patchName = cuttingPatches[i];
793  label cuttingPatchi = srcBm.findPatchID(patchName);
794 
795  if (cuttingPatchi == -1)
796  {
798  << "Unable to find patch '" << patchName
799  << "' in mesh '" << srcRegion_.name() << "'. "
800  << " Available patches include:" << srcBm.names()
801  << exit(FatalError);
802  }
803 
804  cuttingPatches_[i] = cuttingPatchi;
805  }
806 
807  DynamicList<label> srcIDs(patchMap.size());
808  DynamicList<label> tgtIDs(patchMap.size());
809 
810  forAllConstIters(patchMap, iter)
811  {
812  const word& tgtPatchName = iter.key();
813  const word& srcPatchName = iter.val();
814 
815  const polyPatch& srcPatch = srcBm[srcPatchName];
816 
817  // We want to map all the global patches, including constraint
818  // patches (since they might have mappable properties, e.g.
819  // jumpCyclic). We'll fix the value afterwards.
820  if (!isA<processorPolyPatch>(srcPatch))
821  {
822  const polyPatch& tgtPatch = tgtBm[tgtPatchName];
823 
824  srcIDs.append(srcPatch.index());
825  tgtIDs.append(tgtPatch.index());
826  }
827  }
828 
829  srcPatchID_.transfer(srcIDs);
830  tgtPatchID_.transfer(tgtIDs);
831 
832  // calculate volume addressing and weights
833  calculate(methodName, normalise);
834 
835  // calculate patch addressing and weights
836  calculatePatchAMIs(AMIMethodName);
837 }
838 
839 
840 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
841 
842 Foam::meshToMesh::meshToMesh
843 (
844  const polyMesh& src,
845  const polyMesh& tgt,
846  const interpolationMethod& method,
847  const procMapMethod& mapMethod,
848  bool interpAllPatches
849 )
850 :
851  srcRegion_(src),
852  tgtRegion_(tgt),
853  procMapMethod_(mapMethod),
854  srcPatchID_(),
855  tgtPatchID_(),
856  patchAMIs_(),
857  cuttingPatches_(),
858  srcToTgtCellAddr_(),
859  tgtToSrcCellAddr_(),
860  srcToTgtCellWght_(),
861  tgtToSrcCellWght_(),
862  srcToTgtCellVec_(),
863  tgtToSrcCellVec_(),
864  V_(0.0),
865  singleMeshProc_(-1),
866  srcMapPtr_(nullptr),
867  tgtMapPtr_(nullptr)
868 {
869  constructNoCuttingPatches
870  (
871  interpolationMethodNames_[method],
872  interpolationMethodAMI(method),
873  interpAllPatches
874  );
875 }
876 
877 
878 Foam::meshToMesh::meshToMesh
879 (
880  const polyMesh& src,
881  const polyMesh& tgt,
882  const word& methodName,
883  const word& AMIMethodName,
884  const procMapMethod& mapMethod,
885  bool interpAllPatches
886 )
887 :
888  srcRegion_(src),
889  tgtRegion_(tgt),
890  procMapMethod_(mapMethod),
891  srcPatchID_(),
892  tgtPatchID_(),
893  patchAMIs_(),
894  cuttingPatches_(),
895  srcToTgtCellAddr_(),
896  tgtToSrcCellAddr_(),
897  srcToTgtCellWght_(),
898  tgtToSrcCellWght_(),
899  srcToTgtCellVec_(),
900  tgtToSrcCellVec_(),
901  V_(0.0),
902  singleMeshProc_(-1),
903  srcMapPtr_(nullptr),
904  tgtMapPtr_(nullptr)
905 {
906  constructNoCuttingPatches(methodName, AMIMethodName, interpAllPatches);
907 }
908 
909 
910 Foam::meshToMesh::meshToMesh
911 (
912  const polyMesh& src,
913  const polyMesh& tgt,
914  const interpolationMethod& method,
915  const HashTable<word>& patchMap,
916  const wordList& cuttingPatches,
917  const procMapMethod& mapMethod,
918  const bool normalise
919 )
920 :
921  srcRegion_(src),
922  tgtRegion_(tgt),
923  procMapMethod_(mapMethod),
924  srcPatchID_(),
925  tgtPatchID_(),
926  patchAMIs_(),
927  cuttingPatches_(),
928  srcToTgtCellAddr_(),
929  tgtToSrcCellAddr_(),
930  srcToTgtCellWght_(),
931  tgtToSrcCellWght_(),
932  V_(0.0),
933  singleMeshProc_(-1),
934  srcMapPtr_(nullptr),
935  tgtMapPtr_(nullptr)
936 {
937  constructFromCuttingPatches
938  (
939  interpolationMethodNames_[method],
940  interpolationMethodAMI(method),
941  patchMap,
942  cuttingPatches,
943  normalise
944  );
945 }
946 
947 
948 Foam::meshToMesh::meshToMesh
949 (
950  const polyMesh& src,
951  const polyMesh& tgt,
952  const word& methodName, // internal mapping
953  const word& AMIMethodName, // boundary mapping
954  const HashTable<word>& patchMap,
955  const wordList& cuttingPatches,
956  const procMapMethod& mapMethod,
957  const bool normalise
958 )
959 :
960  srcRegion_(src),
961  tgtRegion_(tgt),
962  procMapMethod_(mapMethod),
963  srcPatchID_(),
964  tgtPatchID_(),
965  patchAMIs_(),
966  cuttingPatches_(),
967  srcToTgtCellAddr_(),
968  tgtToSrcCellAddr_(),
969  srcToTgtCellWght_(),
970  tgtToSrcCellWght_(),
971  V_(0.0),
972  singleMeshProc_(-1),
973  srcMapPtr_(nullptr),
974  tgtMapPtr_(nullptr)
975 {
976  constructFromCuttingPatches
977  (
978  methodName,
979  AMIMethodName,
980  patchMap,
981  cuttingPatches,
982  normalise
983  );
984 }
985 
986 
987 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
988 
990 {}
991 
992 
993 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
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:1069
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:643
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:194
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
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 >
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:369
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:989
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
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:346
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:62
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 (stdout output on master, null elsewhere)
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(), const label comm=UPstream::worldComm)
Distribute data. Note:schedule only used for.
Definition: mapDistributeBaseTemplates.C:122
Foam::meshToMesh::procMapMethod
procMapMethod
Enumeration specifying processor parallel map construction method.
Definition: meshToMesh.H:82
Foam::OSstream::name
virtual const fileName & name() const
Get 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:353
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
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::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:450
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:1094
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::UPstream::commsTypes::nonBlocking
Foam::UPstream::msgType
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:540
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:463
Foam::nl
constexpr char nl
Definition: Ostream.H:404
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:80
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
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
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:188
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)