36 namespace cellsToCellss
44 const Foam::scalar Foam::cellsToCellss::intersection::tolerance_ = 1
e-6;
49 bool Foam::cellsToCellss::intersection::intersect
51 const polyMesh& srcMesh,
52 const polyMesh& tgtMesh,
58 tetOverlapVolume().cellCellOverlapMinDecomp
64 treeBoundBox(tgtMesh.points(), tgtMesh.cellPoints()[tgtCelli]),
65 tolerance_*srcMesh.cellVolumes()[srcCelli]
70 Foam::scalar Foam::cellsToCellss::intersection::interVol
72 const polyMesh& srcMesh,
73 const polyMesh& tgtMesh,
79 tetOverlapVolume().cellCellOverlapVolumeMinDecomp
85 treeBoundBox(tgtMesh.points(), tgtMesh.cellPoints()[tgtCelli])
90 bool Foam::cellsToCellss::intersection::findInitialSeeds
92 const polyMesh& srcMesh,
93 const polyMesh& tgtMesh,
96 const label startSeedI,
101 const cellList& srcCells = srcMesh.cells();
102 const faceList& srcFaces = srcMesh.faces();
105 for (
label i = startSeedI; i < srcCellIDs.size(); i++)
107 const label srcI = srcCellIDs[i];
113 tgtMesh.cellTree().findBox
115 treeBoundBox(srcCells[srcI].bb(srcPts, srcFaces))
121 const label tgtI = tgtIDs[j];
123 if (intersect(srcMesh, tgtMesh, srcI, tgtI))
136 Pout<<
"could not find starting seed" <<
endl;
143 Foam::scalar Foam::cellsToCellss::intersection::calculateAddressing
145 const polyMesh& srcMesh,
146 const polyMesh& tgtMesh,
151 const label srcSeedI,
152 const label tgtSeedI,
160 label srcCelli = srcSeedI;
161 label tgtCelli = tgtSeedI;
163 List<DynamicList<label>> srcToTgtAddr(srcMesh.nCells());
164 List<DynamicList<scalar>> srcToTgtWght(srcMesh.nCells());
166 List<DynamicList<label>> tgtToSrcAddr(tgtMesh.nCells());
167 List<DynamicList<scalar>> tgtToSrcWght(tgtMesh.nCells());
170 DynamicList<label> nbrTgtCells(10);
173 DynamicList<label> visitedTgtCells(10);
176 labelList seedCells(srcMesh.nCells(), -1);
177 seedCells[srcCelli] = tgtCelli;
184 visitedTgtCells.clear();
187 nbrTgtCells.append(tgtCelli);
188 appendNbrCells(tgtCelli, tgtMesh, visitedTgtCells, nbrTgtCells);
192 tgtCelli = nbrTgtCells.remove();
193 visitedTgtCells.append(tgtCelli);
195 scalar vol = interVol(srcMesh, tgtMesh, srcCelli, tgtCelli);
198 if (vol/srcVol[srcCelli] > tolerance_)
201 srcToTgtAddr[srcCelli].append(tgtCelli);
202 srcToTgtWght[srcCelli].append(vol);
204 tgtToSrcAddr[tgtCelli].append(srcCelli);
205 tgtToSrcWght[tgtCelli].append(vol);
207 appendNbrCells(tgtCelli, tgtMesh, visitedTgtCells, nbrTgtCells);
213 while (!nbrTgtCells.empty());
215 mapFlag[srcCelli] =
false;
231 while (srcCelli != -1);
234 forAll(srcToTgtCellAddr, i)
236 srcToTgtCellAddr[i].transfer(srcToTgtAddr[i]);
237 srcToTgtCellWght[i].transfer(srcToTgtWght[i]);
240 forAll(tgtToSrcCellAddr, i)
242 tgtToSrcCellAddr[i].transfer(tgtToSrcAddr[i]);
243 tgtToSrcCellWght[i].transfer(tgtToSrcWght[i]);
250 void Foam::cellsToCellss::intersection::setNextCells
252 const polyMesh& srcMesh,
253 const polyMesh& tgtMesh,
259 const DynamicList<label>& visitedCells,
263 const labelList& srcNbrCells = srcMesh.cellCells()[srcCelli];
267 bool valuesSet =
false;
270 label cellS = srcNbrCells[i];
272 if (mapFlag[cellS] && seedCells[cellS] == -1)
276 label cellT = visitedCells[j];
278 if (intersect(srcMesh, tgtMesh, cellS, cellT))
280 seedCells[cellS] = cellT;
301 bool foundNextSeed =
false;
302 for (
label i = startSeedI; i < srcCellIDs.size(); i++)
304 label cellS = srcCellIDs[i];
311 foundNextSeed =
true;
314 if (seedCells[cellS] != -1)
317 tgtCelli = seedCells[cellS];
327 Pout<<
"Advancing front stalled: searching for new "
328 <<
"target cell" <<
endl;
364 initialise(srcMesh, tgtMesh);
367 const labelList srcCellIDs(maskCells(srcMesh, tgtMesh));
376 label startSeedI = 0;
423 forAll(srcToTgtWght, srcCelli)
427 forAll(srcToTgtWght[srcCelli], i)
429 srcToTgtWght[srcCelli][i] /= v;
#define forAll(list, i)
Loop across all elements in list.
Macros for easy insertion into run-time selection tables.
void clear()
Clear the list, i.e. set size to zero.
A List with indirect addressing.
Class to calculate interpolative addressing and weights between the cells of two overlapping meshes.
Intersection-based cells-to-cells interpolation class. Volume conservative.
virtual scalar calculate(const polyMesh &srcMesh, const polyMesh &tgtMesh)
Calculate the addressing and weights.
virtual ~intersection()
Destructor.
virtual void normalise(const polyMesh &mesh, labelListList &localOtherCells, scalarListList &weights) const
Normalise the weights for a given mesh.
intersection()
Construct null.
Mesh consisting of general polyhedral cells.
Calculates the overlap volume of two cells using tetrahedral decomposition.
scalar cellVolumeMinDecomp(const primitiveMesh &mesh, const label celli) const
Calculates the cell volume.
A class for handling words, derived from string.
defineTypeNameAndDebug(intersection, 0)
addToRunTimeSelectionTable(cellsToCells, intersection, word)
const dimensionedScalar e
Elementary charge.
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.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
List< labelList > labelListList
A List of labelList.
prefixOSstream Pout(cout, "Pout")