SuperQuadricParticle.h
Go to the documentation of this file.
1 //Copyright (c) 2013-2023, The MercuryDPM Developers Team. All rights reserved.
2 //For the list of developers, see <http://www.MercuryDPM.org/Team>.
3 //
4 //Redistribution and use in source and binary forms, with or without
5 //modification, are permitted provided that the following conditions are met:
6 // * Redistributions of source code must retain the above copyright
7 // notice, this list of conditions and the following disclaimer.
8 // * Redistributions in binary form must reproduce the above copyright
9 // notice, this list of conditions and the following disclaimer in the
10 // documentation and/or other materials provided with the distribution.
11 // * Neither the name MercuryDPM nor the
12 // names of its contributors may be used to endorse or promote products
13 // derived from this software without specific prior written permission.
14 //
15 //THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND
16 //ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
17 //WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
18 //DISCLAIMED. IN NO EVENT SHALL THE MERCURYDPM DEVELOPERS TEAM BE LIABLE FOR ANY
19 //DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES
20 //(INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
21 //LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND
22 //ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
23 //(INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
24 //SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
25 
26 #ifndef SuperQuadricParticle_H
27 #define SuperQuadricParticle_H
28 
29 #include "Math/SmallVector.h"
30 #include "BaseParticle.h"
31 #include "BaseInteractable.h"
32 #include "NonSphericalParticle.h"
34 
37 
57 {
58 public:
64 
71 
77 
82  ~SuperQuadricParticle() override;
83 
87  SuperQuadricParticle* copy() const override;
88 
92  void write(std::ostream& os) const override;
93 
98  void read(std::istream& is) override;
99 
103  std::string getName() const override;
104 
111  void setAxesAndExponents(const Mdouble& a1, const Mdouble& a2, const Mdouble& a3, const Mdouble& eps1,
112  const Mdouble& eps2);
113 
120  void setAxesAndExponents(const Vec3D& axes, const Mdouble& eps1, const Mdouble& eps2);
121 
127  void setAxes(const Mdouble& a1, const Mdouble& a2, const Mdouble& a3);
128 
134  void setAxes(const Vec3D& axes) override;
135 
141  void setExponents(const Mdouble& eps1, const Mdouble& eps2) override;
142 
145 
150  Vec3D getAxes() const override;
151 
157  Mdouble getExponentEps1() const override;
158 
164  Mdouble getExponentEps2() const override;
165 
169  Mdouble getVolume() const override;
170 
174  void setInertia() override;
175 
176  void setRadius(const Mdouble radius) override;
177 
182  BaseInteraction* getInteractionWith(BaseParticle* P, unsigned timeStamp,
183  InteractionHandler* interactionHandler) override;
184 
189  Mdouble getCurvature(const LabFixedCoordinates& labFixedCoordinates) const override;
190 
194  bool isInContactWith(const BaseParticle* p) const override;
195 
199  SmallVector<3> computeShapeGradientLabFixed(const LabFixedCoordinates& labFixedCoordinates) const;
200 
201 
207  InteractionHandler* interactionHandler);
208 
212  SmallMatrix<3, 3> computeHessianLabFixed(const LabFixedCoordinates& labFixedCoordinates) const;
213 
217  Mdouble computeShape(const LabFixedCoordinates& labFixedCoordinates) const;
218 
224  const SuperQuadricParticle* p1,
225  const SuperQuadricParticle* p2) const;
226 
232  const SuperQuadricParticle* p1,
233  const SuperQuadricParticle* p2) const;
234 
239 
243  Mdouble overlapFromContactPoint(const LabFixedCoordinates& contactPoint, const LabFixedCoordinates& normal) const;
244 
249 
254  SmallVector<4> getContactPointPlanB(const SuperQuadricParticle* pOther, unsigned numberOfSteps) const;
255 
260  bool computeContactPoint(SmallVector<4>& contactPoint, const SuperQuadricParticle* p1,
261  const SuperQuadricParticle* p2) const;
262 
263 
264  void writeDebugMessageStep1(const SuperQuadricParticle* pQuad, const SmallVector<4>& contactPointPlanB) const;
265 
266  void writeDebugMessageStep2(const SuperQuadricParticle* pQuad, const Vec3D& dAxesThis, const Mdouble& dn11,
267  const Mdouble& dn12, const Vec3D& dAxesOther, const Mdouble& dn21,
268  const Mdouble& dn22) const;
269 
270  void writeDebugMessageStep3(const Vec3D& axesThis, const Mdouble& n11, const Mdouble& n12, const Vec3D& axesOther,
271  const Mdouble& n21, const Mdouble& n22) const;
272 
273  void
275  const unsigned int& counter) const;
276 
280  Mdouble getInteractionRadius(const BaseParticle* particle) const;
281 
283  void computeMass(const ParticleSpecies& s) override;
284 
285 private:
286 
290  void setBoundingRadius();
291 
298 
303 };
304 
305 #endif
Vec3D LabFixedCoordinates
Definition: SuperQuadricParticle.h:36
Vec3D BodyFixedCoordinates
Definition: SuperQuadricParticle.h:35
Stores information about interactions between two interactable objects; often particles but could be ...
Definition: BaseInteraction.h:60
Definition: BaseParticle.h:54
Container to store Interaction objects.
Definition: InteractionHandler.h:45
Base class for all non-spherical particle types.
Definition: NonSphericalParticle.h:37
Definition: ParticleSpecies.h:37
Data type for small dense matrix.
Definition: SmallMatrix.h:68
Definition: SmallVector.h:62
Definition: SuperQuadricParticle.h:57
void writeDebugMessageStep1(const SuperQuadricParticle *pQuad, const SmallVector< 4 > &contactPointPlanB) const
Definition: SuperQuadricParticle.cc:953
SmallVector< 3 > computeShapeGradientLabFixed(const LabFixedCoordinates &labFixedCoordinates) const
Compute and get the gradient of the shape-function at the given (lab-fixed) position.
Definition: SuperQuadricParticle.cc:493
void read(std::istream &is) override
Read function: read in the information for this superquadric from the given input-stream,...
Definition: SuperQuadricParticle.cc:122
void setRadius(const Mdouble radius) override
Sets the particle's radius_ (and adjusts the mass_ accordingly, based on the particle's species)
Definition: SuperQuadricParticle.cc:239
Vec3D axes_
Lengths of principal axes (a1, a2, a3).
Definition: SuperQuadricParticle.h:302
void setExponents(const Mdouble &eps1, const Mdouble &eps2) override
Set the exponents to eps1 and eps2 for this superquadric. We use the super-ellipsoid definition state...
Definition: SuperQuadricParticle.cc:169
Vec3D getAxes() const override
Get the axes-lengths of this superquadric. We use the super-ellipsoid definition stated in Chapter 2 ...
Definition: SuperQuadricParticle.cc:182
SmallVector< 4 > getContactPoint(const SuperQuadricParticle *p, BaseInteraction *C) const
Compute the contact point between this and the given superquadric particle.
Definition: SuperQuadricParticle.cc:371
Mdouble getExponentEps2() const override
Get the second exponent of this superquadric. We use the super-ellipsoid definition stated in Chapter...
Definition: SuperQuadricParticle.cc:192
BaseInteraction * getInteractionWithSuperQuad(SuperQuadricParticle *p, unsigned timeStamp, InteractionHandler *interactionHandler)
Checks if this superquadric is in interaction with the given superquadric, and if so,...
Definition: SuperQuadricParticle.cc:322
SuperQuadricParticle * copy() const override
Copy method. It calls to copy constructor of this superquadric, useful for polymorphism.
Definition: SuperQuadricParticle.cc:90
void writeDebugMessageStep3(const Vec3D &axesThis, const Mdouble &n11, const Mdouble &n12, const Vec3D &axesOther, const Mdouble &n21, const Mdouble &n22) const
Definition: SuperQuadricParticle.cc:891
void computeMass(const ParticleSpecies &s) override
Computes the particle's (inverse) mass and inertia.
Definition: SuperQuadricParticle.cc:983
SmallMatrix< 3, 3 > computeHessianLabFixed(const LabFixedCoordinates &labFixedCoordinates) const
Compute and get the hessian ("second derivative") of the shape-function at the given (lab-fixed) posi...
Definition: SuperQuadricParticle.cc:526
Mdouble getCurvature(const LabFixedCoordinates &labFixedCoordinates) const override
Get the mean curvature of this superquadric at the given (lab-fixed) position, see Podlozhyuk et al....
Definition: SuperQuadricParticle.cc:664
void writeDebugMessageMiddleOfLoop(const SuperQuadricParticle &p1, const SuperQuadricParticle &p2, SmallVector< 4 > &contactPointPlanB, const unsigned int &counter) const
Definition: SuperQuadricParticle.cc:869
std::string getName() const override
Returns the name of the class, here "SuperQuadric".
Definition: SuperQuadricParticle.cc:112
void writeDebugMessageStep2(const SuperQuadricParticle *pQuad, const Vec3D &dAxesThis, const Mdouble &dn11, const Mdouble &dn12, const Vec3D &dAxesOther, const Mdouble &dn21, const Mdouble &dn22) const
Definition: SuperQuadricParticle.cc:915
SmallVector< 4 > computeResidualContactDetection(const SmallVector< 4 > &position, const SuperQuadricParticle *p1, const SuperQuadricParticle *p2) const
Objective function for contact detection between the two given superquadrics. See Podlozhyuk et al....
Definition: SuperQuadricParticle.cc:566
Mdouble eps1_
Blockiness parameters.
Definition: SuperQuadricParticle.h:297
void setInertia() override
Compute and set the inertia-tensor for this superquadric. For internal use only.
Definition: SuperQuadricParticle.cc:216
SmallVector< 4 > getContactPointPlanB(const SuperQuadricParticle *pOther, unsigned numberOfSteps) const
If the "normal" procedure fails to find a contact point, use an alternative approach that involves st...
Definition: SuperQuadricParticle.cc:720
bool computeContactPoint(SmallVector< 4 > &contactPoint, const SuperQuadricParticle *p1, const SuperQuadricParticle *p2) const
Perform the actual Newton-iterations to find the contact point. Note, that it is given back as a para...
Definition: SuperQuadricParticle.cc:816
Mdouble getInteractionRadius(const BaseParticle *particle) const
returns the radius plus half the interactionDistance of the mixed species
Definition: SuperQuadricParticle.cc:976
void write(std::ostream &os) const override
Write function: write this SuperQuadric to the given output-stream, for example a restart-file.
Definition: SuperQuadricParticle.cc:100
Mdouble overlapFromContactPoint(const LabFixedCoordinates &contactPoint, const LabFixedCoordinates &normal) const
Compute the distance between the contact-point and surface of this superquadric particle.
Definition: SuperQuadricParticle.cc:392
bool isInContactWith(const BaseParticle *p) const override
Get whether or not this superquadric is in contact with the given particle.
Definition: SuperQuadricParticle.cc:685
SuperQuadricParticle()
Basic Particle constructor, creates a superquadric with axes (1,1,1) and exponents (2,...
Definition: SuperQuadricParticle.cc:38
~SuperQuadricParticle() override
Destructor, needs to be implemented and checked to see if it is the largest or smallest particle curr...
Definition: SuperQuadricParticle.cc:76
Mdouble eps2_
Definition: SuperQuadricParticle.h:297
BaseInteraction * getInteractionWith(BaseParticle *P, unsigned timeStamp, InteractionHandler *interactionHandler) override
Checks if this superquadric is in interaction with the given particle, and if so, returns vector of p...
Definition: SuperQuadricParticle.cc:286
SmallVector< 4 > getInitialGuessForContact(const SuperQuadricParticle *pQuad, BaseInteraction *C) const
Get an initial guess for the contact-point between this particle and the given particle.
Definition: SuperQuadricParticle.cc:435
void setAxes(const Mdouble &a1, const Mdouble &a2, const Mdouble &a3)
Set the axes-lengths to a1, a2 and a3 for this superquadric. We use the super-ellipsoid definition st...
Definition: SuperQuadricParticle.cc:164
Mdouble getExponentEps1() const override
Get the first exponent of this superquadric. We use the super-ellipsoid definition stated in Chapter ...
Definition: SuperQuadricParticle.cc:187
Mdouble computeShape(const LabFixedCoordinates &labFixedCoordinates) const
Compute and get the shape-functiion at the given (lab-fixed) position.
Definition: SuperQuadricParticle.cc:473
void setAxesAndExponents(const Mdouble &a1, const Mdouble &a2, const Mdouble &a3, const Mdouble &eps1, const Mdouble &eps2)
Set the geometrical properties of the superquadrics, namely the axes-lengths a1, a2 and a3,...
Definition: SuperQuadricParticle.cc:133
SmallMatrix< 4, 4 > getJacobianOfContactDetectionObjective(const SmallVector< 4 > &contactPoint, const SuperQuadricParticle *p1, const SuperQuadricParticle *p2) const
Compute and return the derivative of functionThatShouldBecomeZeroForContactDetection,...
Definition: SuperQuadricParticle.cc:623
Mdouble getVolume() const override
Definition: SuperQuadricParticle.cc:204
void setBoundingRadius()
Get the radius of the sphere that fits precisely around the particle.
Definition: SuperQuadricParticle.cc:244
Definition: Vector.h:51
double P
Uniform pressure.
Definition: TwenteMeshGluing.cpp:73