55 Foam::directions::directionTypeNames_;
60 void Foam::directions::writeOBJ(Ostream& os,
const point& pt)
62 os <<
"v " << pt.x() <<
' ' << pt.y() <<
' ' << pt.z() <<
endl;
66 void Foam::directions::writeOBJ
77 os <<
"l " << vertI + 1 <<
' ' << vertI + 2 <<
endl;
83 void Foam::directions::writeOBJ
85 const fileName& fName,
86 const primitiveMesh& mesh,
90 Pout<<
"Writing cell info to " << fName <<
" as vectors at the cellCentres"
93 OFstream xDirStream(fName);
99 const point& ctr = mesh.cellCentres()[celli];
104 const labelList& nbrs = mesh.cellCells()[celli];
113 writeOBJ(xDirStream, ctr, ctr + scale*dirs[celli], vertI);
118 void Foam::directions::check2D
120 const twoDPointCorrector* correct2DPtr,
126 if (
mag(correct2DPtr->planeNormal() & vec) > 1
e-6)
129 <<
"is not normal to plane defined in dynamicMeshDict."
131 <<
"Either make case 3D or adjust vector."
140 const polyMesh& mesh,
149 List<directionInfo> changedFacesInfo(pp.size());
155 label meshFacei = pp.start() + patchFacei;
157 label celli = mesh.faceOwner()[meshFacei];
159 if (!hexMatcher().
isA(mesh, celli))
162 <<
"useHexTopology specified but cell " << celli
163 <<
" on face " << patchFacei <<
" of patch " << pp.name()
167 const vector& cutDir = ppField[patchFacei];
183 changedFaces[patchFacei] = meshFacei;
184 changedFacesInfo[patchFacei] =
196 changedFaces[patchFacei] = pp.start() + patchFacei;
197 changedFacesInfo[patchFacei] =
206 List<directionInfo> faceInfo(mesh.nFaces()), cellInfo(mesh.nCells());
207 FaceCellWave<directionInfo> directionCalc
214 mesh.globalData().nTotalCells() + 1
225 label index = cellInfo[celli].index();
231 <<
"Cell " << celli <<
" never visited to determine "
232 <<
"local coordinate system" <<
endl
233 <<
"Using direction " << defaultDir <<
" instead" <<
endl;
235 dirField[celli] = defaultDir;
239 else if (index == -2)
242 dirField[celli] = cellInfo[celli].n();
246 else if (index == -1)
249 <<
"Illegal index " << index <<
endl
261 reduce(nGeom, sumOp<label>());
262 reduce(nTopo, sumOp<label>());
263 reduce(nUnset, sumOp<label>());
265 Info<<
"Calculated local coords for " << defaultDir
267 <<
" Geometric cut cells : " << nGeom <<
endl
268 <<
" Topological cut cells : " << nTopo <<
endl
269 <<
" Unset cells : " << nUnset <<
endl
288 const word coordSystem(
dict.lookup(
"coordinateSystem"));
295 if (coordSystem !=
"fieldBased")
299 directionType wantedDir = directionTypeNames_[wantedDirs[i]];
305 else if (wantedDir ==
e1)
309 else if (wantedDir ==
e2)
317 if (coordSystem ==
"global")
322 check2D(correct2DPtr,
e1);
325 check2D(correct2DPtr,
e2);
330 Info<<
"Global Coordinate system:" <<
endl
349 else if (coordSystem ==
"patchLocal")
353 const word patchName(patchDict.
lookup(
"patch"));
360 <<
"Cannot find patch "
377 <<
"Discarding user specified e1 since 2D case." <<
endl
378 <<
"Recalculated e1 from face normal and planeNormal as "
382 Switch useTopo(
dict.lookup(
"useHexTopology"));
387 if (wantE3 || wantE2)
405 if (wantE1 || wantE2)
430 else if (coordSystem ==
"fieldBased")
450 <<
"Unknown coordinate system "
451 << coordSystem <<
endl
452 <<
"Known types are global, patchLocal and fieldBased"
#define forAll(list, i)
Loop across all elements in list.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
fileName & instance() const
Return the instance directory, constant, system, <time> etc.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Initialise the NamedEnum HashTable from the static list of names.
const Field< PointType > & faceNormals() const
Return face normals for patch.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
label size() const
Return the number of elements in the UList.
T & operator[](const label)
Return element of UList.
A list of keyword definitions, which are a keyword followed by any number of values (e....
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
static label edgeToFaceIndex(const primitiveMesh &mesh, const label celli, const label facei, const label edgeI)
Given edge on hex cell find corresponding edge on face. Is either.
directions(const polyMesh &mesh, const dictionary &dict, const twoDPointCorrector *correct2DPtr=nullptr)
Construct from mesh and dictionary and optional 2D corrector.
directionType
Enumeration listing the possible coordinate directions.
label findIndex(const word &patchName) const
Find patch index given a name.
Mesh consisting of general polyhedral cells.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
A patch is a list of labels that address the faces in the global face list.
A class for managing temporary objects.
Class applies a two-dimensional correction to mesh motion point field.
const vector & planeNormal() const
Return plane normal.
A class for handling words, derived from string.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
#define WarningInFunction
Report a warning using Foam::Warning.
errorManipArg< error, int > exit(error &err, const int errNo=1)
scalar minDist(const List< pointIndexHit > &hitList)
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.
errorManip< error > abort(error &err)
vector point
Point is a vector.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
IOField< vector > vectorIOField
vectorField with IO.
Vector< scalar > vector
A scalar version of the templated Vector.
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
dimensioned< scalar > mag(const dimensioned< Type > &)
Field< vector > vectorField
Specialisation of Field<T> for vector.
prefixOSstream Pout(cout, "Pout")