35 namespace cellsToCellss
45 bool Foam::cellsToCellss::matching::intersect
47 const polyMesh& srcMesh,
48 const polyMesh& tgtMesh,
53 return tgtMesh.pointInCell
55 srcMesh.cellCentres()[srcCelli],
62 bool Foam::cellsToCellss::matching::findInitialSeeds
64 const polyMesh& srcMesh,
65 const polyMesh& tgtMesh,
68 const label startSeedI,
73 const cellList& srcCells = srcMesh.cells();
74 const faceList& srcFaces = srcMesh.faces();
77 for (
label i = startSeedI; i < srcCellIDs.size(); i++)
79 label srcI = srcCellIDs[i];
83 const point srcCtr(srcCells[srcI].centre(srcPts, srcFaces));
84 label tgtI = tgtMesh.cellTree().findInside(srcCtr);
86 if (tgtI != -1 && intersect(srcMesh, tgtMesh, srcI, tgtI))
98 Pout<<
"could not find starting seed" <<
endl;
105 Foam::scalar Foam::cellsToCellss::matching::calculateAddressing
107 const polyMesh& srcMesh,
108 const polyMesh& tgtMesh,
113 const label srcSeedI,
114 const label tgtSeedI,
123 labelList srcTgtSeed(srcMesh.nCells(), -1);
125 List<DynamicList<label>> srcToTgt(srcMesh.nCells());
126 List<DynamicList<label>> tgtToSrc(tgtMesh.nCells());
128 DynamicList<label> srcSeeds(10);
133 label srcCelli = srcSeedI;
134 label tgtCelli = tgtSeedI;
139 srcToTgt[srcCelli].append(tgtCelli);
140 tgtToSrc[tgtCelli].append(srcCelli);
143 mapFlag[srcCelli] =
false;
146 V += srcVc[srcCelli];
160 while (srcCelli >= 0);
163 forAll(srcToTgtCellAddr, i)
165 srcToTgtCellWght[i] =
scalarList(srcToTgt[i].size(), srcVc[i]);
166 srcToTgtCellAddr[i].transfer(srcToTgt[i]);
169 forAll(tgtToSrcCellAddr, i)
171 tgtToSrcCellWght[i] =
scalarList(tgtToSrc[i].size(), tgtVc[i]);
172 tgtToSrcCellAddr[i].transfer(tgtToSrc[i]);
179 void Foam::cellsToCellss::matching::appendToDirectSeeds
181 const polyMesh& srcMesh,
182 const polyMesh& tgtMesh,
185 DynamicList<label>& srcSeeds,
190 const labelList& srcNbr = srcMesh.cellCells()[srcSeedI];
191 const labelList& tgtNbr = tgtMesh.cellCells()[tgtSeedI];
195 label srcI = srcNbr[i];
197 if (mapFlag[srcI] && (srcTgtSeed[srcI] == -1))
205 label tgtI = tgtNbr[j];
207 if (intersect(srcMesh, tgtMesh, srcI, tgtI))
212 srcTgtSeed[srcI] = tgtI;
213 srcSeeds.append(srcI);
222 mapFlag[srcI] =
false;
229 srcSeedI = srcSeeds.remove();
230 tgtSeedI = srcTgtSeed[srcSeedI];
247 initialise(srcMesh, tgtMesh);
250 const labelList srcCellIDs(maskCells(srcMesh, tgtMesh));
259 label startSeedI = 0;
304 forAll(srcToTgtWght, srcCelli)
306 if (srcToTgtWght[srcCelli].size() > 1)
308 srcToTgtAddr[srcCelli].
resize(1);
309 srcToTgtWght[srcCelli].
resize(1);
310 srcToTgtWght[srcCelli][0] = 1;
#define forAll(list, i)
Loop across all elements in list.
Macros for easy insertion into run-time selection tables.
void resize(const label)
Alias for setSize(const label)
A List with indirect addressing.
Class to calculate interpolative addressing and weights between the cells of two overlapping meshes.
Matching, one-to-one, cells-to-cells interpolation class.
virtual scalar calculate(const polyMesh &srcMesh, const polyMesh &tgtMesh)
Calculate the addressing and weights.
virtual void normalise(const polyMesh &mesh, labelListList &localOtherCells, scalarListList &weights) const
Normalise the weights for a given mesh.
virtual ~matching()
Destructor.
matching()
Construct null.
Mesh consisting of general polyhedral cells.
A class for handling words, derived from string.
defineTypeNameAndDebug(intersection, 0)
addToRunTimeSelectionTable(cellsToCells, intersection, word)
List< scalarList > scalarListList
List< label > labelList
A List of labels.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
List< cell > cellList
list of cells
Ostream & endl(Ostream &os)
Add newline and flush stream.
vectorField pointField
pointField is a vectorField.
List< bool > boolList
Bool container classes.
List< scalar > scalarList
A List of scalars.
vector point
Point is a vector.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
List< labelList > labelListList
A List of labelList.
prefixOSstream Pout(cout, "Pout")