61 Info<<
"Updating vertex markup. curLeft: " 62 << curLeft <<
" curRight: " << curRight <<
endl;
64 tmp<scalarField> tvertexMarkup(
new scalarField(p.size()));
69 if (p[pI].
x() < curLeft - small)
71 vertexMarkup[pI] = -1;
73 else if (p[pI].
x() > curRight + small)
87 void Foam::movingConeTopoFvMesh::addZonesAndModifiers()
96 || topoChanger_.size()
100 <<
"Zones and modifiers already present. Skipping." 106 Info<<
"Time = " << time().timeName() << endl
107 <<
"Adding zones and modifiers to the mesh" <<
endl;
113 boolList flipZone1(fc.size(),
false);
114 label nZoneFaces1 = 0;
117 boolList flipZone2(fc.size(),
false);
118 label nZoneFaces2 = 0;
124 fc[facei].
x() > -0.003501
125 && fc[facei].
x() < -0.003499
128 if ((fa[facei] &
vector(1, 0, 0)) < 0)
130 flipZone1[nZoneFaces1] =
true;
133 zone1[nZoneFaces1] = facei;
134 Info<<
"face " << facei <<
" for zone 1. Flip: " 135 << flipZone1[nZoneFaces1] <<
endl;
140 fc[facei].
x() > -0.00701
141 && fc[facei].
x() < -0.00699
144 zone2[nZoneFaces2] = facei;
146 if ((fa[facei] &
vector(1, 0, 0)) > 0)
148 flipZone2[nZoneFaces2] =
true;
151 Info<<
"face " << facei <<
" for zone 2. Flip: " 152 << flipZone2[nZoneFaces2] <<
endl;
157 zone1.setSize(nZoneFaces1);
158 flipZone1.setSize(nZoneFaces1);
160 zone2.setSize(nZoneFaces2);
161 flipZone2.setSize(nZoneFaces2);
166 List<pointZone*> pz(0);
167 List<faceZone*> fz(2);
168 List<cellZone*> cz(0);
175 "rightExtrusionFaces",
186 "leftExtrusionFaces",
196 Info<<
"Adding mesh zones." <<
endl;
197 addZones(pz, fz, cz);
202 List<polyMeshModifier*> tm(2);
206 new layerAdditionRemoval
211 "rightExtrusionFaces",
212 motionDict_.subDict(
"right").lookup<scalar>(
"minThickness"),
213 motionDict_.subDict(
"right").lookup<scalar>(
"maxThickness")
217 tm[nMods] =
new layerAdditionRemoval
222 "leftExtrusionFaces",
223 motionDict_.subDict(
"left").lookup<scalar>(
"minThickness"),
224 motionDict_.subDict(
"left").lookup<scalar>(
"maxThickness")
229 Info<<
"Adding " << nMods <<
" mesh modifiers" <<
endl;
230 topoChanger_.addTopologyModifiers(tm);
241 motionDict_(dynamicMeshDict().optionalSubDict(typeName +
"Coeffs")),
242 motionVelAmplitude_(motionDict_.
lookup(
"motionVelAmplitude")),
243 motionVelPeriod_(motionDict_.
lookup<scalar>(
"motionVelPeriod")),
246 motionVelAmplitude_*
sin(time().value()*
pi/motionVelPeriod_)
248 leftEdge_(motionDict_.
lookup<scalar>(
"leftEdge")),
249 curLeft_(motionDict_.
lookup<scalar>(
"leftObstacleEdge")),
250 curRight_(motionDict_.
lookup<scalar>(
"rightObstacleEdge"))
253 <<
" Initial curMotionVel_:" << curMotionVel_
256 addZonesAndModifiers();
262 faceZones().findZoneID(
"leftExtrusionFaces")
270 faceZones().findZoneID(
"rightExtrusionFaces")
274 motionMask_ = vertexMarkup
301 motionVelAmplitude_*
sin(
time().value()*
pi/motionVelPeriod_);
303 Pout<<
"time:" <<
time().
value() <<
" curMotionVel_:" << curMotionVel_
304 <<
" curLeft:" << curLeft_ <<
" curRight:" << curRight_
307 if (topoChangeMap.
valid())
309 Info<<
"Topology change. Calculating motion points" <<
endl;
311 if (topoChangeMap().hasMotionPoints())
313 Info<<
"Topology change. Has premotion points" <<
endl;
318 topoChangeMap().preMotionPoints(),
325 topoChangeMap().preMotionPoints()
332 Info<<
"Topology change. Already set mesh points" <<
endl;
352 Info<<
"No topology change" <<
endl;
362 Info <<
"Executing mesh motion" <<
endl;
371 faceZones().findZoneID(
"leftExtrusionFaces")
379 faceZones().findZoneID(
"rightExtrusionFaces")
#define forAll(list, i)
Loop across all elements in list.
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.
bool movePoints()
Do what is necessary if the mesh has moved.
movingConeTopoFvMesh(const IOobject &io)
Construct from database.
Vector< scalar > vector
A scalar version of the templated Vector.
virtual bool update()
Update the mesh for both mesh motion and topology change.
const Time & time() const
Return the top-level database.
autoPtr< mapPolyMesh > changeMesh(const bool inflate, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Macros for easy insertion into run-time selection tables.
virtual const pointField & points() const
Return raw points.
List< bool > boolList
Bool container classes.
vectorField pointField
pointField is a vectorField.
stressControl lookup("compactNormalStress") >> compactNormalStress
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Abstract base class for a topology changing fvMesh.
const Type & value() const
Return const reference to value.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
List< label > labelList
A List of labels.
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set)
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedScalar sin(const dimensionedScalar &ds)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
prefixOSstream Pout(cout, "Pout")
const meshFaceZones & faceZones() const
Return face zones.
dimensioned< scalar > mag(const dimensioned< Type > &)
Field< vector > vectorField
Specialisation of Field<T> for vector.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
A class for managing temporary objects.
polyTopoChanger topoChanger_
virtual ~movingConeTopoFvMesh()
Destructor.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
dimensionedScalar deltaT() const
Return time step.
#define InfoInFunction
Report an information message using Foam::Info.