SineWall Class Reference

#include <SineWall.h>

+ Inheritance diagram for SineWall:

Public Member Functions

 SineWall ()
 Default constructor, sets a chute with default parameters. More...
 
 SineWall (const SineWall &other)
 Copy constructor. More...
 
 SineWall (Mdouble length, Mdouble sw_wavn, Mdouble sw_phshift, Mdouble sw_amp)
 Constructor in which all parameters are set. More...
 
 ~SineWall () override
 
void set (Mdouble length, Mdouble sw_wavn, Mdouble sw_phshift, Mdouble sw_amp)
 
SineWallcopy () const override
 Pure virtual function that can be overwritten in inherited classes in order to copy a BaseWall. More...
 
bool getDistanceAndNormal (const BaseParticle &P, Mdouble &distance, Vec3D &normal_return) const override
 Pure virtual function that computes the distance of a BaseParticle to this wall and returns the normal of this wall if there is a collision. More...
 
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...
 
void read (std::istream &is) override
 
void write (std::ostream &os) const override
 
std::string getName () const override
 
- 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
 
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

Mdouble l_
 
Mdouble sw_amp_
 
Mdouble sw_wavn_
 
Mdouble sw_phshift_
 

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

A half-infinite surface which has a sinusoidal direction in the x direction, and no variation in the y direction.

The equation of the surface is z = sw_amp * sin( sw_wavn * x + sw_phshift ) for x < l Because this is nonlinear, we must solve for the interaction location using Newton's method.

Constructor & Destructor Documentation

◆ SineWall() [1/3]

SineWall::SineWall ( )

Default constructor, sets a chute with default parameters.

10 {
11  l_ = 1.0;
12  sw_wavn_ = 0.0;
13  sw_phshift_ = 0.0;
14  sw_amp_ = 0.0;
15  logger(DEBUG, "SineWall() constructor finished");
16 }
Logger< MERCURYDPM_LOGLEVEL > logger("MercuryKernel")
Definition of different loggers with certain modules. A user can define its own custom logger here.
@ DEBUG
Mdouble l_
Definition: SineWall.h:55
Mdouble sw_wavn_
Definition: SineWall.h:57
Mdouble sw_amp_
Definition: SineWall.h:56
Mdouble sw_phshift_
Definition: SineWall.h:58

References DEBUG, l_, logger, sw_amp_, sw_phshift_, and sw_wavn_.

Referenced by copy().

◆ SineWall() [2/3]

SineWall::SineWall ( const SineWall other)

Copy constructor.

18  : BaseWall(other)
19 {
20  l_ = other.l_;
21  sw_wavn_ = other.sw_wavn_;
22  sw_phshift_ = other.sw_phshift_;
23  sw_amp_ = other.sw_amp_;
24  logger(DEBUG, "SineWall copy constructor finished");
25 }
BaseWall()
Default constructor.
Definition: BaseWall.cc:15

References DEBUG, l_, logger, sw_amp_, sw_phshift_, and sw_wavn_.

◆ SineWall() [3/3]

SineWall::SineWall ( Mdouble  length,
Mdouble  sw_wavn,
Mdouble  sw_phshift,
Mdouble  sw_amp 
)

Constructor in which all parameters are set.

28 {
29  l_ = length;
30  sw_wavn_ = sw_wavn;
31  sw_phshift_ = sw_phshift;
32  sw_amp_ = sw_amp;
33  logger(DEBUG, "SineWall(params) constructor finished");
34 }

References DEBUG, l_, logger, sw_amp_, sw_phshift_, and sw_wavn_.

◆ ~SineWall()

SineWall::~SineWall ( )
override
37 {
38  logger(DEBUG, "SineWall destructor finished");
39 }

References DEBUG, and logger.

Member Function Documentation

◆ copy()

SineWall * SineWall::copy ( ) const
overridevirtual

Pure virtual function that can be overwritten in inherited classes in order to copy a BaseWall.

Returns
A pointer to the new BaseWall.

Implements BaseWall.

50 {
51  return new SineWall(*this);
52 }
SineWall()
Default constructor, sets a chute with default parameters.
Definition: SineWall.cc:9

References SineWall().

◆ getDistanceAndNormal()

bool SineWall::getDistanceAndNormal ( const BaseParticle P,
Mdouble distance,
Vec3D normal_return 
) const
overridevirtual

Pure virtual function that computes the distance of a BaseParticle to this wall and returns the normal of this wall if there is a collision.

Beware, the distance and normal are output parameters, not return values!

Parameters
[in]PReference to the BaseParticle we want to compute the distance to the BaseWall of.
[out]distanceDistance of the BaseParticle to the BaseWall.
[out]normal_returnThe normal of the wall. Is only given if there is a collision.
Returns
A boolean which indicates if there is a collision between the BaseParticle and the wall.

Implements BaseWall.

55 {
56  /* Define some shortcuts */
57  double x0 = p.getPosition().X;
58  double y0 = p.getPosition().Y;
59  double z0 = p.getPosition().Z;
60  double a = p.getRadius();
61 
62  // more shortcuts
63  double A = sw_amp_;
64  double k = sw_wavn_;
65  double ph = sw_phshift_;
66 
67  /* First, check whether the particle is definitely out of contact with the
68  * chute. This will be so if the particle's interaction radius is small and
69  * the particle's z-position is high. */
70  // TODO do this properly
71  if (z0 > A + a)
72  return false;
73 
74  /* Or, if the particle lies under the surface, then it is expelled
75  * perpendicularly. This shouldn't happen very often provided that the step
76  * size is small.*/
77  else if (z0 < A * sin(k * x0 + ph))
78  {
79  fprintf(stderr, "Particle at x0=%f, z0=%f is under SineWall\n", x0, z0);
80  // exit(-1);
81 
82  distance = A * sin(k * x0 + ph) - z0;
83  normal_return = Vec3D(0, 0, 1);
84  return true;
85  }
86 
87  /* Otherwise, find the contact point by using Newton's method to minimise
88  * (half of the squared) distance */
89  else
90  {
91  // Iterate according to Newton-Raphson
92  double q = x0;
93  double correction = 0;
94  // fprintf(stderr, "x0 = %lf, z0 = %lf, q = %lf\n", x0, z0, q);
95  do
96  {
97  double z = A * sin(k * q + ph);
98  double dzdq = A * k * cos(k * q + ph);
99  double d2zdq2 = -A * pow(k, 2) * sin(k * q + ph);
100 
101  // derivative of half the squared distance -- to be zeroed
102  double dRdq = q - x0 + (z - z0) * dzdq;
103  // second derivative of half the squared distance
104  double d2Rdq2 = 1 + (z - z0) * d2zdq2 + pow(dzdq, 2);
105  // The required correction to our current estimate. (Note the sign
106  // convention that this is to be subtracted.)
107  correction = dRdq / d2Rdq2;
108  //fprintf(stderr, "x0 = %.12lf, z0 = %.12lf, q = %.12lf, dRdq = %e, d2Rdq2 = %e, correction = %.e\n",
109  // x0, z0, q, dRdq, d2Rdq2, correction);
110  // N-R update
111  q -= correction;
112  // Perhaps this form is useful for reducing rounding errors?
113  // q = (q*d2Rdq2 - dRdq)/d2Rdq2;
114  } while (fabs(correction) > 1e-14);
115  //fprintf(stderr, "correction = %e\n", fabs(correction));
116 
117  // The standard Newton-Raphson iteration won't work, because the problem is
118  // too oscillatory (and N-R will converge on a near-root minimum, not the
119  // root). Neither does the following, a modified N-R which takes second and
120  // third derivatives into account.
121  /*
122  double q = x0;
123  double correction = 0;
124  do
125  {
126  // derivative of half the squared distance
127  double dRdq = q - x0 - A*pow(k,2)*z0*cos(k*q+ph) + 0.5*pow(A,2)*pow(k,3)*sin(2*k*q+2*ph);
128  // second derivative
129  double d2Rdq2 = 1 + pow(A,2)*pow(k,4)*cos(2*k*q+2*ph) + A*pow(k,3)*z0*sin(k*q+ph);
130  // third derivative
131  double d2Rdq2d = A*pow(k,4)*cos(k*q+ph) * (z0 - 4*A*k*sin(k*q+ph));
132 
133  double discriminant = d2Rdq2*d2Rdq2 - 2*dRdq*d2Rdq2d;
134  double newq = discriminant > 0 ?
135  q - (d2Rdq2 - sqrt(discriminant)) / d2Rdq2d
136  : q - dRdq / d2Rdq2;
137  double correction = newq/q - 1;
138  fprintf(stderr, "x0 = %.12lf, z0 = %.12lf, q = %.12lf, newq = %.12lf, correction = %e\n",
139  x0, z0, q, newq, correction);
140  q = newq;
141  } while (fabs(correction) > 1e-7);
142  fprintf(stderr, "converged, fabs(correction) = %e\n", fabs(correction));
143  */
144 
145  Mdouble distanceSquared
146  = pow(q - x0, 2) + pow(sw_amp_ * sin(sw_wavn_ * q + sw_phshift_) - z0, 2);
147 
148  if (distanceSquared >= mathsFunc::square(p.getWallInteractionRadius(this))
149  && z0 > sw_amp_ * sin(sw_wavn_ * q + sw_phshift_))
150  {
151  return false;
152  }
153  else
154  {
155  Vec3D ContactPoint;
156  distance = sqrt(distanceSquared);
157  ContactPoint.X = q;
158  ContactPoint.Y = y0;
159  ContactPoint.Z = sw_amp_ * sin(sw_wavn_ * q + sw_phshift_);
160  // normal_return = ContactPoint - p.getPosition();
161  normal_return = ContactPoint - p.getPosition();
162  normal_return /= normal_return.getLength();
163  return true;
164  }
165  }
166 }
AnnoyingScalar cos(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:136
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.)
float * p
Definition: Tutorial_Map_using.cpp:9
The matrix class, also used for vectors and row-vectors.
Definition: Eigen/Eigen/src/Core/Matrix.h:186
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
const Scalar * a
Definition: level2_cplx_impl.h:32
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_DEVICE_FUNC const Scalar & q
Definition: SpecialFunctionsImpl.h:2019
Vector< double > x0(2, 0.0)
Real fabs(const Real &a)
Definition: boostmultiprec.cpp:117
T square(const T val)
squares a number
Definition: ExtendedMath.h:86

References a, cos(), e(), boost::multiprecision::fabs(), Vec3D::getLength(), k, p, Eigen::bfloat16_impl::pow(), Eigen::numext::q, sin(), sqrt(), mathsFunc::square(), sw_amp_, sw_phshift_, sw_wavn_, Vec3D::X, Global::x0, Vec3D::Y, and Vec3D::Z.

Referenced by getInteractionWith().

◆ getInteractionWith()

BaseInteraction * SineWall::getInteractionWith ( BaseParticle p,
unsigned  timeStamp,
InteractionHandler interactionHandler 
)
overridevirtual

Returns the interaction between this wall and a given particle, nullptr if there is no interaction.

Parameters
[in]pPointer to the BaseParticle which we want to check the interaction for.
[in]timeStampThe time at which we want to look at the interaction.
[in]interactionHandlerA pointer to the InteractionHandler in which the interaction can be found.
Returns
A pointer to the BaseInteraction that happened between this BaseWall and the BaseParticle at the timeStamp.
Todo:
Hacked please fix

Reimplemented from BaseWall.

170 {
171  Mdouble distance;
172  Vec3D normal;
173  if (getDistanceAndNormal(*p, distance, normal))
174  {
175  BaseInteraction* c = interactionHandler->getInteraction(p, this, timeStamp);
176  c->setNormal(-normal);
177  c->setDistance(distance);
178  c->setOverlap(p->getRadius() - distance);
179  // c->setContactPoint( p->getPosition() - (p->getRadius() - 0.5*c->getOverlap()) * c->getNormal());
180  c->setContactPoint(p->getPosition() + distance * normal);
182  return c;
183  }
184  else
185  return nullptr;
186 }
Stores information about interactions between two interactable objects; often particles but could be ...
Definition: BaseInteraction.h:39
BaseInteraction * getInteraction(BaseInteractable *P, BaseInteractable *I, unsigned timeStamp)
Returns the Interaction between the BaseInteractable's P and I.
Definition: InteractionHandler.cc:126
bool getDistanceAndNormal(const BaseParticle &P, Mdouble &distance, Vec3D &normal_return) const override
Pure virtual function that computes the distance of a BaseParticle to this wall and returns the norma...
Definition: SineWall.cc:54
void normal(const Vector< double > &x, Vector< double > &normal)
Definition: free_surface_rotation.cc:65
int c
Definition: calibrate.py:100

References calibrate::c, getDistanceAndNormal(), InteractionHandler::getInteraction(), WallFunction::normal(), and p.

◆ getName()

std::string SineWall::getName ( ) const
overridevirtual
Returns
The string "SineWall".

Implements BaseObject.

217 {
218  return "SineWall";
219 }

◆ read()

void SineWall::read ( std::istream &  is)
overridevirtual
Parameters
[in,out]isThe input stream from which the SineWall is read.

Reimplemented from BaseWall.

192 {
193  BaseWall::read(is);
194  std::string dummy;
195  is >> dummy >> l_
196  >> dummy >> sw_wavn_
197  >> dummy >> sw_phshift_
198  >> dummy >> sw_amp_;
199 }
void read(std::istream &is) override
Function that reads a BaseWall from an input stream, usually a restart file.
Definition: BaseWall.cc:57
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286

References l_, BaseWall::read(), oomph::Global_string_for_annotation::string(), sw_amp_, sw_phshift_, and sw_wavn_.

◆ set()

void SineWall::set ( Mdouble  length,
Mdouble  sw_wavn,
Mdouble  sw_phshift,
Mdouble  sw_amp 
)
42 {
43  l_ = length;
44  sw_wavn_ = sw_wavn;
45  sw_phshift_ = sw_phshift;
46  sw_amp_ = sw_amp;
47 }

References l_, sw_amp_, sw_phshift_, and sw_wavn_.

◆ write()

void SineWall::write ( std::ostream &  os) const
overridevirtual
Parameters
[in,out]osThe outpus stream to which the SineWall is written.

Reimplemented from BaseWall.

205 {
206  BaseWall::write(os);
207  os << " Length " << l_
208  << " sw_wavn " << sw_wavn_
209  << " sw_phshift " << sw_phshift_
210  << " sw_amp " << sw_amp_;
211 }
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_, sw_amp_, sw_phshift_, sw_wavn_, and BaseWall::write().

Member Data Documentation

◆ l_

Mdouble SineWall::l_
private

Referenced by read(), set(), SineWall(), and write().

◆ sw_amp_

Mdouble SineWall::sw_amp_
private

◆ sw_phshift_

Mdouble SineWall::sw_phshift_
private

◆ sw_wavn_

Mdouble SineWall::sw_wavn_
private

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