MercuryDPM  0.10
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
TangentialSpringParticle.cc
Go to the documentation of this file.
1 //Copyright (c) 2013-2014, 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 
27 #include "../ParticleHandler.h"
28 
30 {
32 }
33 
35 {
36  //std::cout<<"TSP copy constructor"<<std::endl;
38 }
39 
41 {
42 }
43 
46 {
47  //std::cout<<"TangentialSpringParticle copy"<<std::endl;
48  return new TangentialSpringParticle(*this);
49 }
50 
52 {
53  for(std::vector<CTangentialSpring>::iterator it = tangentialSprings.begin(); it!=tangentialSprings.end();it++)
54  {
55  it->reverse();
56  }
57 }
58 
60 
61 void TangentialSpringParticle::print(std::ostream& os) const
62 {
63  os<<"TSP ";
65  os<<" "<<tangentialSprings;
66 }
67 
68 void TangentialSpringParticle::read(std::istream& is)
69 {
70  std::string dummy;
71  is >> dummy;
73  is >> tangentialSprings;
74 }
75 
77 {
79  //Check if some TangentialSpring information has to be moved
80  for(unsigned int i=0;i<get_TangentialSprings().size();)
81  {
83  if (get_TangentialSprings()[i].pParticle>newPos)
84  {
85  int otherId=get_TangentialSprings()[i].pParticle;
86  TangentialSpringParticle* otherTSParticle=dynamic_cast<TangentialSpringParticle*>(getHandler()->getObject(otherId));
87  //Copy Tangentalspring to new location
89  CTS.pParticle=newPos;
90  CTS.reverse();
91  otherTSParticle->get_TangentialSprings().push_back(CTS);
92 
93  //Remove now unused tangentialspring
95  get_TangentialSprings().pop_back();
96  }
97  else
98  {
99  i++;
100  }
101  }
102 }
103 
104 void TangentialSpringParticle::oldRead(std::istream& is,std::string& x)
105 {
107  position.X=atof(x.c_str());
109  is >> position.Y >> position.Z >> velocity >> radius >> angle >> angularVelocity >> invMass >> invInertia >> tangentialSprings;
110  set_Position(position);
111  set_Velocity(velocity);
112  set_Angle(angle);
113  set_AngularVelocity(angularVelocity);
114  set_Radius(radius);
115  if (invMass) set_Mass(1.0/invMass); else set_Mass(1e20);
116  if (invInertia) set_inertia(1.0/invInertia); else set_inertia(1e20);
117 }
118 
Vec3D velocity
Current particle position.
Definition: BaseParticle.h:193
void set_inertia(const Mdouble _new)
TangentialSpringParticle()
Basic Particle contructors, creates an Particle at (0,0,0) with radius, mass and inertia equal to 1...
Mdouble X
Definition: Vector.h:44
virtual ~TangentialSpringParticle()
Particle destructor, needs to be implemented and checked if it removes tangential spring information...
T * getObject(const unsigned int id) const
Gets a pointer to the Object at the specified index in the BaseHandler.
Definition: BaseHandler.h:176
Mdouble radius
Inverse Particle inverse inertia (for computation optimization)
Definition: BaseParticle.h:185
virtual void print(std::ostream &os) const
Particle print function, which accepts an os std::stringstream as input.
int pParticle
A pointer to the particle in contact; NULL if the contact is with a wall (The other particle is the p...
Vec3D angularVelocity
Current angular position.
Definition: BaseParticle.h:191
Vec3D angle
non changing identifier of particle
Definition: BaseParticle.h:190
ParticleHandler * getHandler() const
CTangentialSprings & get_TangentialSprings()
Mdouble invInertia
Particle inertia.
Definition: BaseParticle.h:184
void oldRead(std::istream &is, std::string &x)
virtual void read(std::istream &is)
Particle read function, which accepts an os std::stringstream as input.
double Mdouble
Definition: ExtendedMath.h:33
Member variable of #Particle storing all tangential springs of particle PI with contacting particles...
Stores the tangential spring needed to compute a tangential elastic force between particles PI...
CTangentialSprings tangentialSprings
void reset()
Resets the tangential springs.
Mdouble Y
Definition: Vector.h:44
void set_Radius(const Mdouble _new)
void set_Velocity(const Vec3D &_new)
Mdouble invMass
Particle mass.
Definition: BaseParticle.h:182
Vec3D position
Current angular velocity.
Definition: BaseParticle.h:192
void set_AngularVelocity(const Vec3D &_new)
void set_Angle(const Vec3D &_new)
virtual void moveInHandler(int newPos)
void print(std::ostream &os) const
Particle print function, which accepts an os std::stringstream as input.
Implementation of a 3D vector (by Vitaliy).
Definition: Vector.h:40
virtual TangentialSpringParticle * copy() const
Particle copy method. It calls to copy contrustor of this Particle, usefull for polymorfism.
void read(std::istream &is)
Particle read function, which accepts an os std::stringstream as input.
void set_Position(const Vec3D &_new)
Mdouble Z
Definition: Vector.h:44
void set_Mass(const Mdouble _new)