48void Foam::linearValveFvMesh::addZonesAndModifiers()
61 <<
"Zones and modifiers already present. Skipping."
68 <<
"Adding zones and modifiers to the mesh" <<
endl;
71 List<pointZone*> pz(1);
74 pz[0] =
new pointZone(
"cutPointZone", 0,
pointZones());
79 List<faceZone*> fz(3);
82 const word innerSliderName
84 motionDict_.subDict(
"slider").get<word>(
"inside")
86 const polyPatch& innerSlider =
boundaryMesh()[innerSliderName];
98 const word outerSliderName
100 motionDict_.subDict(
"slider").get<word>(
"outside")
102 const polyPatch& outerSlider =
boundaryMesh()[outerSliderName];
114 fz[2] =
new faceZone(
"cutFaceZone", 2,
faceZones());
116 List<cellZone*> cz(0);
118 Info<<
"Adding point, face and cell zones" <<
endl;
122 Info<<
"Adding topology modifiers" <<
endl;
132 outerSliderName +
"Zone",
133 innerSliderName +
"Zone",
149void Foam::linearValveFvMesh::makeSlidersDead()
158 topoChanges[modI].disable();
163 <<
"Don't know what to do with mesh modifier "
164 << modI <<
" of type " << topoChanges[modI].type()
171void Foam::linearValveFvMesh::makeSlidersLive()
180 topoChanges[modI].enable();
185 <<
"Don't know what to do with mesh modifier "
186 << modI <<
" of type " << topoChanges[modI].type()
193bool Foam::linearValveFvMesh::attached()
const
222 <<
" named " << topoChanges[modI].name()
223 <<
" out of sync: Should be" << result
231 Info<<
"linearValveFvMesh: attached!" <<
endl;
235 Info<<
"linearValveFvMesh: detached!" <<
endl;
245Foam::linearValveFvMesh::linearValveFvMesh(
const IOobject&
io)
259 ).optionalSubDict(
typeName +
"Coeffs")
263 addZonesAndModifiers();
280 Info<<
"Decoupling sliding interfaces" <<
endl;
288 msPtr_->updateMesh();
292 Info<<
"Sliding interfaces decoupled" <<
endl;
300 setMorphTimeIndex(3*time().
timeIndex() + 1);
303 msPtr_->updateMesh();
307 if (topoChangeMap().hasMotionPoints())
309 Info<<
"Topology change; executing pre-motion" <<
endl;
310 movePoints(topoChangeMap().preMotionPoints());
317 movePoints(msPtr_->curPoints());
320 Info<<
"Coupling sliding interfaces" <<
endl;
323 setMorphTimeIndex(3*time().
timeIndex() + 2);
326 Info<<
"Moving points post slider attach" <<
endl;
328 msPtr_->updateMesh();
330 Info<<
"Sliding interfaces coupled: " << attached() <<
endl;
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
label size() const noexcept
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
static dictionary readContents(const IOobject &io)
Read and return contents, testing for "dictionary" type. The IOobject will not be registered.
@ MUST_READ
Reading required.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
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 size() const noexcept
The number of entries in the list.
const dictionary & optionalSubDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary, otherwise return this dictionary.
static autoPtr< dynamicFvMesh > New(const IOobject &io)
Select, construct and return the dynamicFvMesh.
const Time & time() const
Return the top-level database.
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
virtual void movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
virtual ~linearValveFvMesh()
Destructor.
virtual bool update()
Update the mesh for both mesh motion and topology change.
Virtual base class for mesh motion solver.
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.
void addZones(PtrList< pointZone > &&pz, PtrList< faceZone > &&fz, PtrList< cellZone > &&cz)
Add mesh zones.
const pointZoneMesh & pointZones() const noexcept
Return point zone mesh.
List of mesh modifiers defining the mesh dynamics.
Abstract base class for a topology changing fvMesh.
polyTopoChanger topoChanger_
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
#define InfoInFunction
Report an information message using Foam::Info.
Different types of constants.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
messageStream Info
Information stream (stdout output on master, null elsewhere).
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
errorManip< error > abort(error &err)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
#define forAll(list, i)
Loop across all elements in list.