Go to the documentation of this file.
46 const label startSeedI,
53 for (
label i = startSeedI; i < srcCellIDs.
size(); i++)
55 label srcI = srcCellIDs[i];
59 const point& srcCc = srcCcs[srcI];
61 tgt_.cellTree().findNearest(srcCc, GREAT);
66 tgtSeedI = hit.
index();
73 <<
"Unable to find nearest target cell"
74 <<
" for source cell " << srcI
75 <<
" with centre " << srcCc
85 Pout<<
"could not find starting seed" <<
endl;
112 label srcCellI = srcSeedI;
113 label tgtCellI = tgtSeedI;
118 findNearestCell(src_, tgt_, srcCellI, tgtCellI);
121 srcToTgt[srcCellI].
append(tgtCellI);
122 tgtToSrc[tgtCellI].
append(srcCellI);
125 mapFlag[srcCellI] =
false;
128 V_ += srcVc[srcCellI];
140 while (srcCellI >= 0);
148 forAll(tgtToSrc, targetCellI)
150 if (tgtToSrc[targetCellI].size() > 1)
152 const vector& tgtC = tgtCc[targetCellI];
156 label srcCellI = srcCells[0];
157 scalar d =
magSqr(tgtC - srcCc[srcCellI]);
159 for (
label i = 1; i < srcCells.size(); i++)
161 label srcI = srcCells[i];
162 scalar dNew =
magSqr(tgtC - srcCc[srcI]);
171 srcCells.
append(srcCellI);
177 forAll(tgtToSrc, tgtCellI)
179 if (tgtToSrc[tgtCellI].empty())
181 label srcCellI = findMappedSrcCell(tgtCellI, tgtToSrc);
183 findNearestCell(tgt_, src_, tgtCellI, srcCellI);
185 tgtToSrc[tgtCellI].
append(srcCellI);
190 forAll(srcToTgtCellAddr, i)
192 srcToTgtCellWght[i] =
scalarList(srcToTgt[i].size(), srcVc[i]);
193 srcToTgtCellAddr[i].
transfer(srcToTgt[i]);
196 forAll(tgtToSrcCellAddr, i)
198 tgtToSrcCellWght[i] =
scalarList(tgtToSrc[i].size(), tgtVc[i]);
199 tgtToSrcCellAddr[i].
transfer(tgtToSrc[i]);
215 const vector& p1 = Cc1[cell1];
229 scalar dTest =
magSqr(Cc2[
c2] - p1);
234 appendNbrCells(cell2, mesh2, visitedCells, cells2);
237 }
while (cells2.size() > 0);
250 const labelList& srcNbr = src_.cellCells()[srcCellI];
255 label cellI = srcNbr[i];
263 for (
label i = startSeedI; i < srcCellIDs.
size(); i++)
265 label cellI = srcCellIDs[i];
273 (void)findInitialSeeds
286 const label tgtCellI,
293 testCells.
append(tgtCellI);
302 visitedCells.
append(tgtI);
304 if (tgtToSrc[tgtI].size())
306 return tgtToSrc[tgtI][0];
310 const labelList& nbrCells = tgt_.cellCells()[tgtI];
314 if (
findIndex(visitedCells, nbrCells[i]) == -1)
316 testCells.
append(nbrCells[i]);
321 }
while (testCells.size());
375 boolList mapFlag(src_.nCells(),
false);
381 label startSeedI = 0;
List< scalar > scalarList
A List of scalars.
label index() const
Return index.
virtual void setNextNearestCells(label &startSeedI, label &srcCellI, label &tgtCellI, boolList &mapFlag, const labelList &srcCellIDs) const
Set the next cells for the marching front algorithm.
virtual void calculate(labelListList &srcToTgtAddr, scalarListList &srcToTgtWght, pointListList &srcToTgtVec, labelListList &tgtToSrcAddr, scalarListList &tgtToSrcWght, pointListList &tgtToSrcVec)
Calculate addressing and weights and optionally offset vectors.
#define forAll(list, i)
Loop across all elements in list.
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
bool hit() const
Is there a hit.
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Mesh consisting of general polyhedral cells.
virtual ~mapNearestMethod()
Destructor.
Base class for mesh-to-mesh calculation methods.
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
virtual bool findInitialSeeds(const labelList &srcCellIDs, const boolList &mapFlag, const label startSeedI, label &srcSeedI, label &tgtSeedI) const
Find indices of overlapping cells in src and tgt meshes - returns.
void append(const T &)
Append an element at the end of the list.
virtual label findMappedSrcCell(const label tgtCellI, const List< DynamicList< label > > &tgtToSrc) const
Find a source cell mapped to target cell tgtCellI.
void clear()
Clear the addressed list, i.e. set the size to zero.
virtual void calculateAddressing(labelListList &srcToTgtCellAddr, scalarListList &srcToTgtCellWght, labelListList &tgtToSrcCellAddr, scalarListList &tgtToSrcCellWght, const label srcSeedI, const label tgtSeedI, const labelList &srcCellIDs, boolList &mapFlag, label &startSeedI)
Calculate the mesh-to-mesh addressing and weights.
Macros for easy insertion into run-time selection tables.
mapNearestMethod(const mapNearestMethod &)
Disallow default bitwise copy construct.
virtual void findNearestCell(const polyMesh &mesh1, const polyMesh &mesh2, const label cell1, label &cell2) const
Find the nearest cell on mesh2 for cell1 on mesh1.
errorManip< error > abort(error &err)
const dimensionedScalar c2
Second radiation constant: default SI units: [m.K].
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
const vectorField & cellCentres() const
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
prefixOSstream Pout(cout, "Pout")
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
T remove()
Remove and return the top element.
A List with indirect addressing.
void size(const label)
Override size to be inconsistent with allocated storage.
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< scalar > magSqr(const dimensioned< Type > &)