Go to the documentation of this file.
70 Info<<
"Updating vertex markup. curLeft: "
71 << curLeft <<
" curRight: " << curRight <<
endl;
73 tmp<scalarField> tvertexMarkup(
new scalarField(
p.size()));
78 if (
p[pI].
x() < curLeft - SMALL)
80 vertexMarkup[pI] = -1;
82 else if (
p[pI].
x() > curRight + SMALL)
96 void Foam::movingConeTopoFvMesh::addZonesAndModifiers()
103 || faceZones().size()
104 || cellZones().size()
105 || topoChanger_.size()
109 <<
"Zones and modifiers already present. Skipping."
115 Info<<
"Time = " << time().timeName() <<
endl
116 <<
"Adding zones and modifiers to the mesh" <<
endl;
122 boolList flipZone1(fc.size(),
false);
123 label nZoneFaces1 = 0;
126 boolList flipZone2(fc.size(),
false);
127 label nZoneFaces2 = 0;
133 fc[facei].
x() > -0.003501
134 && fc[facei].
x() < -0.003499
137 if ((fa[facei] &
vector(1, 0, 0)) < 0)
139 flipZone1[nZoneFaces1] =
true;
142 zone1[nZoneFaces1] = facei;
143 Info<<
"face " << facei <<
" for zone 1. Flip: "
144 << flipZone1[nZoneFaces1] <<
endl;
149 fc[facei].
x() > -0.00701
150 && fc[facei].
x() < -0.00699
153 zone2[nZoneFaces2] = facei;
155 if ((fa[facei] &
vector(1, 0, 0)) > 0)
157 flipZone2[nZoneFaces2] =
true;
160 Info<<
"face " << facei <<
" for zone 2. Flip: "
161 << flipZone2[nZoneFaces2] <<
endl;
166 zone1.setSize(nZoneFaces1);
167 flipZone1.setSize(nZoneFaces1);
169 zone2.setSize(nZoneFaces2);
170 flipZone2.setSize(nZoneFaces2);
175 List<pointZone*> pz(0);
176 List<faceZone*> fz(2);
177 List<cellZone*> cz(0);
184 "rightExtrusionFaces",
186 std::move(flipZone1),
195 "leftExtrusionFaces",
197 std::move(flipZone2),
205 Info<<
"Adding mesh zones." <<
endl;
206 addZones(pz, fz, cz);
211 List<polyMeshModifier*> tm(2);
215 new layerAdditionRemoval
220 "rightExtrusionFaces",
221 motionDict_.subDict(
"right").get<scalar>(
"minThickness"),
222 motionDict_.subDict(
"right").get<scalar>(
"maxThickness")
226 tm[nMods] =
new layerAdditionRemoval
231 "leftExtrusionFaces",
232 motionDict_.subDict(
"left").get<scalar>(
"minThickness"),
233 motionDict_.subDict(
"left").get<scalar>(
"maxThickness")
238 Info<<
"Adding " << nMods <<
" mesh modifiers" <<
endl;
239 topoChanger_.addTopologyModifiers(tm);
247 Foam::movingConeTopoFvMesh::movingConeTopoFvMesh
267 ).optionalSubDict(typeName +
"Coeffs")
284 motionVelAmplitude_ = motionDict_.get<
vector>(
"motionVelAmplitude");
285 motionVelPeriod_ = motionDict_.get<scalar>(
"motionVelPeriod");
287 motionVelAmplitude_*
sin(time().value()*
pi/motionVelPeriod_);
288 leftEdge_ = motionDict_.get<scalar>(
"leftEdge");
289 curLeft_ = motionDict_.get<scalar>(
"leftObstacleEdge");
290 curRight_ = motionDict_.get<scalar>(
"rightObstacleEdge");
292 Pout<<
"Initial time:" << time().value()
293 <<
" Initial curMotionVel_:" << curMotionVel_
296 addZonesAndModifiers();
302 faceZones().findZoneID(
"leftExtrusionFaces")
310 faceZones().findZoneID(
"rightExtrusionFaces")
314 motionMask_ = vertexMarkup
344 motionVelAmplitude_*
sin(time().value()*
pi/motionVelPeriod_);
346 Pout<<
"time:" << time().value() <<
" curMotionVel_:" << curMotionVel_
347 <<
" curLeft:" << curLeft_ <<
" curRight:" << curRight_
352 Info<<
"Topology change. Calculating motion points" <<
endl;
354 if (topoChangeMap().hasMotionPoints())
356 Info<<
"Topology change. Has premotion points" <<
endl;
361 topoChangeMap().preMotionPoints(),
368 topoChangeMap().preMotionPoints()
371 )*curMotionVel_*time().deltaT().value();
375 Info<<
"Topology change. Already set mesh points" <<
endl;
390 )*curMotionVel_*time().deltaT().value();
395 Info<<
"No topology change" <<
endl;
401 )*curMotionVel_*time().deltaT().value();
405 Info <<
"Executing mesh motion" <<
endl;
406 movePoints(newPoints);
414 faceZones().findZoneID(
"leftExtrusionFaces")
422 faceZones().findZoneID(
"rightExtrusionFaces")
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
List< label > labelList
A List of labels.
vectorField pointField
pointField is a vectorField.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
#define InfoInFunction
Report an information message using Foam::Info.
virtual ~movingConeTopoFvMesh()
Destructor.
A class for managing temporary objects.
dimensionedScalar sin(const dimensionedScalar &ds)
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
List< bool > boolList
A List of bools.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar pos0(const dimensionedScalar &ds)
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
virtual bool update()
Update the mesh for both mesh motion and topology change.
#define forAll(list, i)
Loop across all elements in list.
Field< vector > vectorField
Specialisation of Field<T> for vector.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Macros for easy insertion into run-time selection tables.
Vector< scalar > vector
A scalar version of the templated Vector.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
constexpr scalar pi(M_PI)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Abstract base class for a topology changing fvMesh.
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
constant condensation/saturation model.
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)