NurbsWall Class Reference

This function defines a wall via a NurbsSurface. More...

#include <NurbsWall.h>

+ Inheritance diagram for NurbsWall:

Public Member Functions

 NurbsWall ()
 Default constructor: make a wall with default parameters. More...
 
 NurbsWall (const NurbsWall &other)
 Copy constructor, copies another wall. More...
 
 NurbsWall (const NurbsSurface &nurbsSurface)
 Constructor in which all parameters of the wall are set. More...
 
 ~NurbsWall ()
 Default destructor. More...
 
NurbsWallcopy () const override
 Copy this wall and return a pointer to the copy. More...
 
void set (const NurbsSurface &nurbsSurface)
 Defines a wall, given a NurbsSurface. More...
 
bool getDistanceAndNormal (const BaseParticle &P, Mdouble &distance, Vec3D &normal_return) const final
 Compute the distance from the Screw 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 read (std::istream &is) override
 Reads this wall from an input stream, for example a restart file. More...
 
void write (std::ostream &os) const override
 Writes this wall to an output stream, for example a restart file. More...
 
std::string getName () const override
 Returns the name of the object, here the string "Screw". More...
 
void writeVTK (VTKContainer &vtk) const override
 
void writeWallDetailsVTK (VTKData &data) 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)
 
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 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
 

Protected Attributes

NurbsSurface nurbsSurface_
 

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 function defines a wall via a NurbsSurface.

Constructor & Destructor Documentation

◆ NurbsWall() [1/3]

NurbsWall::NurbsWall ( )

Default constructor: make a wall with default parameters.

Make a default NurbsWall which is centered in the origin, has a length of 1, one revelation, a radius of 1, turns with 1 revelation per second, is infinitely thin and starts at its normal initial point.

40 {
41  logger(DEBUG, "NurbsWall() constructor finished.");
42 }
Logger< MERCURYDPM_LOGLEVEL > logger("MercuryKernel")
Definition of different loggers with certain modules. A user can define its own custom logger here.
@ DEBUG

References DEBUG, and logger.

Referenced by copy().

◆ NurbsWall() [2/3]

NurbsWall::NurbsWall ( const NurbsWall other)

Copy constructor, copies another wall.

Parameters
[in]otherThe NurbsWall that has to be copied.
48  : BaseWall(other), nurbsSurface_(other.nurbsSurface_)
49 {
50  logger(DEBUG, "NurbsWall(const NurbsWall&) copy constructor finished.");
51 }
BaseWall()
Default constructor.
Definition: BaseWall.cc:36
NurbsSurface nurbsSurface_
Definition: NurbsWall.h:95

References DEBUG, and logger.

◆ NurbsWall() [3/3]

NurbsWall::NurbsWall ( const NurbsSurface nurbsSurface)

Constructor in which all parameters of the wall are set.

54  : nurbsSurface_(nurbsSurface)
55 {
56  logger(DEBUG, "NurbsWall(const NurbsSurface&) constructor finished.");
57 }

References DEBUG, and logger.

◆ ~NurbsWall()

NurbsWall::~NurbsWall ( )

Default destructor.

60 {
61  logger(DEBUG, "~NurbsWall() finished, destroyed the wall.");
62 }

References DEBUG, and logger.

Member Function Documentation

◆ copy()

NurbsWall * NurbsWall::copy ( ) const
overridevirtual

Copy this wall and return a pointer to the copy.

Implements BaseWall.

Reimplemented in WearableNurbsWall.

65 {
66  return new NurbsWall(*this);
67 }
NurbsWall()
Default constructor: make a wall with default parameters.
Definition: NurbsWall.cc:39

References NurbsWall().

◆ getDistanceAndNormal()

bool NurbsWall::getDistanceAndNormal ( const BaseParticle P,
Mdouble distance,
Vec3D normal_return 
) const
finalvirtual

Compute the distance from the Screw 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.

Implements BaseWall.

75 {
76  //transform coordinates into position-orientation frame
77  Vec3D position = p.getPosition() - getPosition();
78  getOrientation().rotateBack(position);
80  if (nurbsSurface_.getDistance(position, p.getRadius()+s->getInteractionDistance(), distance, normal_return)) {
81  getOrientation().rotate(normal_return);
82  return true;
83  } else {
84  return false;
85  }
86 }
DPMBase * getDPMBase()
Gets the problem that is solved using this handler.
Definition: BaseHandler.h:725
const Quaternion & getOrientation() const
Returns the orientation of this BaseInteractable.
Definition: BaseInteractable.h:230
const ParticleSpecies * getSpecies() const
Returns a pointer to the species of this BaseInteractable.
Definition: BaseInteractable.h:108
const Vec3D & getPosition() const
Returns the position of this BaseInteractable.
Definition: BaseInteractable.h:218
BaseSpecies is the class from which all other species are derived.
Definition: BaseSpecies.h:50
Mdouble getInteractionDistance() const
returns the largest separation distance at which adhesive short-range forces can occur.
Definition: BaseSpecies.h:146
WallHandler * getHandler() const
A function which returns the WallHandler that handles this BaseWall.
Definition: BaseWall.cc:134
SpeciesHandler speciesHandler
A handler to that stores the species type i.e. LinearViscoelasticSpecies, etc.
Definition: DPMBase.h:1427
bool getDistance(Vec3D P, double radius, double &distance, Vec3D &normal) const
Definition: NurbsSurface.cc:198
void rotate(Vec3D &position) const
Definition: Quaternion.cc:563
void rotateBack(Vec3D &position) const
Definition: Quaternion.cc:610
std::enable_if<!std::is_pointer< typename U::MixedSpeciesType >::value, typename U::MixedSpeciesType * >::type getMixedObject(const U *S, const U *T)
Definition: SpeciesHandler.h:74
Definition: Vector.h:51

References NurbsSurface::getDistance(), BaseHandler< T >::getDPMBase(), BaseWall::getHandler(), BaseSpecies::getInteractionDistance(), SpeciesHandler::getMixedObject(), BaseInteractable::getOrientation(), BaseInteractable::getPosition(), BaseParticle::getRadius(), BaseInteractable::getSpecies(), nurbsSurface_, Quaternion::rotate(), Quaternion::rotateBack(), and DPMBase::speciesHandler.

◆ getName()

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

Returns the name of the object, here the string "Screw".

Returns
The string "NurbsWall".

Implements BaseObject.

Reimplemented in WearableNurbsWall.

112 {
113  return "NurbsWall";
114 }

◆ read()

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

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

Parameters
[in,out]isInput stream from which the wall must be read.

Reimplemented from BaseWall.

Reimplemented in WearableNurbsWall.

92 {
93  BaseWall::read(is);
94  std::string dummy;
95  is >> dummy;
96  is >> nurbsSurface_;
97 }
void read(std::istream &is) override
Function that reads a BaseWall from an input stream, usually a restart file.
Definition: BaseWall.cc:78

References nurbsSurface_, and BaseWall::read().

Referenced by WearableNurbsWall::read().

◆ set()

void NurbsWall::set ( const NurbsSurface nurbsSurface)

Defines a wall, given a NurbsSurface.

69  {
70  nurbsSurface_ = nurbsSurface;
71 }

References nurbsSurface_.

Referenced by main(), and Nurbs::setupInitialConditions().

◆ write()

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

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

Parameters
[in,out]osOutput stream to which the wall must be written.

Reimplemented from BaseWall.

Reimplemented in WearableNurbsWall.

103 {
104  BaseWall::write(os);
105  os << " NurbsSurface " << nurbsSurface_;
106 }
void write(std::ostream &os) const override
Function that writes a BaseWall to an output stream, usually a restart file.
Definition: BaseWall.cc:100

References nurbsSurface_, and BaseWall::write().

Referenced by WearableNurbsWall::write().

◆ writeVTK()

void NurbsWall::writeVTK ( VTKContainer vtk) const
overridevirtual

adds extra information to the points and triangleStrips vectors needed to plot the wall in vtk format

Parameters
pointsCoordinates of the vertices of the triangulated surfaces (in the VTK file this is called POINTS)
triangleStripsIndices of three vertices forming one triangulated surface (in the VTK file this is called CELL)

Reimplemented from BaseWall.

117 {
118  unsigned nu = 5 * nurbsSurface_.getControlPoints().size();
119  unsigned nv = 5 * nurbsSurface_.getControlPoints()[0].size();
120 
121  const double lbU = nurbsSurface_.getLowerBoundU();
122  const double ubU = nurbsSurface_.getUpperBoundU();
123  const double lbV = nurbsSurface_.getLowerBoundV();
124  const double ubV = nurbsSurface_.getUpperBoundV();
125 
126  //create points
127  size_t nPoints = vtk.points.size();
128  for (double u = 0; u <= nu; u++) {
129  for (double v = 0; v <= nv; v++) {
130  Vec3D p = nurbsSurface_.evaluate(lbU+(ubU-lbU)*u/nu, lbV+(ubV-lbV)*v/nv);
131  getOrientation().rotate(p);
132  p += getPosition();
133  vtk.points.push_back(p);
134  }
135  }
136 
137  //create connectivity matrix
138  //vtk.triangleStrips.reserve(vtk.triangleStrips.size()+nu);
139  for (unsigned i = 0; i < nu; ++i) {
140  std::vector<double> cell;
141  cell.reserve(2*nv);
142  for (unsigned j = 0; j <= nv; ++j) {
143  cell.push_back(nPoints + j + i * (nv+1));
144  cell.push_back(nPoints + j + (i+1) * (nv+1));
145  }
146  vtk.triangleStrips.push_back(cell);
147  }
148 }
Vec3D evaluate(double u, double v) const
Definition: NurbsSurface.cc:170
const std::vector< std::vector< Vec3D > > & getControlPoints() const
Definition: NurbsSurface.h:141
double getLowerBoundV() const
Definition: NurbsSurface.h:131
double getUpperBoundU() const
Definition: NurbsSurface.h:125
double getLowerBoundU() const
Definition: NurbsSurface.h:119
double getUpperBoundV() const
Definition: NurbsSurface.h:137
const std::complex< Mdouble > i
Definition: ExtendedMath.h:51
std::vector< std::vector< double > > triangleStrips
Definition: BaseWall.h:40
std::vector< Vec3D > points
Definition: BaseWall.h:39

References NurbsSurface::evaluate(), NurbsSurface::getControlPoints(), NurbsSurface::getLowerBoundU(), NurbsSurface::getLowerBoundV(), BaseInteractable::getOrientation(), BaseInteractable::getPosition(), NurbsSurface::getUpperBoundU(), NurbsSurface::getUpperBoundV(), constants::i, nurbsSurface_, VTKContainer::points, Quaternion::rotate(), and VTKContainer::triangleStrips.

◆ writeWallDetailsVTK()

void NurbsWall::writeWallDetailsVTK ( VTKData data) const
overridevirtual

Reimplemented from BaseWall.

Reimplemented in WearableNurbsWall.

151 {
152  const std::vector<std::vector<Vec3D>>& controlPoints = nurbsSurface_.getControlPoints();
153  const std::vector<std::vector<double>>& weights = nurbsSurface_.getWeights();
154 
155  // Reserve memory for number of points and cells about to be added
156  const unsigned int np = controlPoints.size() * controlPoints[0].size();
157  const unsigned int nc = (controlPoints.size() - 1) * (controlPoints[0].size() - 1);
158  data.reservePoints(np, { "Weight", "ID" });
159  data.reserveCells(nc);
160 
161  // Number of points in v-direction
162  size_t nv = controlPoints[0].size();
163  // Point index offset, the point indices have to be offset by the number of points already present at the start (from previously added data)
164  size_t pio = data.getPoints().size();
165  // Get last id only when possible and increment it, otherwise start with 0
166  double id = data.getPointData()["ID"].empty() ? 0 : data.getPointData()["ID"].back() + 1;
167 
168  for (int i = 0; i < controlPoints.size(); i++)
169  {
170  for (int j = 0; j < controlPoints[i].size(); j++)
171  {
172  Vec3D p = controlPoints[i][j];
173  getOrientation().rotate(p);
174  p += getPosition();
175  data.addToPoints(p);
176  data.addToPointData("Weight", weights[i][j]);
177  data.addToPointData("ID", id);
178 
179  if (i > 0 && j > 0)
180  {
181  // Basic 2D/1D mapping for indexing: nv * i + j (+ point index offset)
182  // 4 points to form a rectangle: point, point to the left, point down, point to the left and down
183  data.addToConnectivity({ nv*i+j+pio, nv*(i-1)+j+pio, nv*i+j-1+pio, nv*(i-1)+j-1+pio });
184  data.addToTypes(8);
185  }
186  }
187  }
188 }
const std::vector< std::vector< Mdouble > > & getWeights() const
Definition: NurbsSurface.h:144
void addToTypes(int type)
Adds a type to the types vector.
Definition: VTKData.cc:43
void addToPoints(Vec3D point)
Adds a point to the points vector.
Definition: VTKData.cc:28
std::unordered_map< std::string, std::vector< Mdouble > > getPointData() const
Definition: VTKData.cc:53
void reserveCells(unsigned int n)
Reserves additional memory for connectivity and types vectors.
Definition: VTKData.cc:76
void addToPointData(const std::string &key, Mdouble value)
Adds a value to the pointData values vector corresponding to the given key.
Definition: VTKData.cc:33
void reservePoints(unsigned int n, const std::vector< std::string > &keys={})
Reserves additional memory for the points vector and optionally for the values vectors of the pointDa...
Definition: VTKData.cc:68
std::vector< Vec3D > getPoints() const
Definition: VTKData.cc:48
void addToConnectivity(const std::vector< size_t > &indices)
Adds a vector of indices to the connectivity vector.
Definition: VTKData.cc:38

References VTKData::addToConnectivity(), VTKData::addToPointData(), VTKData::addToPoints(), VTKData::addToTypes(), NurbsSurface::getControlPoints(), BaseInteractable::getOrientation(), VTKData::getPointData(), VTKData::getPoints(), BaseInteractable::getPosition(), NurbsSurface::getWeights(), constants::i, nurbsSurface_, VTKData::reserveCells(), VTKData::reservePoints(), and Quaternion::rotate().

Member Data Documentation

◆ nurbsSurface_


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