Coil Class Reference

This class defines a coil in the z-direction from a (constant) starting point, a (constant) length L, a (constant) radius r, a (constant) number or revelations N and a (constant) rotation speed (rev/s) More...

#include <Coil.h>

+ Inheritance diagram for Coil:

Public Member Functions

 Coil ()
 Default constructor, sets a coil with default parameters. More...
 
 Coil (const Coil &other)
 Copy constructor, makes a coil with the same properties as the input Coil. More...
 
 Coil (Vec3D Start, Mdouble L, Mdouble r, Mdouble N, Mdouble omega, Mdouble thickness)
 Constructor in which all parameters are set. More...
 
 ~Coil () override
 Default destructor. More...
 
void set (Vec3D Start, Mdouble length, Mdouble radius, Mdouble numberOfRevelations, Mdouble omega, Mdouble thickness)
 Set all parameters of this Coil. More...
 
Coilcopy () const override
 Copy this Coil and return a pointer to the copy, useful for polymorphism. More...
 
bool getDistanceAndNormal (const BaseParticle &P, Mdouble &distance, Vec3D &normal_return) const override
 Compute the distance from the Coil for a given BaseParticle and return if there is a collision. If there is a collision, also return the normal vector of the interaction point. More...
 
void move_time (Mdouble dt)
 Rotate the Coil for a period dt, so that the offset_ changes with omega_*dt. More...
 
void read (std::istream &is) override
 Reads a Coil from an input stream, for example a restart file. More...
 
MERCURYDPM_DEPRECATED void oldRead (std::istream &is)
 Reads an old-style Coil from an input stream, for example an old restart file. More...
 
void write (std::ostream &os) const override
 Writes a Coil to an output stream, for example a restart file. More...
 
std::string getName () const override
 Returns the name of the object, in this case the string "Coil". More...
 
- Public Member Functions inherited from BaseWall
 BaseWall ()
 Default constructor. More...
 
 BaseWall (const BaseWall &w)
 Copy constructor. More...
 
 ~BaseWall () override
 Default destructor. More...
 
virtual bool getDistanceNormalOverlap (const BaseParticle &P, Mdouble &distance, Vec3D &normal_return, Mdouble &overlap) const
 
virtual bool getDistanceNormalOverlapSuperquadric (const SuperQuadricParticle &p, Mdouble &distance, Vec3D &normal_return, Mdouble &overlap) const
 
virtual Vec3D getFurthestPointSuperQuadric (const Vec3D &normalBodyFixed, const Vec3D &axes, Mdouble eps1, Mdouble eps2) const
 
virtual void setHandler (WallHandler *handler)
 A function which sets the WallHandler for this BaseWall. More...
 
WallHandlergetHandler () const
 A function which returns the WallHandler that handles this BaseWall. More...
 
void setIndSpecies (unsigned int indSpecies) override
 Define the species of this wall using the index of the species in the SpeciesHandler in this DPMBase. More...
 
void setSpecies (const ParticleSpecies *species)
 Defines the species of the current wall. More...
 
bool isFixed () const override
 
void setForceControl (Vec3D forceGoal, Vec3D gainFactor, Vec3D baseVelocity={0, 0, 0})
 Slowly adjusts the force on a wall towards a specified goal, by adjusting (prescribing) the velocity of the wall. More...
 
virtual bool isLocal (Vec3D &min, Vec3D &max) const
 
bool getLinePlaneIntersect (Vec3D &intersect, const Vec3D &p0, const Vec3D &p1, const Vec3D &n, const Vec3D &p)
 
bool isInsideWallVTK (const Vec3D &point, const Vec3D &normal, const Vec3D &position) const
 
void projectOntoWallVTK (Vec3D &point0, const Vec3D &point1, const Vec3D &normal, const Vec3D &position) const
 
void intersectVTK (std::vector< Vec3D > &points, Vec3D normal, Vec3D position) const
 
virtual BaseInteractiongetInteractionWithSuperQuad (SuperQuadricParticle *p, unsigned timeStamp, InteractionHandler *interactionHandler)
 
virtual void writeVTK (VTKContainer &vtk) const
 
void getVTK (std::vector< Vec3D > &points, std::vector< std::vector< double >> &triangleStrips)
 
const Vec3D getAxis () const
 
BaseInteractiongetInteractionWith (BaseParticle *p, unsigned timeStamp, InteractionHandler *interactionHandler) override
 Returns the interaction between this wall and a given particle, nullptr if there is no interaction. More...
 
virtual void actionsOnRestart ()
 No implementation but can be overidden in its derived classes. More...
 
virtual void actionsAfterParticleGhostUpdate ()
 No implementation but can be overidden in its derived classes. More...
 
virtual void handleParticleAddition (unsigned int id, BaseParticle *p)
 Handles the addition of particles to the particleHandler. More...
 
virtual void handleParticleRemoval (unsigned int id)
 Handles the addition of particles to the particleHandler. More...
 
virtual void checkInteractions (InteractionHandler *interactionHandler, unsigned int timeStamp)
 Check if all interactions are valid. More...
 
bool getVTKVisibility () const
 
void setVTKVisibility (bool vtkVisibility)
 
void addRenderedWall (BaseWall *w)
 
BaseWallgetRenderedWall (size_t i) const
 
std::vector< BaseWall * > getRenderedWalls () const
 
void removeRenderedWalls ()
 
void renderWall (VTKContainer &vtk)
 
void addParticlesAtWall (unsigned numElements=50)
 
void setVelocityControl (Vec3D forceGoal, Vec3D gainFactor, Vec3D baseVelocity)
 
virtual void writeWallDetailsVTK (VTKData &data) const
 
virtual void computeWear ()
 
- Public Member Functions inherited from BaseInteractable
 BaseInteractable ()
 Default BaseInteractable constructor. More...
 
 BaseInteractable (const BaseInteractable &p)
 Copy constructor. More...
 
 ~BaseInteractable () override
 Destructor, it simply destructs the BaseInteractable and all the objects it contains. More...
 
unsigned int getIndSpecies () const
 Returns the index of the species associated with the interactable object. More...
 
const ParticleSpeciesgetSpecies () const
 Returns a pointer to the species of this BaseInteractable. More...
 
void setSpecies (const ParticleSpecies *species)
 Sets the species of this BaseInteractable. More...
 
const Vec3DgetForce () const
 Returns the force on this BaseInteractable. More...
 
const Vec3DgetTorque () const
 Returns the torque on this BaseInteractable. More...
 
void setForce (const Vec3D &force)
 Sets the force on this BaseInteractable. More...
 
void setTorque (const Vec3D &torque)
 Sets the torque on this BaseInteractable. More...
 
void addForce (const Vec3D &addForce)
 Adds an amount to the force on this BaseInteractable. More...
 
void addTorque (const Vec3D &addTorque)
 Adds an amount to the torque on this BaseInteractable. More...
 
virtual void resetForceTorque (int numberOfOMPthreads)
 
void sumForceTorqueOMP ()
 
const Vec3DgetPosition () const
 Returns the position of this BaseInteractable. More...
 
const QuaterniongetOrientation () const
 Returns the orientation of this BaseInteractable. More...
 
virtual void setPosition (const Vec3D &position)
 Sets the position of this BaseInteractable. More...
 
void setOrientationViaNormal (Vec3D normal)
 Sets the orientation of this BaseInteractable by defining the vector that results from the rotation of the (1,0,0) vector. More...
 
void setOrientationViaEuler (Vec3D eulerAngle)
 Sets the orientation of this BaseInteractable by defining the euler angles. More...
 
virtual void setOrientation (const Quaternion &orientation)
 Sets the orientation of this BaseInteractable. More...
 
virtual void move (const Vec3D &move)
 Moves this BaseInteractable by adding an amount to the position. More...
 
virtual void rotate (const Vec3D &angularVelocityDt)
 Rotates this BaseInteractable. More...
 
const std::vector< BaseInteraction * > & getInteractions () const
 Returns a list of interactions which belong to this interactable. More...
 
void addInteraction (BaseInteraction *I)
 Adds an interaction to this BaseInteractable. More...
 
bool removeInteraction (BaseInteraction *I)
 Removes an interaction from this BaseInteractable. More...
 
void copyInteractionsForPeriodicParticles (const BaseInteractable &p)
 Copies interactions to this BaseInteractable whenever a periodic copy made. More...
 
void setVelocity (const Vec3D &velocity)
 set the velocity of the BaseInteractable. More...
 
void setAngularVelocity (const Vec3D &angularVelocity)
 set the angular velocity of the BaseInteractble. More...
 
void addVelocity (const Vec3D &velocity)
 adds an increment to the velocity. More...
 
void addAngularVelocity (const Vec3D &angularVelocity)
 add an increment to the angular velocity. More...
 
virtual const Vec3DgetVelocity () const
 Returns the velocity of this interactable. More...
 
virtual const Vec3DgetAngularVelocity () const
 Returns the angular velocity of this interactable. More...
 
void setPrescribedPosition (const std::function< Vec3D(double)> &prescribedPosition)
 Allows the position of an infinite mass interactable to be prescribed. More...
 
void applyPrescribedPosition (double time)
 Computes the position from the user defined prescribed position function. More...
 
void setPrescribedVelocity (const std::function< Vec3D(double)> &prescribedVelocity)
 Allows the velocity of an infinite mass interactable to be prescribed. More...
 
void applyPrescribedVelocity (double time)
 Computes the velocity from the user defined prescribed velocity function. More...
 
void setPrescribedOrientation (const std::function< Quaternion(double)> &prescribedOrientation)
 Allows the orientation of the infinite mass interactbale to be prescribed. More...
 
void applyPrescribedOrientation (double time)
 Computes the orientation from the user defined prescribed orientation function. More...
 
void setPrescribedAngularVelocity (const std::function< Vec3D(double)> &prescribedAngularVelocity)
 Allows the angular velocity of the infinite mass interactable to be prescribed. More...
 
void applyPrescribedAngularVelocity (double time)
 Computes the angular velocity from the user defined prescribed angular velocity. More...
 
virtual const Vec3D getVelocityAtContact (const Vec3D &contact) const
 Returns the velocity at the contact point, use by many force laws. More...
 
void integrateBeforeForceComputation (double time, double timeStep)
 This is part of integrate routine for objects with infinite mass. More...
 
void integrateAfterForceComputation (double time, double timeStep)
 This is part of the integration routine for objects with infinite mass. More...
 
virtual Mdouble getInvMass () const
 
virtual Mdouble getCurvature (const Vec3D &labFixedCoordinates) const
 
virtual bool isFaceContact (const Vec3D &normal) const
 
- Public Member Functions inherited from BaseObject
 BaseObject ()=default
 Default constructor. More...
 
 BaseObject (const BaseObject &p)=default
 Copy constructor, copies all the objects BaseObject contains. More...
 
virtual ~BaseObject ()=default
 virtual destructor More...
 
virtual void moveInHandler (unsigned int index)
 Except that it is virtual, it does the same thing as setIndex() does. More...
 
void setIndex (unsigned int index)
 Allows one to assign an index to an object in the handler/container. More...
 
void setId (unsigned long id)
 Assigns a unique identifier to each object in the handler (container) which remains constant even after the object is deleted from the container/handler. More...
 
unsigned int getIndex () const
 Returns the index of the object in the handler. More...
 
unsigned int getId () const
 Returns the unique identifier of any particular object. More...
 
void setGroupId (unsigned groupId)
 
unsigned getGroupId () const
 

Private Attributes

Vec3D start_
 The centre of the lower end of the Coil. More...
 
Mdouble l_
 The length of the Coil. More...
 
Mdouble r_
 
Mdouble n_
 
Mdouble omega_
 
Mdouble offset_
 
Mdouble thickness_
 

Additional Inherited Members

- Static Public Member Functions inherited from BaseWall
static void addToVTK (const std::vector< Vec3D > &points, VTKContainer &vtk)
 Takes the points provided and adds a triangle strip connecting these points to the vtk container. More...
 

Detailed Description

This class defines a coil in the z-direction from a (constant) starting point, a (constant) length L, a (constant) radius r, a (constant) number or revelations N and a (constant) rotation speed (rev/s)

q is a new coordinate going from 0 to 1 and t is the time, x=xs+r*cos(2*pi(offset+N*q)), y=ys+r*sin(2*pi*(offset+N*q)) and z=zs+q*L

Todo:

IFCD: Can someone look at the details of the documentation of class Coil? I can't make sense of them.

Coil is now fixed in Z-direction, centered around the Z-axis. Consider converting to more general parameters, with direction of choise and central axis of choice.

Constructor & Destructor Documentation

◆ Coil() [1/3]

Coil::Coil ( )

Default constructor, sets a coil with default parameters.

Default constructor, make a Coil which is centered in the origin, has length 1, radius 1, 1 revelation, rotates with 1 rev/s, has thickness 0 and has not rotated.

15 {
16  start_.setZero();
17  l_ = 1.0;
18  r_ = 1.0;
19  n_ = 1.0;
20  omega_ = 1.0;
21  offset_ = 0.0;
22  thickness_ = 0.0;
23  logger(DEBUG, "Coil() constructor finished");
24 }
Logger< MERCURYDPM_LOGLEVEL > logger("MercuryKernel")
Definition of different loggers with certain modules. A user can define its own custom logger here.
@ DEBUG
Mdouble n_
Definition: Coil.h:105
Mdouble thickness_
Definition: Coil.h:120
Mdouble omega_
Definition: Coil.h:110
Mdouble l_
The length of the Coil.
Definition: Coil.h:95
Mdouble r_
Definition: Coil.h:100
Vec3D start_
The centre of the lower end of the Coil.
Definition: Coil.h:90
Mdouble offset_
Definition: Coil.h:115
void setZero()
Sets all elements to zero.
Definition: Vector.cc:23

References DEBUG, l_, logger, n_, offset_, omega_, r_, Vec3D::setZero(), start_, and thickness_.

Referenced by copy().

◆ Coil() [2/3]

Coil::Coil ( const Coil other)

Copy constructor, makes a coil with the same properties as the input Coil.

Parameters
[in]otherCoil to be copied into this coil
29  : BaseWall(other)
30 {
31  start_ = other.start_;
32  l_ = other.l_;
33  r_ = other.r_;
34  n_ = other.n_;
35  omega_ = other.omega_;
36  offset_ = other.offset_;
37  thickness_ = other.thickness_;
38  logger(DEBUG, "Coil(const Coil&) copy constructor finished");
39 }
BaseWall()
Default constructor.
Definition: BaseWall.cc:15

References DEBUG, l_, logger, n_, offset_, omega_, r_, start_, and thickness_.

◆ Coil() [3/3]

Coil::Coil ( Vec3D  start,
Mdouble  l,
Mdouble  r,
Mdouble  n,
Mdouble  omega,
Mdouble  thickness 
)

Constructor in which all parameters are set.

Parameters
[in]startA Vec3D which denotes the centre of the lower end of the Coil.
[in]lThe length of the Coil, must be positive.
[in]rThe radius of the Coil, must be positive.
[in]nThe number of revelations of the Coil, must be positive.
[in]omegaThe rotation speed of the Coil in rev/s.
[in]thicknessThe thickness of the "spiral" of the Coil, must be non-negative.

Make a Coil by assigning all input parameters to the data-members of this class, and setting the offset_ to 0.

52 {
53  start_ = start;
54  l_ = l;
55  r_ = r;
56  n_ = n;
57  omega_ = omega;
58  thickness_ = thickness;
59  offset_ = 0.0;
60  logger(DEBUG, "Coil(params) constructor with parameters finished.");
61 }
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Vector::Scalar omega(const Vector &t, const Vector &s, RealScalar angle)
Definition: IDRS.h:36
r
Definition: UniformPSDSelfTest.py:20
void start(const unsigned &i)
(Re-)start i-th timer
Definition: oomph_utilities.cc:243

References DEBUG, l_, logger, n, n_, offset_, Eigen::internal::omega(), omega_, UniformPSDSelfTest::r, r_, oomph::CumulativeTimings::start(), start_, and thickness_.

◆ ~Coil()

Coil::~Coil ( )
override

Default destructor.

64 {
65  logger(DEBUG, "~Coil() destructor finished.");
66 }

References DEBUG, and logger.

Member Function Documentation

◆ copy()

Coil * Coil::copy ( ) const
overridevirtual

Copy this Coil and return a pointer to the copy, useful for polymorphism.

Returns
A pointer to a copy of this Coil.

Implements BaseWall.

96 {
97  return new Coil(*this);
98 }
Coil()
Default constructor, sets a coil with default parameters.
Definition: Coil.cc:14

References Coil().

◆ getDistanceAndNormal()

bool Coil::getDistanceAndNormal ( const BaseParticle p,
Mdouble distance,
Vec3D normal_return 
) const
overridevirtual

Compute the distance from the Coil for a given BaseParticle and return if there is a collision. If there is a collision, also return the normal vector of the interaction point.

Parameters
[in]pBaseParticle we want to calculate the distance and whether it collided of.
[out]distanceThe distance of the BaseParticle to this Coil.
[out]normal_returnIf there was a collision, the normal vector to this Coil will be placed here.
Returns
A boolean which says whether or not there was a collision.

This function computes whether or not there is a collision between a given BaseParticle and this Coil. If there is a collision, this function also computes the distance between the BaseParticle and Coil and the normal of the IntersectionOfWalls at the intersection point.

Todo:
Make this function readable and explain the steps in the details.

Implements BaseWall.

112 {
113  const Mdouble interactionRadius = p.getWallInteractionRadius(this) + thickness_;
114  Mdouble Rsqr = pow(p.getPosition().X - start_.X, 2) + pow(p.getPosition().Y - start_.Y, 2);
115  if (Rsqr > pow(r_ + interactionRadius, 2) ||
116  Rsqr < pow(r_ - interactionRadius, 2) ||
117  p.getPosition().Z > l_ + start_.Z + interactionRadius ||
118  p.getPosition().Z < start_.Z - interactionRadius)
119  {
120  //std::cout<<"Particle is out of first bound checking"<<std::endl;
121  //std::cout<<"Rule 1: "<<pow(r-P.getRadius()-thickness_,2)<<"<"<<Rsqr<<"<"<<pow(r+P.getRadius()+thickness_,2)<<std::endl;
122  //std::cout<<"Rule 2: "<<Start.Z-P.getRadius()-thickness_<<"<"<<P.getPosition().Z<<"<"<<L+Start.Z+P.getRadius()+thickness_<<std::endl;
123  return false;
124  }
125  Mdouble R = sqrt(Rsqr);
126  Mdouble alpha = atan2(p.getPosition().Y - start_.Y, p.getPosition().X - start_.X);
127  Mdouble dz = p.getPosition().Z - start_.Z;
128 
129  //To find the contact point we have to minimize (with respect to q)
130  //Distance^2=(x-x0-r*cos(2*Pi*(offset+N*q)))^2+(y-y0-r*sin(2*Pi*(offset+N*q)))^2+(z-z0-q*L)^2
131  //Using polar coordinated (i.e. x-x0=R*cos(alpha), y-y0=R*sin(alpha) and dz=z-z0)
132  //Distance^2=R^2+r^2-2*R*r*cos(alpha-2*Pi*(offset+N*q))+(dz-q*L)^2
133 
134  //So we have to minimize:
135  //Distance^2=R^2+r^2-2*R*r*cos(alpha-2*Pi*(offset+N*q))+(dz-q*L)^2
136  //For this we use the Euler algoritm
137 
138  Mdouble q; //Current guess
139  Mdouble dd; //Derivative at current guess
140  Mdouble ddd; //Second derivative at current guess
141  Mdouble q0 = dz / l_; //Minimum of the parabolic part
142 
143  //The initial guess will be in the maximum of the cos closest to the minimum of the parabolic part
144  //Minima of the cos are at
145  //alpha-2*Pi*(offset+N*q)=2*k*Pi (k=integer)
146  //q=alpha/(2*Pi*N)-k/N-offset/N (k=integer)
147 
148  Mdouble k = round(alpha / 2.0 / constants::pi - (offset_ + n_ * q0));
149  q = alpha / (2 * constants::pi * n_) - k / n_ - offset_ / n_;
150 
151  //Now apply Newton's method
152  do
153  {
154  dd = -4.0 * R * r_ * constants::pi * n_ * sin(alpha - 2.0 * constants::pi * (n_ * q + offset_)) -
155  2.0 * l_ * (dz - q * l_);
156  ddd = 8.0 * R * r_ * constants::sqr_pi * n_ * n_ * cos(alpha - 2.0 * constants::pi * (n_ * q + offset_)) +
157  2.0 * l_ * l_;
158  q -= dd / ddd;
159  } while (fabs(dd / ddd) > 1e-14);
160 
161  //Check if the location is actually on the coil, otherwise a point collision with the end of the coil calculated
162  if (q < 0) //Left boundary of the coil
163  {
164  q = 0;
165  }
166  else if (q > 1) //right boundary of the coil
167  {
168  q = 1;
169  }
170 
171  Mdouble distanceSquared =
172  R * R + r_ * r_ - 2 * R * r_ * cos(alpha - 2 * constants::pi * (offset_ + n_ * q)) + pow(dz - q * l_, 2);
173  //If distance is too large there is no contact
174  if (distanceSquared >= mathsFunc::square(p.getWallInteractionRadius(this) + thickness_))
175  {
176  //std::cout<<"Particle is out of second bound checking, distance^2="<<distance<<" max="<<(P.getRadius()+thickness_)*(P.getRadius()+thickness_)<<std::endl;
177  return false;
178  }
179 
180  Vec3D ContactPoint;
181  distance = sqrt(distanceSquared) - thickness_;
182  ContactPoint.X = start_.X + r_ * cos(2.0 * constants::pi * (offset_ + n_ * q));
183  ContactPoint.Y = start_.Y + r_ * sin(2.0 * constants::pi * (offset_ + n_ * q));
184  ContactPoint.Z = start_.Z + q * l_;
185  normal_return = ContactPoint - p.getPosition();
186  normal_return /= normal_return.getLength();
187  return true;
188 }
AnnoyingScalar cos(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:136
AnnoyingScalar atan2(const AnnoyingScalar &y, const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:139
AnnoyingScalar sin(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:137
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
Array< double, 1, 3 > e(1./3., 0.5, 2.)
@ R
Definition: StatisticsVector.h:21
float * p
Definition: Tutorial_Map_using.cpp:9
Definition: Kernel/Math/Vector.h:30
Mdouble Y
Definition: Kernel/Math/Vector.h:45
Mdouble Z
Definition: Kernel/Math/Vector.h:45
Mdouble X
the vector components
Definition: Kernel/Math/Vector.h:45
static Mdouble getLength(const Vec3D &a)
Calculates the length of a Vec3D: .
Definition: Vector.cc:350
RealScalar alpha
Definition: level1_cplx_impl.h:151
char char char int int * k
Definition: level2_impl.h:374
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16 pow(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:625
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16 round(const bfloat16 &a)
Definition: BFloat16.h:646
EIGEN_DEVICE_FUNC const Scalar & q
Definition: SpecialFunctionsImpl.h:2019
Real fabs(const Real &a)
Definition: boostmultiprec.cpp:117
const Mdouble sqr_pi
Definition: ExtendedMath.h:25
const Mdouble pi
Definition: ExtendedMath.h:23
T square(const T val)
squares a number
Definition: ExtendedMath.h:86

References alpha, atan2(), cos(), e(), boost::multiprecision::fabs(), Vec3D::getLength(), k, l_, n_, offset_, p, constants::pi, Eigen::bfloat16_impl::pow(), Eigen::numext::q, R, r_, Eigen::bfloat16_impl::round(), sin(), constants::sqr_pi, sqrt(), mathsFunc::square(), start_, thickness_, Vec3D::X, Vec3D::Y, and Vec3D::Z.

◆ getName()

std::string Coil::getName ( ) const
overridevirtual

Returns the name of the object, in this case the string "Coil".

Returns
The string "Coil".

Implements BaseObject.

246 {
247  return "Coil";
248 }

◆ move_time()

void Coil::move_time ( Mdouble  dt)

Rotate the Coil for a period dt, so that the offset_ changes with omega_*dt.

Parameters
[in]dtThe time for which the Screw has to be turned.
194 {
195  offset_ += omega_ * dt;
196 }

References offset_, and omega_.

Referenced by CoilSelfTest::actionsBeforeTimeStep().

◆ oldRead()

void Coil::oldRead ( std::istream &  is)

Reads an old-style Coil from an input stream, for example an old restart file.

Deprecated:
If you have old restart files, please convert them to the current version of the restart files and use read(std::istream&) instead of oldRead(std::istream&).
Parameters
[in,out]isThe input stream from which the Coil is read.
218 {
219  std::string dummy;
220  is >> dummy >> start_
221  >> dummy >> l_
222  >> dummy >> r_ >> dummy >> n_
223  >> dummy >> omega_
224  >> dummy >> offset_;
225 }
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286

References l_, n_, offset_, omega_, r_, start_, and oomph::Global_string_for_annotation::string().

◆ read()

void Coil::read ( std::istream &  is)
overridevirtual

Reads a Coil from an input stream, for example a restart file.

Parameters
[in,out]isThe input stream from which the Coil is read.

Reimplemented from BaseWall.

202 {
203  BaseWall::read(is);
204  std::string dummy;
205  is >> dummy >> start_
206  >> dummy >> l_
207  >> dummy >> r_
208  >> dummy >> n_
209  >> dummy >> omega_
210  >> dummy >> thickness_
211  >> dummy >> offset_;
212 }
void read(std::istream &is) override
Function that reads a BaseWall from an input stream, usually a restart file.
Definition: BaseWall.cc:57

References l_, n_, offset_, omega_, r_, BaseWall::read(), start_, oomph::Global_string_for_annotation::string(), and thickness_.

◆ set()

void Coil::set ( Vec3D  start,
Mdouble  length,
Mdouble  radius,
Mdouble  numberOfRevelations,
Mdouble  omega,
Mdouble  thickness 
)

Set all parameters of this Coil.

Set the parameters of this Coil by assigning all input parameters to the data-members of this class, and setting the offset_ to 0.

Parameters
[in]startA Vec3D which denotes the centre of the lower end of the Coil.
[in]lengthThe length of the Coil, must be positive.
[in]radiusThe radius of the Coil, must be positive.
[in]numberOfRevelationsThe number of revelations of the Coil, must be positive.
[in]omegaThe rotation speed of the Coil in rev/s.
[in]thicknessThe thickness of the "spiral" of the Coil, must be non-negative.
82 {
83  start_ = start;
84  l_ = length;
85  r_ = radius;
86  n_ = numberOfRevelations;
87  omega_ = omega;
88  thickness_ = thickness;
89  offset_ = 0.0;
90 }
radius
Definition: UniformPSDSelfTest.py:15

References l_, n_, offset_, Eigen::internal::omega(), omega_, r_, UniformPSDSelfTest::radius, oomph::CumulativeTimings::start(), start_, and thickness_.

Referenced by CoilSelfTest::setupInitialConditions().

◆ write()

void Coil::write ( std::ostream &  os) const
overridevirtual

Writes a Coil to an output stream, for example a restart file.

Parameters
[in,out]osThe outpus stream to which the Coil is written.

Reimplemented from BaseWall.

231 {
232  BaseWall::write(os);
233  os << " Start " << start_
234  << " Length " << l_
235  << " Radius " << r_
236  << " Revolutions " << n_
237  << " Omega " << omega_
238  << " Thickness " << thickness_
239  << " Offset " << offset_;
240 }
void write(std::ostream &os) const override
Function that writes a BaseWall to an output stream, usually a restart file.
Definition: BaseWall.cc:79

References l_, n_, offset_, omega_, r_, start_, thickness_, and BaseWall::write().

Member Data Documentation

◆ l_

Mdouble Coil::l_
private

The length of the Coil.

Referenced by Coil(), getDistanceAndNormal(), oldRead(), read(), set(), and write().

◆ n_

Mdouble Coil::n_
private

The number of revelations of the Coil.

Referenced by Coil(), getDistanceAndNormal(), oldRead(), read(), set(), and write().

◆ offset_

Mdouble Coil::offset_
private

The amount the Coil has rotated, this is a number between 0 and 1.

Referenced by Coil(), getDistanceAndNormal(), move_time(), oldRead(), read(), set(), and write().

◆ omega_

Mdouble Coil::omega_
private

The rotation speed of the Coil, in rev/s.

Referenced by Coil(), move_time(), oldRead(), read(), set(), and write().

◆ r_

Mdouble Coil::r_
private

The radius of the Coil.

Referenced by Coil(), getDistanceAndNormal(), oldRead(), read(), set(), and write().

◆ start_

Vec3D Coil::start_
private

The centre of the lower end of the Coil.

Referenced by Coil(), getDistanceAndNormal(), oldRead(), read(), set(), and write().

◆ thickness_

Mdouble Coil::thickness_
private

The thickness of the "spiral" of the Coil.

Referenced by Coil(), getDistanceAndNormal(), read(), set(), and write().


The documentation for this class was generated from the following files: