Go to the documentation of this file.
70 Pout<<
"participating source mesh cells: " <<
cells.size() <<
endl;
83 scalar threshold = tolerance_*src_.cellVolumes()[srcCellI];
87 treeBoundBox bbTgtCell(tgt_.points(), tgt_.cellPoints()[tgtCellI]);
103 const label srcCellI,
109 treeBoundBox bbTgtCell(tgt_.points(), tgt_.cellPoints()[tgtCellI]);
127 const label srcCellI,
133 treeBoundBox bbTgtCell(tgt_.points(), tgt_.cellPoints()[tgtCellI]);
146 if (volAndInertia.
first() <= ROOTVSMALL)
148 volAndInertia.
first() = 0.0;
156 return volAndInertia;
173 label nbrCellI = nbrCells[i];
177 (
findIndex(visitedCells, nbrCellI) == -1)
178 && (
findIndex(nbrCellIDs, nbrCellI) == -1)
181 nbrCellIDs.
append(nbrCellI);
195 srcToTgtAddr.
setSize(src_.nCells());
196 srcToTgtWght.
setSize(src_.nCells());
197 tgtToSrcAddr.
setSize(tgt_.nCells());
198 tgtToSrcWght.
setSize(tgt_.nCells());
204 else if (!tgt_.nCells())
208 Pout<<
"mesh interpolation: hhave " << src_.nCells() <<
" source "
209 <<
" cells but no target cells" <<
endl;
231 if (!src_.nCells() || !tgt_.nCells())
235 Pout<<
"mesh interpolation: cells not on processor: Source cells = "
236 << src_.nCells() <<
", target cells = " << tgt_.nCells()
261 word fName(
"addressing_" + mesh1.
name() +
"_to_" + mesh2.
name());
268 OFstream os(src_.time().path()/fName +
".obj");
271 forAll(mesh1ToMesh2Addr, i)
273 const labelList& addr = mesh1ToMesh2Addr[i];
276 label cellI = addr[j];
284 os <<
"v " <<
p.x() <<
' ' <<
p.y() <<
' ' <<
p.z() <<
nl;
286 os <<
"v " << c0.
x() <<
' ' << c0.
y() <<
' ' << c0.
z()
289 os <<
"l " << vertI - 1 <<
' ' << vertI <<
nl;
virtual const pointField & points() const
Return raw points.
A class for handling words, derived from string.
const point & max() const
Maximum describing the bounding box.
#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,.
virtual ~meshToMeshMethod()
Destructor.
const Type2 & second() const
Return second.
scalar cellCellOverlapVolumeMinDecomp(const primitiveMesh &meshA, const label cellAI, const primitiveMesh &meshB, const label cellBI, const treeBoundBox &cellBbB) const
Calculates the overlap volume.
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
void inflate(const scalar s)
Inflate box by factor*mag(span) in all dimensions.
Standard boundBox + extra functionality for use in octree.
virtual Tuple2< scalar, point > interVolAndCentroid(const label srcCellI, const label tgtCellI)
Return the intersection volume and centroid between two cells.
static bool & parRun()
Is this a parallel run?
const Type1 & first() const
Return first.
const cellList & cells() const
virtual scalar interVol(const label srcCellI, const label tgtCellI) const
Return the intersection volume between two cells.
static scalar tolerance_
Tolerance used in volume overlap calculations.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Mesh consisting of general polyhedral cells.
void writeConnectivity(const polyMesh &mesh1, const polyMesh &mesh2, const labelListList &mesh1ToMesh2Addr) const
Write the connectivity (debugging)
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.
bool overlaps(const boundBox &) const
Overlaps/touches boundingBox?
const polyMesh & tgt_
Reference to the target mesh.
const labelListList & cellCells() const
const word & name() const
Return name.
const point & min() const
Minimum describing the bounding box.
static void cellCellOverlapMinDecomp(const primitiveMesh &meshA, const label cellAI, const primitiveMesh &meshB, const label cellBI, const treeBoundBox &cellBbB, tetsOp &combineTetsOp)
Cell overlap calculation.
const polyMesh & src_
Reference to the source mesh.
Tuple2< scalar, point > cellCellOverlapMomentMinDecomp(const primitiveMesh &meshA, const label cellAI, const primitiveMesh &meshB, const label cellBI, const treeBoundBox &cellBbB) const
Calculates the overlap volume and moment.
virtual void appendNbrCells(const label tgtCellI, const polyMesh &mesh, const DynamicList< label > &visitedTgtCells, DynamicList< label > &nbrTgtCellIDs) const
Append target cell neihgbour cells to cellIDs list.
labelList maskCells() const
Return src cell IDs for the overlap region.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
const boundBox & bounds() const
Return mesh bounding box.
void setSize(const label)
Reset size of List.
const vectorField & cellCentres() const
virtual const faceList & faces() const
Return raw faces.
virtual bool initialise(labelListList &srcToTgtAddr, scalarListList &srcToTgtWght, labelListList &tgtToTgtAddr, scalarListList &tgtToTgtWght) const
prefixOSstream Pout(cout, "Pout")
meshToMeshMethod(const meshToMeshMethod &)
Disallow default bitwise copy construct.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
const dimensionedScalar e
Elementary charge.
A bounding box defined in terms of the points at its extremities.
const dimensionedScalar c
Speed of light in a vacuum.
A 2-tuple for storing two objects of different types.
Calculates the overlap volume of two cells using tetrahedral decomposition.
virtual bool intersect(const label srcCellI, const label tgtCellI) const
Return the true if cells intersect.
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
A cell is defined as a list of faces with extra functionality.
word name(const complex &)
Return a string representation of a complex.