48 bool Foam::hexCellLooper::walkHex
51 const label startFacei,
52 const label startEdgeI,
58 label facei = startFacei;
60 label edgeI = startEdgeI;
68 Pout<<
" walkHex : inserting cut onto edge:" << edgeI
74 loopWeights[cutI] = 0.5;
84 if (edgeI == startEdgeI)
94 Pout<<
"hexCellLooper::walkHex" <<
"Problem : cell:" << celli
95 <<
" collected loop:";
97 Pout<<
"loopWeights:" << loopWeights <<
endl;
109 void Foam::hexCellLooper::makeFace
118 facePoints.setSize(loop.size());
119 faceVerts.setSize(loop.size());
123 label cut = loop[cutI];
127 label edgeI = getEdge(cut);
129 const edge&
e = mesh().edges()[edgeI];
131 const point& v0 = mesh().points()[
e.start()];
132 const point& v1 = mesh().points()[
e.end()];
135 loopWeights[cutI]*v1 + (1.0-loopWeights[cutI])*v0;
139 label vertI = getVertex(cut);
141 facePoints[cutI] = mesh().points()[vertI];
144 faceVerts[cutI] = cutI;
181 if (mesh().
cellShapes()[celli].model() == hex_)
196 success = walkHex(celli, face0, edgeI, loop, loopWeights);
218 <<
"could not cut cell " << celli <<
endl;
220 fileName cutsFile(
"hexCellLooper_" +
name(celli) +
".obj");
222 Pout<<
"hexCellLooper : writing cell to " << cutsFile <<
endl;
243 label elem = loop[elemI];
245 if (loopSet.
found(elem))
257 makeFace(loop, loopWeights, faceVerts, facePoints);
259 if ((faceVerts.
mag(facePoints) < small) || (loop.
size() < 3))
262 <<
" on points:" << facePoints <<
endl
274 const plane& cutPlane,
Various functions to operate on Lists.
#define forAll(list, i)
Loop across all elements in list.
Macros for easy insertion into run-time selection tables.
bool insert(const Key &key)
Insert a new entry.
bool found(const Key &) const
Return true if hashedEntry is found in table.
void size(const label)
Override size to be inconsistent with allocated storage.
void setSize(const label)
Reset size of List.
A List with indirect addressing.
bool empty() const
Return true if the UList is empty (ie, size() is zero)
A static collection of cell models, and a means of looking them up.
Ostream & writeCuts(Ostream &os, const labelList &, const scalarField &) const
Write cut descriptions to Ostream.
static label edgeToEVert(const primitiveMesh &mesh, const label edgeI)
Convert edgeI to eVert.
const polyMesh & mesh() const
A face is a list of labels corresponding to mesh vertices.
scalar mag(const pointField &) const
Return scalar magnitude.
A class for handling file names.
Implementation of cellLooper. Does pure geometric cut through cell.
virtual bool cut(const vector &refDir, const label celli, const boolList &vertIsCut, const boolList &edgeIsCut, const scalarField &edgeWeight, labelList &loop, scalarField &loopWeights) const
Create cut along circumference of celli. Gets current mesh cuts.
Implementation of cellLooper.
hexCellLooper(const polyMesh &mesh)
Construct from components.
virtual ~hexCellLooper()
Destructor.
virtual bool cut(const vector &refDir, const label celli, const boolList &vertIsCut, const boolList &edgeIsCut, const scalarField &edgeWeight, labelList &loop, scalarField &loopWeights) const
Create cut along circumference of celli. Gets current mesh cuts.
Geometric class that creates a 2D plane and can return the intersection point between a line and the ...
Mesh consisting of general polyhedral cells.
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
const cellShapeList & cellShapes
volScalarField scalarField(fieldObject, mesh)
#define WarningInFunction
Report a warning using Foam::Warning.
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.
Ostream & endl(Ostream &os)
Add newline and flush stream.
word name(const bool)
Return a word representation of a bool.
errorManip< error > abort(error &err)
vectorField pointField
pointField is a vectorField.
vector point
Point is a vector.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
defineTypeNameAndDebug(combustionModel, 0)
prefixOSstream Pout(cout, "Pout")