70 Info<<
"Updating vertex markup. curLeft: "
71 << curLeft <<
" curRight: " << curRight <<
endl;
74 auto& vertexMarkup = tvertexMarkup.ref();
78 if (
p[pI].
x() < curLeft - SMALL)
80 vertexMarkup[pI] = -1;
82 else if (
p[pI].
x() > curRight + SMALL)
96void 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;
123 label nZoneFaces1 = 0;
127 label nZoneFaces2 = 0;
133 fc[facei].
x() > -0.003501
134 && fc[facei].
x() < -0.003499
137 if (
fa[facei].
x() < 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].
x() > 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);
184 "rightExtrusionFaces",
186 std::move(flipZone1),
195 "leftExtrusionFaces",
197 std::move(flipZone2),
205 Info<<
"Adding mesh zones." <<
endl;
206 addZones(pz, fz, cz);
220 "rightExtrusionFaces",
221 motionDict_.subDict(
"right").get<scalar>(
"minThickness"),
222 motionDict_.subDict(
"right").get<scalar>(
"maxThickness")
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);
247Foam::movingConeTopoFvMesh::movingConeTopoFvMesh
265 ).optionalSubDict(
typeName +
"Coeffs")
282 motionVelAmplitude_ = motionDict_.get<
vector>(
"motionVelAmplitude");
283 motionVelPeriod_ = motionDict_.get<scalar>(
"motionVelPeriod");
285 motionVelAmplitude_*
sin(time().value()*
pi/motionVelPeriod_);
286 leftEdge_ = motionDict_.get<scalar>(
"leftEdge");
287 curLeft_ = motionDict_.get<scalar>(
"leftObstacleEdge");
288 curRight_ = motionDict_.get<scalar>(
"rightObstacleEdge");
290 Pout<<
"Initial time:" << time().value()
291 <<
" Initial curMotionVel_:" << curMotionVel_
294 addZonesAndModifiers();
298 faceZones()[
"leftExtrusionFaces"].
patch().localPoints()
303 faceZones()[
"rightExtrusionFaces"].
patch().localPoints()
306 motionMask_ = vertexMarkup
336 motionVelAmplitude_*
sin(time().value()*
pi/motionVelPeriod_);
338 Pout<<
"time:" << time().value() <<
" curMotionVel_:" << curMotionVel_
339 <<
" curLeft:" << curLeft_ <<
" curRight:" << curRight_
344 Info<<
"Topology change. Calculating motion points" <<
endl;
346 if (topoChangeMap().hasMotionPoints())
348 Info<<
"Topology change. Has premotion points" <<
endl;
353 topoChangeMap().preMotionPoints(),
360 topoChangeMap().preMotionPoints()
363 )*curMotionVel_*time().deltaTValue();
367 Info<<
"Topology change. Already set mesh points" <<
endl;
382 )*curMotionVel_*time().deltaTValue();
387 Info<<
"No topology change" <<
endl;
393 )*curMotionVel_*time().deltaTValue();
397 Info <<
"Executing mesh motion" <<
endl;
398 movePoints(newPoints);
404 faceZones()[
"leftExtrusionFaces"].
patch().localPoints()
409 faceZones()[
"rightExtrusionFaces"].
patch().localPoints()
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
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.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
void size(const label n)
Older name for setAddressableSize.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
const dictionary & optionalSubDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary, otherwise return this dictionary.
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
A subset of mesh faces organised as a primitive patch.
const Time & time() const
Return the top-level database.
virtual void movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Cell layer addition mesh modifier.
Sample topoChangerFvMesh that moves an object in x direction and introduces/removes layers.
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
virtual bool update()
Update the mesh for both mesh motion and topology change.
virtual ~movingConeTopoFvMesh()
Destructor.
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
A class for managing temporary objects.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Abstract base class for a topology changing fvMesh.
polyTopoChanger topoChanger_
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define InfoInFunction
Report an information message using Foam::Info.
constexpr scalar pi(M_PI)
Different types of constants.
Namespace for finite-area.
const std::string patch
OpenFOAM patch number as a std::string.
dimensionedScalar pos0(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
dimensionedScalar sin(const dimensionedScalar &ds)
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.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< bool > boolList
A List of bools.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &f1, const label comm)
#define forAll(list, i)
Loop across all elements in list.