64Foam::displacementLayeredMotionMotionSolver::getCylindrical
66 const label cellZoneI,
70 auto iter = cylSystems_.cfind(cellZoneI);
77 cylSystems_.set(cellZoneI,
new coordSystem::cylindrical(zoneDict));
79 return *cylSystems_[cellZoneI];
83void Foam::displacementLayeredMotionMotionSolver::calcZoneMask
85 const label cellZoneI,
91 isZoneEdge.resize(
mesh().nEdges());
107 for (
const label celli : cz)
109 isZonePoint.set(
mesh().cellPoints(celli));
110 isZoneEdge.set(
mesh().cellEdges(celli));
130 <<
"On cellZone " << cz.name() <<
" marked "
137void Foam::displacementLayeredMotionMotionSolver::walkStructured
139 const label cellZoneI,
140 const bitSet& isZonePoint,
153 const label pointi = seedPoints[i];
171 for (
const label pointi : isZonePoint)
187 for (
const label edgei : isZoneEdge)
207 mesh().globalData().nTotalPoints()
211 for (
const label pointi : isZonePoint)
213 const auto& pointInfo = allPointInfo[pointi];
215 distance[pointi] = pointInfo.dist();
216 data[pointi] = pointInfo.data();
219 if (patchPoints.size())
221 patchPoints[pointi] = pointInfo.index();
228Foam::tmp<Foam::vectorField>
229Foam::displacementLayeredMotionMotionSolver::faceZoneEvaluate
239 auto&
fld = tfld.ref();
243 if (
type ==
"fixedValue")
247 else if (
type ==
"timeVaryingUniformFixedValue")
251 fld = timeSeries(
mesh().time().timeOutputValue());
253 else if (
type ==
"slip")
255 if ((patchi % 2) != 1)
258 <<
"FaceZone:" << fz.name()
264 else if (
type ==
"follow")
269 else if (
type ==
"uniformFollow")
277 pointDisplacement_.boundaryField()[
patchID].patchInternalField()
284 <<
"Unknown faceZonePatch type " <<
type
285 <<
" for faceZone " << fz.name() <<
nl
293void Foam::displacementLayeredMotionMotionSolver::cellZoneSolve
295 const label cellZoneI,
299 bitSet isZonePoint, isZoneEdge;
300 calcZoneMask(cellZoneI, isZonePoint, isZoneEdge);
304 if (patchesDict.size() != 2)
307 <<
"Two faceZones (patches) must be specified per cellZone. "
308 <<
" cellZone:" << cellZoneI
309 <<
" patches:" << patchesDict.toc()
319 for (
const entry& dEntry : patchesDict)
321 const word& faceZoneName = dEntry.keyword();
326 <<
"Cannot find faceZone " << faceZoneName
359 pointDisplacement_.correctBoundaryConditions();
362 for (
const entry& dEntry : patchesDict)
364 const word& faceZoneName = dEntry.keyword();
365 const dictionary& faceZoneDict = dEntry.dict();
369 const labelList& fzMeshPoints = fz().meshPoints();
373 if (isZonePoint[fzMeshPoints[i]])
375 meshPoints.append(fzMeshPoints[i]);
395 <<
"For cellZone:" << cellZoneI
396 <<
" for faceZone:" << fz.name()
397 <<
" nPoints:" << tseed().size()
398 <<
" have patchField:"
399 <<
" min:" <<
limits.min()
400 <<
" max:" <<
limits.max()
421 patchDisp[patchi].correctBoundaryConditions();
437 mesh().cellZones()[cellZoneI].
name(),
446 for (
const label pointi : isZonePoint)
448 const scalar d1 = patchDist[0][pointi];
449 const scalar d2 = patchDist[1][pointi];
463 const word interpolationScheme(zoneDict.get<
word>(
"interpolationScheme"));
465 if (interpolationScheme ==
"oneSided")
467 for (
const label pointi : isZonePoint)
469 pointDisplacement_[pointi] = patchDisp[0][pointi];
472 else if (interpolationScheme ==
"linear")
474 for (
const label pointi : isZonePoint)
476 const scalar d1 = patchDist[0][pointi];
477 const scalar d2 = patchDist[1][pointi];
478 const scalar
s = d1/(d1 + d2 + VSMALL);
480 const vector& pd1 = patchDisp[0][pointi];
481 const vector& pd2 = patchDisp[1][pointi];
483 pointDisplacement_[pointi] = (1 -
s)*pd1 +
s*pd2;
486 else if (interpolationScheme ==
"cylindrical")
489 this->getCylindrical(cellZoneI, zoneDict);
495 <<
"cylindrical interpolation not yet available" <<
nl
496 <<
"coordinate system " << cs <<
nl
502 <<
"Invalid interpolationScheme: " << interpolationScheme
503 <<
". Valid schemes: (oneSided linear cylindrical)" <<
nl
511Foam::displacementLayeredMotionMotionSolver::
512displacementLayeredMotionMotionSolver
514 const polyMesh&
mesh,
522Foam::displacementLayeredMotionMotionSolver::
523displacementLayeredMotionMotionSolver
555 pointDisplacement_.boundaryFieldRef().updateCoeffs();
558 for (
const entry& dEntry : coeffDict().subDict(
"regions"))
560 const word& cellZoneName = dEntry.keyword();
565 Info<<
"solving for zone: " << cellZoneName <<
endl;
570 <<
"Cannot find cellZone " << cellZoneName
575 cellZoneSolve(zoneI, regionDict);
592 const vectorField displacement(this->newPoints() - points0_);
596 const label oldPointi = mpm.
pointMap()[pointi];
602 if ((masterPointi != pointi))
608 points0_[pointi] -= displacement[pointi];
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
label size() const noexcept
The number of elements in list.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
static const char *const typeName
Typename for Field.
static tmp< GeometricField< scalar, pointPatchField, pointMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=pointPatchField< scalar >::calculatedType())
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_REGISTER
Do not request registration (bool: false).
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
static FOAM_NO_DANGLING_REFERENCE const pointMesh & New(const polyMesh &mesh, Args &&... args)
Wave propagation of information through grid. Every iteration information goes through one layer of e...
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
wordList names() const
A list of the zone names.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
A cylindrical coordinate system (r-theta-z). The coordinate system angle theta is always in radians.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Mesh motion solver for an (multi-block) extruded fvMesh. Gets given the structure of the mesh blocks ...
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
virtual void updateMesh(const mapPolyMesh &mpm)
Update topology.
virtual void solve()
Solve for motion.
Virtual base class for displacement motion solver.
pointVectorField & pointDisplacement() noexcept
Return reference to the point motion displacement field.
displacementMotionSolver(const displacementMotionSolver &)=delete
No copy construct.
pointVectorField pointDisplacement_
Point motion field.
A keyword and a list of tokens is an 'entry'.
A subset of mesh faces organised as a primitive patch.
const Time & time() const
Return the top-level database.
const word & name() const
Return reference to name.
An interpolation/look-up table of scalar vs <Type> values. The reference scalar values must be monoto...
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const labelList & reversePointMap() const noexcept
Reverse point map.
const labelList & pointMap() const noexcept
Old point map.
Virtual base class for mesh motion solver.
const polyMesh & mesh() const
Return reference to mesh.
virtual tmp< pointField > newPoints()
Provide new points for motion. Solves for motion.
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Application of (multi-)patch point constraints.
void constrainDisplacement(pointVectorField &displacement, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints),.
Determines length of string of edges walked to point.
pointField & points0() noexcept
Return reference to the reference ('0') pointField.
pointIOField points0_
Starting points.
virtual void movePoints(const pointField &)
Update local data for geometry changes.
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
Mesh consisting of general polyhedral cells.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
A class for managing temporary objects.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
#define DebugInfo
Report an information message using Foam::Info.
Namespace for handling debugging switches.
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
scalar distance(const vector &p1, const vector &p2)
List< label > labelList
A List of labels.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
vectorIOField pointIOField
pointIOField is a vectorIOField.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Field< vector > vectorField
Specialisation of Field<T> for vector.
Field< label > labelField
Specialisation of Field<T> for label.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))