Go to the documentation of this file.
58 Info<<
"Updating vertex markup. curLeft: "
59 << curLeft <<
" curRight: " << curRight <<
endl;
66 if (
p[pI].
x() < curLeft - SMALL)
68 vertexMarkup[pI] = -1;
70 else if (
p[pI].
x() > curRight + SMALL)
96 Info<<
"void movingConeTopoFvMesh::addZonesAndModifiers() : "
97 <<
"Zones and modifiers already present. Skipping."
104 <<
"Adding zones and modifiers to the mesh" <<
endl;
110 boolList flipZone1(fc.size(),
false);
111 label nZoneFaces1 = 0;
114 boolList flipZone2(fc.size(),
false);
115 label nZoneFaces2 = 0;
121 fc[faceI].
x() > -0.003501
122 && fc[faceI].
x() < -0.003499
125 if ((fa[faceI] &
vector(1, 0, 0)) < 0)
127 flipZone1[nZoneFaces1] =
true;
130 zone1[nZoneFaces1] = faceI;
131 Info<<
"face " << faceI <<
" for zone 1. Flip: "
132 << flipZone1[nZoneFaces1] <<
endl;
137 fc[faceI].
x() > -0.00701
138 && fc[faceI].
x() < -0.00699
141 zone2[nZoneFaces2] = faceI;
143 if ((fa[faceI] &
vector(1, 0, 0)) > 0)
145 flipZone2[nZoneFaces2] =
true;
148 Info<<
"face " << faceI <<
" for zone 2. Flip: "
149 << flipZone2[nZoneFaces2] <<
endl;
155 flipZone1.
setSize(nZoneFaces1);
158 flipZone2.
setSize(nZoneFaces2);
172 "rightExtrusionFaces",
183 "leftExtrusionFaces",
193 Info<<
"Adding mesh zones." <<
endl;
208 "rightExtrusionFaces",
225 "leftExtrusionFaces",
238 Info<<
"Adding " << nMods <<
" mesh modifiers" <<
endl;
264 ).subDict(typeName +
"Coeffs")
266 motionVelAmplitude_(motionDict_.
lookup(
"motionVelAmplitude")),
271 Foam::
sin(time().value()*M_PI/motionVelPeriod_)
287 faceZones().findZoneID(
"leftExtrusionFaces")
295 faceZones().findZoneID(
"rightExtrusionFaces")
327 Foam::sin(time().value()*M_PI/motionVelPeriod_);
329 Pout<<
"time:" << time().value() <<
" curMotionVel_:" << curMotionVel_
330 <<
" curLeft:" << curLeft_ <<
" curRight:" << curRight_
333 if (topoChangeMap.
valid())
335 Info<<
"Topology change. Calculating motion points" <<
endl;
337 if (topoChangeMap().hasMotionPoints())
339 Info<<
"Topology change. Has premotion points" <<
endl;
386 topoChangeMap().preMotionPoints(),
393 topoChangeMap().preMotionPoints()
395 pos(0.5 -
mag(motionMask_))
396 )*curMotionVel_*time().deltaT().value();
400 Info<<
"Topology change. Already set mesh points" <<
endl;
414 pos(0.5 -
mag(motionMask_))
415 )*curMotionVel_*time().deltaT().value();
420 Info<<
"No topology change" <<
endl;
425 pos(0.5 -
mag(motionMask_))
426 )*curMotionVel_*time().deltaT().value();
430 Info <<
"Executing mesh motion" <<
endl;
431 movePoints(newPoints);
438 faceZones().findZoneID(
"leftExtrusionFaces")
446 faceZones().findZoneID(
"rightExtrusionFaces")
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
virtual const pointField & points() const
Return raw points.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
virtual ~movingConeTopoFvMesh()
Destructor.
Cell layer addition mesh modifier.
movingConeTopoFvMesh(const movingConeTopoFvMesh &)
Disallow default bitwise copy construct.
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
vector curMotionVel_
Motion velocity period.
dimensionedScalar sin(const dimensionedScalar &ds)
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Type & value() const
Return const reference to value.
void addTopologyModifiers(const List< polyMeshModifier * > &tm)
Add given set of topology modifiers to the topoChanger.
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual bool update()
Update the mesh for both mesh motion and topology change.
Constant dispersed-phase particle diameter model.
const faceZoneMesh & faceZones() const
Return face zone mesh.
virtual bool write() const
Write mesh using IO settings from time.
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.
A subset of mesh faces organised as a primitive patch.
const pointZoneMesh & pointZones() const
Return point zone mesh.
scalarField motionMask_
Vertex motion mask.
label size() const
Return number of elements in table.
scalar curLeft_
Current left obstacle position.
Macros for easy insertion into run-time selection tables.
void addZonesAndModifiers()
Add mixer zones and modifiers.
Vector< scalar > vector
A scalar version of the templated Vector.
void setSize(const label)
Reset size of List.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
tmp< scalarField > vertexMarkup(const pointField &p, const scalar curLeft, const scalar curRight) const
Markup motion vertices.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
prefixOSstream Pout(cout, "Pout")
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
volScalarField scalarField(fieldObject, mesh)
const vectorField & faceCentres() const
Abstract base class for a topology changing fvMesh.
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set).
scalar curRight_
Current right obstacle position.
label size() const
Return the number of elements in the PtrList.
const Time & time() const
Return the top-level database.
dictionary motionDict_
Motion dictionary.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
defineTypeNameAndDebug(combustionModel, 0)
void addZones(const List< pointZone * > &pz, const List< faceZone * > &fz, const List< cellZone * > &cz)
Add mesh zones.
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
polyTopoChanger topoChanger_
stressControl lookup("compactNormalStress") >> compactNormalStress
dimensionedScalar pos(const dimensionedScalar &ds)
const vectorField & faceAreas() const