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 | Copyright (C) 2012-2014 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "meshToMesh.H"
27 #include "Time.H"
28 #include "globalIndex.H"
29 #include "meshToMeshMethod.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  defineTypeNameAndDebug(meshToMesh, 0);
36 
37  template<>
38  const char* Foam::NamedEnum
39  <
41  4
42  >::names[] =
43  {
44  "direct",
45  "mapNearest",
46  "cellVolumeWeight",
47  "correctedCellVolumeWeight"
48  };
49 
52 }
53 
54 
55 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
56 
57 template<>
59 (
61  const plusEqOp<sphericalTensor>& cop,
63  const bool secondOrder
64 ) const
65 {
66  mapSrcToTgt(field, cop, result.internalField());
67 }
68 
69 
70 template<>
72 (
74  const minusEqOp<sphericalTensor>& cop,
76  const bool secondOrder
77 ) const
78 {
79  mapSrcToTgt(field, cop, result.internalField());
80 }
81 
82 
83 template<>
85 (
87  const plusEqOp<symmTensor>& cop,
89  const bool secondOrder
90 ) const
91 {
92  mapSrcToTgt(field, cop, result.internalField());
93 }
94 
95 
96 template<>
98 (
100  const minusEqOp<symmTensor>& cop,
102  const bool secondOrder
103 ) const
104 {
105  mapSrcToTgt(field, cop, result.internalField());
106 }
107 
108 
109 template<>
111 (
113  const plusEqOp<tensor>& cop,
115  const bool secondOrder
116 ) const
117 {
118  mapSrcToTgt(field, cop, result.internalField());
119 }
120 
121 
122 template<>
124 (
126  const minusEqOp<tensor>& cop,
128  const bool secondOrder
129 ) const
130 {
131  mapSrcToTgt(field, cop, result.internalField());
132 }
133 
134 
135 template<>
137 (
139  const plusEqOp<sphericalTensor>& cop,
141  const bool secondOrder
142 ) const
143 {
144  mapSrcToTgt(field, cop, result.internalField());
145 }
146 
147 
148 template<>
150 (
152  const minusEqOp<sphericalTensor>& cop,
154  const bool secondOrder
155 ) const
156 {
157  mapSrcToTgt(field, cop, result.internalField());
158 }
159 
160 
161 template<>
163 (
165  const plusEqOp<symmTensor>& cop,
167  const bool secondOrder
168 ) const
169 {
170  mapSrcToTgt(field, cop, result.internalField());
171 }
172 
173 
174 template<>
176 (
178  const minusEqOp<symmTensor>& cop,
180  const bool secondOrder
181 ) const
182 {
183  mapSrcToTgt(field, cop, result.internalField());
184 }
185 
186 
187 template<>
189 (
191  const plusEqOp<tensor>& cop,
193  const bool secondOrder
194 ) const
195 {
196  mapSrcToTgt(field, cop, result.internalField());
197 }
198 
199 
200 template<>
202 (
204  const minusEqOp<tensor>& cop,
206  const bool secondOrder
207 ) const
208 {
209  mapSrcToTgt(field, cop, result.internalField());
210 }
211 
212 
213 template<>
215 (
216  const AMIPatchToPatchInterpolation& AMI,
217  const Field<scalar>& srcField,
218  Field<scalar>& tgtField,
219  const plusEqOp<scalar>& cop
220 ) const
221 {}
222 
223 
224 template<>
226 (
227  const AMIPatchToPatchInterpolation& AMI,
228  const Field<vector>& srcField,
229  Field<vector>& tgtField,
230  const plusEqOp<vector>& cop
231 ) const
232 {}
233 
234 
235 template<>
237 (
238  const AMIPatchToPatchInterpolation& AMI,
239  const Field<sphericalTensor>& srcField,
240  Field<sphericalTensor>& tgtField,
241  const plusEqOp<sphericalTensor>& cop
242 ) const
243 {}
244 
245 
246 template<>
248 (
249  const AMIPatchToPatchInterpolation& AMI,
250  const Field<symmTensor>& srcField,
251  Field<symmTensor>& tgtField,
252  const plusEqOp<symmTensor>& cop
253 ) const
254 {}
255 
256 
257 template<>
259 (
260  const AMIPatchToPatchInterpolation& AMI,
261  const Field<tensor>& srcField,
262  Field<tensor>& tgtField,
263  const plusEqOp<tensor>& cop
264 ) const
265 {}
266 
267 
268 template<>
270 (
271  const AMIPatchToPatchInterpolation& AMI,
272  Field<scalar>& srcField,
273  const Field<scalar>& tgtField,
274  const plusEqOp<scalar>& cop
275 ) const
276 {}
277 
278 
279 template<>
281 (
282  const AMIPatchToPatchInterpolation& AMI,
283  Field<vector>& srcField,
284  const Field<vector>& tgtField,
285  const plusEqOp<vector>& cop
286 ) const
287 {}
288 
289 
290 template<>
292 (
293  const AMIPatchToPatchInterpolation& AMI,
294  Field<sphericalTensor>& srcField,
295  const Field<sphericalTensor>& tgtField,
296  const plusEqOp<sphericalTensor>& cop
297 ) const
298 {}
299 
300 
301 template<>
303 (
304  const AMIPatchToPatchInterpolation& AMI,
305  Field<symmTensor>& srcField,
306  const Field<symmTensor>& tgtField,
307  const plusEqOp<symmTensor>& cop
308 ) const
309 {}
310 
311 
312 template<>
314 (
315  const AMIPatchToPatchInterpolation& AMI,
316  Field<tensor>& srcField,
317  const Field<tensor>& tgtField,
318  const plusEqOp<tensor>& cop
319 ) const
320 {}
321 
322 
324 (
325  const polyMesh& src,
326  const polyMesh& tgt
327 ) const
328 {
329  boundBox intersectBb
330  (
331  max(src.bounds().min(), tgt.bounds().min()),
332  min(src.bounds().max(), tgt.bounds().max())
333  );
334 
335  intersectBb.inflate(0.01);
336 
337  const cellList& srcCells = src.cells();
338  const faceList& srcFaces = src.faces();
339  const pointField& srcPts = src.points();
340 
342  forAll(srcCells, srcI)
343  {
344  boundBox cellBb(srcCells[srcI].points(srcFaces, srcPts), false);
345  if (intersectBb.overlaps(cellBb))
346  {
347  cells.append(srcI);
348  }
349  }
350 
351  if (debug)
352  {
353  Pout<< "participating source mesh cells: " << cells.size() << endl;
354  }
355 
356  return cells;
357 }
358 
359 
361 (
362  const word& descriptor,
363  const labelListList& addr,
364  scalarListList& wght
365 ) const
366 {
367  const label nCell = returnReduce(wght.size(), sumOp<label>());
368 
369  if (nCell > 0)
370  {
371  forAll(wght, cellI)
372  {
373  scalarList& w = wght[cellI];
374  scalar s = sum(w);
375 
376  forAll(w, i)
377  {
378  // note: normalise by s instead of cell volume since
379  // 1-to-1 methods duplicate contributions in parallel
380  w[i] /= s;
381  }
382  }
383  }
384 }
385 
386 
388 (
389  const word& methodName,
390  const polyMesh& src,
391  const polyMesh& tgt
392 )
393 {
394  autoPtr<meshToMeshMethod> methodPtr
395  (
397  (
398  methodName,
399  src,
400  tgt
401  )
402  );
403 
404  methodPtr->calculate
405  (
406  srcToTgtCellAddr_,
407  srcToTgtCellWght_,
408  srcToTgtCellVec_,
409  tgtToSrcCellAddr_,
410  tgtToSrcCellWght_,
411  tgtToSrcCellVec_
412  );
413 
414  V_ = methodPtr->V();
415 
416  if (debug)
417  {
418  methodPtr->writeConnectivity(src, tgt, srcToTgtCellAddr_);
419  }
420 }
421 
422 
423 void Foam::meshToMesh::calculate(const word& methodName)
424 {
425  Info<< "Creating mesh-to-mesh addressing for " << srcRegion_.name()
426  << " and " << tgtRegion_.name() << " regions using "
427  << methodName << endl;
428 
430 
431  if (singleMeshProc_ == -1)
432  {
433  // create global indexing for src and tgt meshes
434  globalIndex globalSrcCells(srcRegion_.nCells());
435  globalIndex globalTgtCells(tgtRegion_.nCells());
436 
437  // Create processor map of overlapping cells. This map gets
438  // (possibly remote) cells from the tgt mesh such that they (together)
439  // cover all of the src mesh
441  const mapDistribute& map = mapPtr();
442 
443  pointField newTgtPoints;
444  faceList newTgtFaces;
445  labelList newTgtFaceOwners;
446  labelList newTgtFaceNeighbours;
447  labelList newTgtCellIDs;
448 
450  (
451  map,
452  tgtRegion_,
453  globalTgtCells,
454  newTgtPoints,
455  newTgtFaces,
456  newTgtFaceOwners,
457  newTgtFaceNeighbours,
458  newTgtCellIDs
459  );
460 
461 
462  // create a new target mesh
463  polyMesh newTgt
464  (
465  IOobject
466  (
467  "newTgt." + Foam::name(Pstream::myProcNo()),
469  tgtRegion_.time(),
471  ),
472  xferMove(newTgtPoints),
473  xferMove(newTgtFaces),
474  xferMove(newTgtFaceOwners),
475  xferMove(newTgtFaceNeighbours),
476  false // no parallel comms
477  );
478 
479  // create some dummy patch info
481  patches[0] = new polyPatch
482  (
483  "defaultFaces",
484  newTgt.nFaces() - newTgt.nInternalFaces(),
485  newTgt.nInternalFaces(),
486  0,
487  newTgt.boundaryMesh(),
488  word::null
489  );
490 
491  newTgt.addPatches(patches);
492 
493  // force calculation of tet-base points used for point-in-cell
494  (void)newTgt.tetBasePtIs();
495 
496  // force construction of cell tree
497 // (void)newTgt.cellTree();
498 
499  if (debug)
500  {
501  Pout<< "Created newTgt mesh:" << nl
502  << " old cells = " << tgtRegion_.nCells()
503  << ", new cells = " << newTgt.nCells() << nl
504  << " old faces = " << tgtRegion_.nFaces()
505  << ", new faces = " << newTgt.nFaces() << endl;
506 
507  if (debug > 1)
508  {
509  Pout<< "Writing newTgt mesh: " << newTgt.name() << endl;
510  newTgt.write();
511  }
512  }
513 
514  calcAddressing(methodName, srcRegion_, newTgt);
515 
516  // per source cell the target cell address in newTgt mesh
518  {
519  labelList& addressing = srcToTgtCellAddr_[i];
520  forAll(addressing, addrI)
521  {
522  addressing[addrI] = newTgtCellIDs[addressing[addrI]];
523  }
524  }
525 
526  // convert target addresses in newTgtMesh into global cell numbering
528  {
529  labelList& addressing = tgtToSrcCellAddr_[i];
530  forAll(addressing, addrI)
531  {
532  addressing[addrI] = globalSrcCells.toGlobal(addressing[addrI]);
533  }
534  }
535 
536  // set up as a reverse distribute
538  (
540  List<labelPair>(),
541  tgtRegion_.nCells(),
542  map.constructMap(),
543  false,
544  map.subMap(),
545  false,
548  flipOp(),
549  labelList()
550  );
551 
552  // set up as a reverse distribute
554  (
556  List<labelPair>(),
557  tgtRegion_.nCells(),
558  map.constructMap(),
559  false,
560  map.subMap(),
561  false,
564  flipOp(),
565  scalarList()
566  );
567 
568  // weights normalisation
570  (
571  "source",
574  );
575 
577  (
578  "target",
581  );
582 
583  // cache maps and reset addresses
584  List<Map<label> > cMap;
585  srcMapPtr_.reset
586  (
587  new mapDistribute(globalSrcCells, tgtToSrcCellAddr_, cMap)
588  );
589  tgtMapPtr_.reset
590  (
591  new mapDistribute(globalTgtCells, srcToTgtCellAddr_, cMap)
592  );
593 
594  // collect volume intersection contributions
595  reduce(V_, sumOp<scalar>());
596  }
597  else
598  {
599  calcAddressing(methodName, srcRegion_, tgtRegion_);
600 
602  (
603  "source",
606  );
607 
609  (
610  "target",
613  );
614  }
615 
616  Info<< " Overlap volume: " << V_ << endl;
617 }
618 
619 
622 {
623  switch (method)
624  {
625  case imDirect:
626  {
628  break;
629  }
630  case imMapNearest:
631  {
633  break;
634  }
635  case imCellVolumeWeight:
636  case imCorrectedCellVolumeWeight:
637  {
639  break;
640  }
641  default:
642  {
644  << "Unhandled enumeration " << method
645  << abort(FatalError);
646  }
647  }
648 
650 }
651 
652 
653 void Foam::meshToMesh::calculatePatchAMIs(const word& AMIMethodName)
654 {
655  if (!patchAMIs_.empty())
656  {
658  << "patch AMI already calculated"
659  << exit(FatalError);
660  }
661 
662  patchAMIs_.setSize(srcPatchID_.size());
663 
664  forAll(srcPatchID_, i)
665  {
666  label srcPatchI = srcPatchID_[i];
667  label tgtPatchI = tgtPatchID_[i];
668 
669  const polyPatch& srcPP = srcRegion_.boundaryMesh()[srcPatchI];
670  const polyPatch& tgtPP = tgtRegion_.boundaryMesh()[tgtPatchI];
671 
672  Info<< "Creating AMI between source patch " << srcPP.name()
673  << " and target patch " << tgtPP.name()
674  << " using " << AMIMethodName
675  << endl;
676 
677  Info<< incrIndent;
678 
679  patchAMIs_.set
680  (
681  i,
683  (
684  srcPP,
685  tgtPP,
687  false,
688  AMIMethodName,
689  -1,
690  true // flip target patch since patch normals are aligned
691  )
692  );
693 
694  Info<< decrIndent;
695  }
696 }
697 
698 
700 (
701  const word& methodName,
702  const word& AMIMethodName,
703  const bool interpAllPatches
704 )
705 {
706  if (interpAllPatches)
707  {
708  const polyBoundaryMesh& srcBM = srcRegion_.boundaryMesh();
709  const polyBoundaryMesh& tgtBM = tgtRegion_.boundaryMesh();
710 
711  DynamicList<label> srcPatchID(srcBM.size());
712  DynamicList<label> tgtPatchID(tgtBM.size());
713  forAll(srcBM, patchI)
714  {
715  const polyPatch& pp = srcBM[patchI];
716  if (!polyPatch::constraintType(pp.type()))
717  {
718  srcPatchID.append(pp.index());
719 
720  label tgtPatchI = tgtBM.findPatchID(pp.name());
721 
722  if (tgtPatchI != -1)
723  {
724  tgtPatchID.append(tgtPatchI);
725  }
726  else
727  {
729  << "Source patch " << pp.name()
730  << " not found in target mesh. "
731  << "Available target patches are " << tgtBM.names()
732  << exit(FatalError);
733  }
734  }
735  }
736 
737  srcPatchID_.transfer(srcPatchID);
738  tgtPatchID_.transfer(tgtPatchID);
739  }
740 
741  // calculate volume addressing and weights
742  calculate(methodName);
743 
744  // calculate patch addressing and weights
745  calculatePatchAMIs(AMIMethodName);
746 }
747 
748 
750 (
751  const word& methodName,
752  const word& AMIMethodName,
753  const HashTable<word>& patchMap,
754  const wordList& cuttingPatches
755 )
756 {
757  DynamicList<label> srcIDs(patchMap.size());
758  DynamicList<label> tgtIDs(patchMap.size());
759 
760  forAllConstIter(HashTable<word>, patchMap, iter)
761  {
762  const word& tgtPatchName = iter.key();
763  const word& srcPatchName = iter();
764 
765  const polyPatch& srcPatch = srcRegion_.boundaryMesh()[srcPatchName];
766 
767  if (!polyPatch::constraintType(srcPatch.type()))
768  {
769  const polyPatch& tgtPatch = tgtRegion_.boundaryMesh()[tgtPatchName];
770 
771  srcIDs.append(srcPatch.index());
772  tgtIDs.append(tgtPatch.index());
773  }
774  }
775 
776  srcPatchID_.transfer(srcIDs);
777  tgtPatchID_.transfer(tgtIDs);
778 
779  // calculate volume addressing and weights
780  calculate(methodName);
781 
782  // calculate patch addressing and weights
783  calculatePatchAMIs(AMIMethodName);
784 
785  // set IDs of cutting patches on target mesh
786  cuttingPatches_.setSize(cuttingPatches.size());
787  forAll(cuttingPatches_, i)
788  {
789  const word& patchName = cuttingPatches[i];
790  cuttingPatches_[i] = tgtRegion_.boundaryMesh().findPatchID(patchName);
791  }
792 }
793 
794 
795 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
796 
798 (
799  const polyMesh& src,
800  const polyMesh& tgt,
801  const interpolationMethod& method,
802  bool interpAllPatches
803 )
804 :
805  srcRegion_(src),
806  tgtRegion_(tgt),
807  srcPatchID_(),
808  tgtPatchID_(),
809  patchAMIs_(),
810  cuttingPatches_(),
811  srcToTgtCellAddr_(),
812  tgtToSrcCellAddr_(),
813  srcToTgtCellWght_(),
814  tgtToSrcCellWght_(),
815  srcToTgtCellVec_(),
816  tgtToSrcCellVec_(),
817  V_(0.0),
818  singleMeshProc_(-1),
819  srcMapPtr_(NULL),
820  tgtMapPtr_(NULL)
821 {
822  constructNoCuttingPatches
823  (
824  interpolationMethodNames_[method],
826  (
827  interpolationMethodAMI(method)
828  ),
829  interpAllPatches
830  );
831 }
832 
833 
835 (
836  const polyMesh& src,
837  const polyMesh& tgt,
838  const word& methodName,
839  const word& AMIMethodName,
840  bool interpAllPatches
841 )
842 :
843  srcRegion_(src),
844  tgtRegion_(tgt),
845  srcPatchID_(),
846  tgtPatchID_(),
847  patchAMIs_(),
848  cuttingPatches_(),
849  srcToTgtCellAddr_(),
850  tgtToSrcCellAddr_(),
851  srcToTgtCellWght_(),
852  tgtToSrcCellWght_(),
853  srcToTgtCellVec_(),
854  tgtToSrcCellVec_(),
855  V_(0.0),
856  singleMeshProc_(-1),
857  srcMapPtr_(NULL),
858  tgtMapPtr_(NULL)
859 {
860  constructNoCuttingPatches(methodName, AMIMethodName, interpAllPatches);
861 }
862 
863 
865 (
866  const polyMesh& src,
867  const polyMesh& tgt,
868  const interpolationMethod& method,
869  const HashTable<word>& patchMap,
870  const wordList& cuttingPatches
871 )
872 :
873  srcRegion_(src),
874  tgtRegion_(tgt),
875  srcPatchID_(),
876  tgtPatchID_(),
877  patchAMIs_(),
878  cuttingPatches_(),
879  srcToTgtCellAddr_(),
880  tgtToSrcCellAddr_(),
881  srcToTgtCellWght_(),
882  tgtToSrcCellWght_(),
883  V_(0.0),
884  singleMeshProc_(-1),
885  srcMapPtr_(NULL),
886  tgtMapPtr_(NULL)
887 {
888  constructFromCuttingPatches
889  (
890  interpolationMethodNames_[method],
892  (
893  interpolationMethodAMI(method)
894  ),
895  patchMap,
896  cuttingPatches
897  );
898 }
899 
900 
902 (
903  const polyMesh& src,
904  const polyMesh& tgt,
905  const word& methodName, // internal mapping
906  const word& AMIMethodName, // boundary mapping
907  const HashTable<word>& patchMap,
908  const wordList& cuttingPatches
909 )
910 :
911  srcRegion_(src),
912  tgtRegion_(tgt),
913  srcPatchID_(),
914  tgtPatchID_(),
915  patchAMIs_(),
916  cuttingPatches_(),
917  srcToTgtCellAddr_(),
918  tgtToSrcCellAddr_(),
919  srcToTgtCellWght_(),
920  tgtToSrcCellWght_(),
921  V_(0.0),
922  singleMeshProc_(-1),
923  srcMapPtr_(NULL),
924  tgtMapPtr_(NULL)
925 {
926  constructFromCuttingPatches
927  (
928  methodName,
929  AMIMethodName,
930  patchMap,
931  cuttingPatches
932  );
933 }
934 
935 
936 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
937 
939 {}
940 
941 
942 // ************************************************************************* //
Foam::AMIInterpolation::interpolationMethod
interpolationMethod
Enumeration specifying interpolation method.
Definition: AMIInterpolation.H:86
Foam::mapDistributeBase::subMap
const labelListList & subMap() const
From subsetted data back to original data.
Definition: mapDistributeBase.H:256
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:979
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::scalarList
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
Foam::AMIInterpolation::interpolationMethodToWord
static word interpolationMethodToWord(const interpolationMethod &method)
Convert interpolationMethod to word representation.
Definition: AMIInterpolation.C:37
w
volScalarField w(IOobject("w", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE), mesh, dimensionedScalar("w", dimensionSet(0, 0, 0, 0, 0, 0, 0), 0.0))
Foam::meshToMesh::singleMeshProc_
label singleMeshProc_
Index of processor that holds all of both sides. -1 in all other.
Definition: meshToMesh.H:127
Foam::polyBoundaryMesh::findPatchID
label findPatchID(const word &patchName) const
Find patch index given a name.
Definition: polyBoundaryMesh.C:678
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::meshToMesh::meshToMesh
meshToMesh(const meshToMesh &)
Disallow default bitwise copy construct.
Foam::labelList
List< label > labelList
A List of labels.
Definition: labelList.H:56
Foam::boundBox::max
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:60
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
Foam::minusEqOp
Definition: ops.H:72
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::meshToMesh::distributeAndMergeCells
void distributeAndMergeCells(const mapDistribute &map, const polyMesh &tgt, const globalIndex &globalI, pointField &tgtPoints, faceList &tgtFaces, labelList &tgtFaceOwners, labelList &tgtFaceNeighbours, labelList &tgtCellIDs) const
Collect pieces of tgt mesh from other procssors and restructure.
Definition: meshToMeshParallelOps.C:504
Foam::meshToMesh::calculatePatchAMIs
void calculatePatchAMIs(const word &amiMethodName)
Calculate patch overlap.
Definition: meshToMesh.C:653
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 Pstream::scheduled.
Definition: mapDistributeBaseTemplates.C:119
Foam::meshToMesh::mapAndOpSrcToTgt
void mapAndOpSrcToTgt(const AMIPatchToPatchInterpolation &AMI, const Field< Type > &srcField, Field< Type > &tgtField, const CombineOp &cop) const
Helper function to interpolate patch field. Template.
Definition: meshToMeshTemplates.C:527
Foam::ListPlusEqOp
Plus op for FixedList<scalar>
Definition: regionSizeDistribution.C:41
Foam::polyMesh::tetBasePtIs
const labelList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:818
Foam::meshToMesh::constructNoCuttingPatches
void constructNoCuttingPatches(const word &methodName, const word &AMIMethodName, const bool interpAllPatches)
Constructor helper.
Definition: meshToMesh.C:700
Foam::AMIInterpolation::imFaceAreaWeight
@ imFaceAreaWeight
Definition: AMIInterpolation.H:90
Foam::boundBox::inflate
void inflate(const scalar s)
Inflate box by factor*mag(span) in all dimensions.
Definition: boundBox.C:210
globalIndex.H
Foam::objectRegistry::time
const Time & time() const
Return time.
Definition: objectRegistry.H:117
meshToMesh.H
Foam::flipOp
Class containing functor to negate primitives. Dummy for all other types.
Definition: flipOp.H:50
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:136
Foam::regIOobject::write
virtual bool write() const
Write using setting from DB.
Definition: regIOobjectWrite.C:126
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::meshToMesh::~meshToMesh
virtual ~meshToMesh()
Destructor.
Definition: meshToMesh.C:938
Foam::incrIndent
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:228
Foam::meshToMesh::calcDistribution
label calcDistribution(const polyMesh &src, const polyMesh &tgt) const
Determine whether the meshes are split across multiple pocessors.
Definition: meshToMeshParallelOps.C:38
Foam::meshToMesh::tgtToSrcCellWght_
scalarListList tgtToSrcCellWght_
Target to source cell interpolation weights.
Definition: meshToMesh.H:111
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::polyBoundaryMesh::names
wordList names() const
Return a list of patch names.
Definition: polyBoundaryMesh.C:527
Foam::meshToMesh::srcToTgtCellWght_
scalarListList srcToTgtCellWght_
Source to target cell interplation weights.
Definition: meshToMesh.H:108
Foam::meshToMesh::srcRegion_
const polyMesh & srcRegion_
Reference to the source mesh.
Definition: meshToMesh.H:84
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::primitiveMesh::nCells
label nCells() const
Definition: primitiveMeshI.H:64
Foam::meshToMesh::calcProcMap
autoPtr< mapDistribute > calcProcMap(const polyMesh &src, const polyMesh &tgt) const
Calculate the mapping between processors.
Definition: meshToMeshParallelOps.C:117
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
Foam::meshToMesh::interpolationMethodNames_
static const NamedEnum< interpolationMethod, 4 > interpolationMethodNames_
Definition: meshToMesh.H:77
Foam::meshToMesh::maskCells
labelList maskCells(const polyMesh &src, const polyMesh &tgt) const
Return src cell IDs for the overlap region.
Definition: meshToMesh.C:324
Foam::UPstream::nonBlocking
@ nonBlocking
Definition: UPstream.H:68
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:59
Foam::Field< scalar >
Foam::polyPatch::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:302
Foam::meshToMesh::mapInternalSrcToTgt
void mapInternalSrcToTgt(const GeometricField< Type, fvPatchField, volMesh > &field, const CombineOp &cop, GeometricField< Type, fvPatchField, volMesh > &result, const bool secondOrder) const
Helper function to interpolate internal field. Optionally uses.
Definition: meshToMeshTemplates.C:501
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::boundBox::overlaps
bool overlaps(const boundBox &) const
Overlaps/touches boundingBox?
Definition: boundBoxI.H:120
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::meshToMesh::tgtToSrcCellAddr_
labelListList tgtToSrcCellAddr_
Target to source cell addressing.
Definition: meshToMesh.H:105
Foam::mapDistribute
Class containing processor-to-processor mapping information.
Definition: mapDistribute.H:152
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
Foam::meshToMesh::calcAddressing
void calcAddressing(const word &methodName, const polyMesh &src, const polyMesh &tgt)
Calculate the addressing between overlapping regions of src and tgt.
Definition: meshToMesh.C:388
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobject.H:273
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
Foam::polyPatch::constraintType
static bool constraintType(const word &pt)
Return true if the given type is a constraint type.
Definition: polyPatch.C:270
Foam::meshToMesh::tgtMapPtr_
autoPtr< mapDistribute > tgtMapPtr_
Target map pointer - parallel running only.
Definition: meshToMesh.H:133
Foam::FatalError
error FatalError
Foam::HashTable::size
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const
Definition: primitiveMeshI.H:52
Foam::plusEqOp
Definition: ops.H:71
Foam::meshToMesh::constructFromCuttingPatches
void constructFromCuttingPatches(const word &methodName, const word &AMIMethodName, const HashTable< word > &patchMap, const wordList &cuttingPatches)
Constructor helper.
Definition: meshToMesh.C:750
Foam::meshToMesh::V_
scalar V_
Cell total volume in overlap region [m3].
Definition: meshToMesh.H:123
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::decrIndent
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
s
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().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))
meshToMeshMethod.H
Foam::meshToMesh::tgtRegion_
const polyMesh & tgtRegion_
Reference to the target mesh.
Definition: meshToMesh.H:87
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::meshToMesh::srcToTgtCellAddr_
labelListList srcToTgtCellAddr_
Source to target cell addressing.
Definition: meshToMesh.H:102
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:405
Foam::HashTable
An STL-conforming hash table.
Definition: HashTable.H:61
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::polyMesh::bounds
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:427
Foam::AMIInterpolation::imDirect
@ imDirect
Definition: AMIInterpolation.H:88
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1004
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::xferMove
Xfer< T > xferMove(T &)
Construct by transferring the contents of the arg.
Foam::primitiveMesh::nFaces
label nFaces() const
Definition: primitiveMeshI.H:58
Foam::meshToMesh::calculate
void calculate(const word &methodName)
Calculate - main driver function.
Definition: meshToMesh.C:423
Foam::sumOp
Definition: ops.H:162
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Foam::meshToMesh::interpolationMethod
interpolationMethod
Enumeration specifying interpolation method.
Definition: meshToMesh.H:68
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
Foam::polyMesh::addPatches
void addPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:878
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::AMIInterpolation
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
Definition: AMIInterpolation.H:77
Foam::meshToMesh::interpolationMethodAMI
static AMIPatchToPatchInterpolation::interpolationMethod interpolationMethodAMI(const interpolationMethod method)
Conversion between mesh and patch interpolation methods.
Definition: meshToMesh.C:621
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::PtrList::size
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
Foam::meshToMesh::srcMapPtr_
autoPtr< mapDistribute > srcMapPtr_
Source map pointer - parallel running only.
Definition: meshToMesh.H:130
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:333
Foam::word::null
static const word null
An empty word.
Definition: word.H:77
Foam::AMIInterpolation::imMapNearest
@ imMapNearest
Definition: AMIInterpolation.H:89
Foam::boundBox
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:55
patches
patches[0]
Definition: createSingleCellMesh.H:36
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::meshToMesh::mapAndOpTgtToSrc
void mapAndOpTgtToSrc(const AMIPatchToPatchInterpolation &AMI, Field< Type > &srcField, const Field< Type > &tgtField, const CombineOp &cop) const
Helper function to interpolate patch field. Template.
Definition: meshToMeshTemplates.C:766
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::mapDistributeBase::constructMap
const labelListList & constructMap() const
From subsetted data to new reconstructed data.
Definition: mapDistributeBase.H:268
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::patchIdentifier::name
const word & name() const
Return name.
Definition: patchIdentifier.H:109
Foam::meshToMesh::normaliseWeights
void normaliseWeights(const word &descriptor, const labelListList &addr, scalarListList &wght) const
Normalise the interpolation weights.
Definition: meshToMesh.C:361
Foam::patchIdentifier::index
label index() const
Return the index of this patch in the boundaryMesh.
Definition: patchIdentifier.H:133
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::NamedEnum
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:52
Foam::meshToMesh::mapInternalTgtToSrc
void mapInternalTgtToSrc(const GeometricField< Type, fvPatchField, volMesh > &field, const CombineOp &cop, GeometricField< Type, fvPatchField, volMesh > &result, const bool secondOrder) const
Helper function to interpolate internal field. Optionally uses.
Definition: meshToMeshTemplates.C:740
Foam::faceAreaIntersect::tmMesh
@ tmMesh
Definition: faceAreaIntersect.H:62
Foam::globalIndex::toGlobal
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82