Membrane.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 MEMBRANE_H
27 #define MEMBRANE_H
28 
29 #include <Mercury3D.h>
30 #include <vector>
31 #include <Math/Vector.h>
32 #include <array>
33 
34 #include <Particles/BaseParticle.h>
35 #include <Walls/MeshTriangle.h>
36 
37 
79 class Membrane : public BaseObject
80 {
81 public:
82 
86  struct Edge
87  {
88  // 3 and 4
92  std::array<unsigned int, 2> vertexId;
93 
97  std::array<BaseParticle*, 2> vertex = {{nullptr}};
98 
99  // 1 and 2
103  std::array<unsigned int, 2> faceVertexId;
104 
108  std::array<BaseParticle*, 2> faceVertex = {{nullptr}};
109 
113  std::array<MeshTriangle*, 2> face = {{nullptr}};
114 
118  std::array<Mdouble, 2> faceInitialArea;
119 
123  std::array<Mdouble, 8> uPre;
124 
129 
134 
139 
144 
149 
154 
158  void applyForce(Mdouble Kn, Mdouble critDampCoeff, Mdouble Ke, Mdouble Kd, bool bendingAreaConstant);
159 
163  void applyStretchForce(Mdouble Kn, Mdouble critDampCoeff);
164 
168  void applyBendForce(Mdouble Kn, Mdouble Kd, bool bendingAreaConstant);
169 
173  void calculateUPre(Vec3D& state, Mdouble& length, std::array<Mdouble, 2>& faceArea);
174 
179 
183  void handleParticleRemoval(unsigned int id);
184 
188  void handleParticleAddition(unsigned int id, BaseParticle* p);
189 
193  void checkActive();
194  };
195 
199  Membrane();
200 
201  /*
202  * \brief Constructor to create a membrane with two species
203  */
204  Membrane(ParticleSpecies* membraneSpecies, ParticleSpecies* membraneParticleSpecies);
205 
209  // Membrane(const Membrane& other);
210 
211 
215  ~Membrane() override;
216 
220  // Membrane& operator=(const Membrane& other);
221 
225  Membrane* copy() const;
226 
230  void read(std::istream& is) override;
231 
235  void write(std::ostream& os) const override;
236 
240  std::string getName() const override;
241 
245  void setKnAndCrittDampCoeff(Mdouble Kn, Mdouble critDampCoeff);
246 
253 
257  void setSpringConstant(Mdouble k);
258 
263 
267  void setKeAndKd(Mdouble Ke, Mdouble Kd);
268 
273  void setBendingAreaConstant(bool areaConstant);
274 
278  void setThickness(Mdouble thickness);
279 
280 
281  Mdouble getKn() { return Kn_; }
283  Mdouble getKe() { return Ke_; }
284  Mdouble getKd() { return Kd_; }
287 
291  void setParticleRadius(Mdouble radius);
292 
297 
301  void setDPMBase(DPMBase* dpm) { DPMBase_ = dpm; }
302 
306  DPMBase* getDPMBase() { return DPMBase_; }
307 
311  std::vector<BaseParticle*> getVertexParticles(){ return vertexParticle_; };
312 
316  std::vector<MeshTriangle*> getFaces(){ return face_; };
317 
321  void saveVertexPositions(std::ostream& os);
322 
327  void loadVertexPositions(std::istream& is);
328 
332  void saveAsOFF(unsigned int d);
333 
337  void saveAsSTL(std::string fileName);
338 
342  void loadFromSTL(BaseParticle& p0, std::string fileName);
343 
347  void loadFromSTL(BaseParticle& p0, std::string fileName, Mdouble eps);
348 
352  void buildMesh(BaseParticle& p0, std::vector<Vec3D> vertexPositions, std::vector<unsigned int> edgeVertices, std::vector<unsigned int> faceVertices);
353 
359 
363  void updateEdgeMass();
364 
369 
373  void applyPressure(Mdouble pressure);
374 
378  void handleParticleRemoval(unsigned int id);
379 
383  void handleParticleAddition(unsigned int id, BaseParticle* p);
384 
388  Mdouble getVolume();
389 
390 private:
391 
395  void createVertexParticles(BaseParticle& p0, std::vector<Vec3D> vertexPositions);
396 
401 
405  void updateFaceNeighbors();
406 
410  unsigned int addVertex(std::vector<Vec3D>& vertices, Vec3D pos, Mdouble eps);
411 
415  bool addEdge(std::vector<unsigned int>& edges, std::map<std::pair<unsigned int, unsigned int>, int>& map, unsigned int v0, unsigned int v1);
416 
420  std::vector<BaseParticle*> vertexParticle_;
421 
425  std::vector<unsigned int> vertexParticleId_;
426 
431  unsigned int vertexInitId_;
432 
436  std::vector<MeshTriangle*> face_;
437 
442  std::vector<unsigned int> faceVertices_;
443 
447  std::vector<Edge> edge_;
448 
453 
458 
463 
468 
473 
478 
483 
488 
493 
497  DPMBase* DPMBase_ = nullptr;
498 };
499 
500 
501 #endif
double Mdouble
Definition: GeneralDefine.h:34
It is an abstract base class due to the purely virtual functions declared below. Even if the function...
Definition: BaseObject.h:51
Definition: BaseParticle.h:54
The DPMBase header includes quite a few header files, defining all the handlers, which are essential....
Definition: DPMBase.h:77
A Membrane consists of masses connected by springs.
Definition: Membrane.h:80
void write(std::ostream &os) const override
Writes a Membrane to an output stream, for example a restart file.
Definition: Membrane.cc:204
Mdouble getKe()
Definition: Membrane.h:283
std::vector< MeshTriangle * > getFaces()
Returns a vecter with pointers to the mesh triangles.
Definition: Membrane.h:316
void read(std::istream &is) override
Reads a Membrane from an input stream, for example a restart file.
Definition: Membrane.cc:89
std::vector< unsigned int > vertexParticleId_
Definition: Membrane.h:425
void updateEdgeMass()
Set the correct edge mass by taking the mass from the conencted vertices.
Definition: Membrane.cc:940
void setKnAndCrittDampCoeff(Mdouble Kn, Mdouble critDampCoeff)
Set the parameters needed for the stretching forces.
Definition: Membrane.cc:286
void computeAdditionalForces()
Compute the forces due to the mass spring system.
Definition: Membrane.cc:966
DPMBase * DPMBase_
Definition: Membrane.h:497
unsigned int addVertex(std::vector< Vec3D > &vertices, Vec3D pos, Mdouble eps)
Helper function to check if a given vertex already exists.
Definition: Membrane.cc:563
void updateFaceNeighbors()
Update the faces to have the correct neighbors.
Definition: Membrane.cc:789
Mdouble thickness_
Definition: Membrane.h:477
Mdouble critDampCoeff_
Definition: Membrane.h:462
Mdouble getKn()
Definition: Membrane.h:281
void saveAsOFF(unsigned int d)
Save the Membrane as a .off file.
Definition: Membrane.cc:416
void loadVertexPositions(std::istream &is)
Load the positions of the vertex particles from a stream and apply them to existing particles.
Definition: Membrane.cc:383
void setSpringConstant(Mdouble k)
Set the spring constant of the membrane.
Definition: Membrane.cc:297
ParticleSpecies * membraneSpecies_
Definition: Membrane.h:487
Mdouble Ke_
Definition: Membrane.h:467
ParticleSpecies * membraneParticleSpecies_
Definition: Membrane.h:492
std::vector< BaseParticle * > vertexParticle_
Definition: Membrane.h:420
Mdouble getThickness()
Definition: Membrane.h:285
void setElasticModulusAndThickness(Mdouble E, Mdouble thickness)
Set the elastic modulus and thickness of the membrane \deltails The supplied values are used to calcu...
Definition: Membrane.cc:291
void loadFromSTL(BaseParticle &p0, std::string fileName)
Load the Membrane geometry from a .stl file.
Definition: Membrane.cc:486
std::vector< MeshTriangle * > face_
Definition: Membrane.h:436
void createVertexParticles(BaseParticle &p0, std::vector< Vec3D > vertexPositions)
Handles the actual creation of vertex particles.
Definition: Membrane.cc:705
unsigned int vertexInitId_
Definition: Membrane.h:431
void setThickness(Mdouble thickness)
Set the thickness of the membrane.
Definition: Membrane.cc:317
void handleParticleRemoval(unsigned int id)
Handles the removal of vertex particles from the particle handler.
Definition: Membrane.cc:999
void adjustVertexParticleSize()
Calculates an updated radius for every vertex particle in order to account for the correct mass.
Definition: Membrane.cc:906
Mdouble getKd()
Definition: Membrane.h:284
void setParticleRadius(Mdouble radius)
Set the radius of the vertex particles.
Definition: Membrane.cc:328
std::string getName() const override
Returns the name of the object, here the string "Membrane".
Definition: Membrane.cc:281
std::vector< Edge > edge_
Definition: Membrane.h:447
void buildMesh(BaseParticle &p0, std::vector< Vec3D > vertexPositions, std::vector< unsigned int > edgeVertices, std::vector< unsigned int > faceVertices)
Build the geometry from specified positions and their connectivity.
Definition: Membrane.cc:613
bool bendingAreaConstant_
Definition: Membrane.h:482
Mdouble Kn_
Definition: Membrane.h:457
Membrane * copy() const
Copy method. It calls the copy constructor of this Object.
Definition: Membrane.cc:80
Mdouble particleRadius_
Definition: Membrane.h:452
bool addEdge(std::vector< unsigned int > &edges, std::map< std::pair< unsigned int, unsigned int >, int > &map, unsigned int v0, unsigned int v1)
Helper function to check if a given edge already exists.
Definition: Membrane.cc:588
void saveVertexPositions(std::ostream &os)
save the current positions of the vertex particles to a stream
Definition: Membrane.cc:334
void saveAsSTL(std::string fileName)
Save the Membrane as a .stl file.
Definition: Membrane.cc:453
void initializeEdgeBendingQuantities()
Compute the forces due to the mass spring system.
Definition: Membrane.cc:727
DPMBase * getDPMBase()
Get the stored pointer to DPMBase.
Definition: Membrane.h:306
std::vector< unsigned int > faceVertices_
Definition: Membrane.h:442
Mdouble getBendingAreaConstant()
Definition: Membrane.h:286
Mdouble Kd_
Definition: Membrane.h:472
std::vector< BaseParticle * > getVertexParticles()
Returns a vector with pointers to the vertex particles.
Definition: Membrane.h:311
~Membrane() override
Copy constructor.
Membrane()
Default constructor.
Definition: Membrane.cc:39
Mdouble getCriticalDampingCoefficient()
Definition: Membrane.h:282
Mdouble getParticleRadius()
Returns the radius of the vertex particles.
Definition: Membrane.h:296
Mdouble getVolume()
Calculate the volume of the membrane.
Definition: Membrane.cc:1066
void handleParticleAddition(unsigned int id, BaseParticle *p)
Handles the addition of vertex particles to the particle handler.
Definition: Membrane.cc:1033
void setCriticalDampingCoefficient(Mdouble coeff)
Set damping coefficient for the distance springs.
Definition: Membrane.cc:303
void setDPMBase(DPMBase *dpm)
Set a pointer to DPMBase.
Definition: Membrane.h:301
void setKeAndKd(Mdouble Ke, Mdouble Kd)
Set the parameters needed for the bending forces.
Definition: Membrane.cc:311
void setBendingAreaConstant(bool areaConstant)
Definition: Membrane.cc:323
void applyPressure(Mdouble pressure)
Apply a surface pressure to the membrane.
Definition: Membrane.cc:978
Definition: ParticleSpecies.h:37
Definition: Vector.h:51
double E
Elastic modulus.
Definition: TwenteMeshGluing.cpp:64
Definition: Membrane.h:87
void calculateUPre(Vec3D &state, Mdouble &length, std::array< Mdouble, 2 > &faceArea)
Calculate some prefactors for the bending penalty.
Definition: Membrane.cc:1200
std::array< Mdouble, 8 > uPre
Definition: Membrane.h:123
bool isBendActive
Definition: Membrane.h:153
Mdouble effectiveMass
Definition: Membrane.h:143
std::array< BaseParticle *, 2 > vertex
Definition: Membrane.h:97
Mdouble initialLength
Definition: Membrane.h:133
void applyForce(Mdouble Kn, Mdouble critDampCoeff, Mdouble Ke, Mdouble Kd, bool bendingAreaConstant)
Calculates and applies all neccesary forces.
Definition: Membrane.cc:1088
void handleParticleAddition(unsigned int id, BaseParticle *p)
Handle the particle addition.
Definition: Membrane.cc:1280
void applyBendForce(Mdouble Kn, Mdouble Kd, bool bendingAreaConstant)
Apply a force due to bending.
Definition: Membrane.cc:1136
Mdouble initialSineHalfTheta
Definition: Membrane.h:138
std::array< BaseParticle *, 2 > faceVertex
Definition: Membrane.h:108
std::array< Mdouble, 2 > faceInitialArea
Definition: Membrane.h:118
void applyStretchForce(Mdouble Kn, Mdouble critDampCoeff)
Apply the force due to stretching only.
Definition: Membrane.cc:1099
std::array< unsigned int, 2 > faceVertexId
Definition: Membrane.h:103
Mdouble getSineHalfTheta()
Calculate the sine of half the bending angle.
Definition: Membrane.cc:1227
void checkActive()
check if the edge should calculate bending or stretch forces
Definition: Membrane.cc:1302
std::array< MeshTriangle *, 2 > face
Definition: Membrane.h:113
bool isStretchActive
Definition: Membrane.h:148
std::array< unsigned int, 2 > vertexId
Definition: Membrane.h:92
void handleParticleRemoval(unsigned int id)
handle the partical removal
Definition: Membrane.cc:1254
Vec3D initialState
Definition: Membrane.h:128