Membrane Class Reference

A Membrane consists of masses connected by springs. More...

#include <Membrane.h>

+ Inheritance diagram for Membrane:

Classes

struct  Edge
 

Public Member Functions

 Membrane ()
 Default constructor. More...
 
 Membrane (ParticleSpecies *membraneSpecies, ParticleSpecies *membraneParticleSpecies)
 
 ~Membrane () override
 Copy constructor. More...
 
Membranecopy () const
 Copy method. It calls the copy constructor of this Object. More...
 
void read (std::istream &is) override
 Reads a Membrane from an input stream, for example a restart file. More...
 
void write (std::ostream &os) const override
 Writes a Membrane to an output stream, for example a restart file. More...
 
std::string getName () const override
 Returns the name of the object, here the string "Membrane". More...
 
void setKnAndCrittDampCoeff (Mdouble Kn, Mdouble critDampCoeff)
 Set the parameters needed for the stretching forces. More...
 
void setElasticModulusAndThickness (Mdouble E, Mdouble thickness)
 Set the elastic modulus and thickness of the membrane \deltails The supplied values are used to calculate the spring constant of the membrane. More...
 
void setSpringConstant (Mdouble k)
 Set the spring constant of the membrane. More...
 
void setCriticalDampingCoefficient (Mdouble coeff)
 Set damping coefficient for the distance springs. More...
 
void setKeAndKd (Mdouble Ke, Mdouble Kd)
 Set the parameters needed for the bending forces. More...
 
void setBendingAreaConstant (bool areaConstant)
 
void setThickness (Mdouble thickness)
 Set the thickness of the membrane. More...
 
Mdouble getKn ()
 
Mdouble getCriticalDampingCoefficient ()
 
Mdouble getKe ()
 
Mdouble getKd ()
 
Mdouble getThickness ()
 
Mdouble getBendingAreaConstant ()
 
void setParticleRadius (Mdouble radius)
 Set the radius of the vertex particles. More...
 
Mdouble getParticleRadius ()
 Returns the radius of the vertex particles. More...
 
void setDPMBase (DPMBase *dpm)
 Set a pointer to DPMBase. More...
 
DPMBasegetDPMBase ()
 Get the stored pointer to DPMBase. More...
 
std::vector< BaseParticle * > getVertexParticles ()
 Returns a vector with pointers to the vertex particles. More...
 
std::vector< MeshTriangle * > getFaces ()
 Returns a vecter with pointers to the mesh triangles. More...
 
void saveVertexPositions (std::ostream &os)
 save the current positions of the vertex particles to a stream More...
 
void loadVertexPositions (std::istream &is)
 Load the positions of the vertex particles from a stream and apply them to existing particles. More...
 
void saveAsOFF (unsigned int d)
 Save the Membrane as a .off file. More...
 
void saveAsSTL (std::string fileName)
 Save the Membrane as a .stl file. More...
 
void loadFromSTL (BaseParticle &p0, std::string fileName)
 Load the Membrane geometry from a .stl file. More...
 
void loadFromSTL (BaseParticle &p0, std::string fileName, Mdouble eps)
 Load the Membrane geometry from a .stl file. More...
 
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. More...
 
void adjustVertexParticleSize ()
 Calculates an updated radius for every vertex particle in order to account for the correct mass. More...
 
void updateEdgeMass ()
 Set the correct edge mass by taking the mass from the conencted vertices. More...
 
void computeAdditionalForces ()
 Compute the forces due to the mass spring system. More...
 
void applyPressure (Mdouble pressure)
 Apply a surface pressure to the membrane. More...
 
void handleParticleRemoval (unsigned int id)
 Handles the removal of vertex particles from the particle handler. More...
 
void handleParticleAddition (unsigned int id, BaseParticle *p)
 Handles the addition of vertex particles to the particle handler. More...
 
Mdouble getVolume ()
 Calculate the volume of the membrane. More...
 
- 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 Member Functions

void createVertexParticles (BaseParticle &p0, std::vector< Vec3D > vertexPositions)
 Handles the actual creation of vertex particles. More...
 
void initializeEdgeBendingQuantities ()
 Compute the forces due to the mass spring system. More...
 
void updateFaceNeighbors ()
 Update the faces to have the correct neighbors. More...
 
unsigned int addVertex (std::vector< Vec3D > &vertices, Vec3D pos, Mdouble eps)
 Helper function to check if a given vertex already exists. More...
 
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. More...
 

Private Attributes

std::vector< BaseParticle * > vertexParticle_
 
std::vector< unsigned int > vertexParticleId_
 
unsigned int vertexInitId_
 
std::vector< MeshTriangle * > face_
 
std::vector< unsigned int > faceVertices_
 
std::vector< Edgeedge_
 
Mdouble particleRadius_
 
Mdouble Kn_
 
Mdouble critDampCoeff_
 
Mdouble Ke_
 
Mdouble Kd_
 
Mdouble thickness_
 
bool bendingAreaConstant_
 
ParticleSpeciesmembraneSpecies_ = nullptr
 
ParticleSpeciesmembraneParticleSpecies_ = nullptr
 
DPMBaseDPMBase_ = nullptr
 

Detailed Description

A Membrane consists of masses connected by springs.

This class assumes, that the masses are connected in a triangular mesh. Additionally, for correct calculation of the spring constant for the distance spring, a hexagonal unit cell is assumed.

This class creates the particles required for the mass spring system. For contact detections, further triangular wall elements are created

Please make sure to provide appropriate a species for the particles as well as the wall elements (e.g giving the particles a species so that they do not interact with other particles in the simulation).

An instance of the membrane may be initialized with the following lines, e.g. by using a STL file:

// Initialize the membrane and set neccesary quantities
membrane.setDPMBase(this);
membrane.setKeAndKd(Ke_, Kd_);
// Set the (approximate) radius of the vertex particles
membrane.setParticleRadius(membraneParticleRadius_);
// Create an initial particle for the vertex particles
p0.setVelocity(Vec3D(0.0, 0.0, 0.0));
// Create the system using the specified particle and stl file
membrane_.readFromSTL(p0, "geometry.stl");
void setVelocity(const Vec3D &velocity)
set the velocity of the BaseInteractable.
Definition: BaseInteractable.cc:350
virtual void setRadius(Mdouble radius)
Sets the particle's radius_ (and adjusts the mass_ accordingly, based on the particle's species)
Definition: BaseParticle.cc:553
void setSpecies(const ParticleSpecies *species)
Definition: BaseParticle.cc:818
A Membrane consists of masses connected by springs.
Definition: Membrane.h:80
void setKnAndCrittDampCoeff(Mdouble Kn, Mdouble critDampCoeff)
Set the parameters needed for the stretching forces.
Definition: Membrane.cc:286
Mdouble critDampCoeff_
Definition: Membrane.h:462
ParticleSpecies * membraneSpecies_
Definition: Membrane.h:487
Mdouble Ke_
Definition: Membrane.h:467
ParticleSpecies * membraneParticleSpecies_
Definition: Membrane.h:492
void setParticleRadius(Mdouble radius)
Set the radius of the vertex particles.
Definition: Membrane.cc:328
Mdouble Kn_
Definition: Membrane.h:457
Mdouble Kd_
Definition: Membrane.h:472
Membrane()
Default constructor.
Definition: Membrane.cc:39
Mdouble getParticleRadius()
Returns the radius of the vertex particles.
Definition: Membrane.h:296
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
A spherical particle is the most simple particle used in MercuryDPM.
Definition: SphericalParticle.h:37
Definition: Vector.h:51

Note that the method computeAdditionalForces needs to be called from the overridden DPMBase::computeAdditionalForces every timestep for the Membrane to work.

As an example please also have a look at the Driver code MembraneDemo

Constructor & Destructor Documentation

◆ Membrane() [1/2]

Membrane::Membrane ( )

Default constructor.

40 {
41  logger(DEBUG, "Membrane() constructed.");
42  vertexInitId_ = 0;
43  particleRadius_ = 0;
44  Kn_ = 0;
45  critDampCoeff_ = 0;
46  Ke_ = 0;
47  Kd_ = 0;
48  thickness_ = 0;
50 }
Logger< MERCURYDPM_LOGLEVEL > logger("MercuryKernel")
Definition of different loggers with certain modules. A user can define its own custom logger here.
@ DEBUG
Mdouble thickness_
Definition: Membrane.h:477
unsigned int vertexInitId_
Definition: Membrane.h:431
bool bendingAreaConstant_
Definition: Membrane.h:482
Mdouble particleRadius_
Definition: Membrane.h:452

References bendingAreaConstant_, critDampCoeff_, DEBUG, Kd_, Ke_, Kn_, logger, particleRadius_, thickness_, and vertexInitId_.

Referenced by copy().

◆ Membrane() [2/2]

Membrane::Membrane ( ParticleSpecies membraneSpecies,
ParticleSpecies membraneParticleSpecies 
)
52  : Membrane()
53 {
54  membraneSpecies_ = membraneSpecies;
55  membraneParticleSpecies_ = membraneParticleSpecies;
56 }

References membraneParticleSpecies_, and membraneSpecies_.

◆ ~Membrane()

Membrane::~Membrane ( )
override

Copy constructor.

Destructor.

Member Function Documentation

◆ addEdge()

bool Membrane::addEdge ( std::vector< unsigned int > &  edges,
std::map< std::pair< unsigned int, unsigned int >, int > &  map,
unsigned int  v0,
unsigned int  v1 
)
private

Helper function to check if a given edge already exists.

Parameters
[in]edgesThe vector, where the new edge should be added
[in]mapA map used to dtermin, if the edge already exists
[in]v0The index of the first point of the edge
[in]v0The index of the second point of the edge
Returns
True if a new edge was added, false if not.

This function adds a given edge to the vector edges, if it is not already contained.

589 {
590  if ( map.find(std::make_pair(v0, v1)) == map.end() ){
591  map[std::make_pair(v0,v1)] = 1;
592  map[std::make_pair(v1,v0)] = 1;
593 
594  edges.push_back(v0);
595  edges.push_back(v1);
596  return true;
597  }
598  return false;
599 }

Referenced by loadFromSTL().

◆ addVertex()

unsigned int Membrane::addVertex ( std::vector< Vec3D > &  vertices,
Vec3D  pos,
Mdouble  eps 
)
private

Helper function to check if a given vertex already exists.

Parameters
[in]verticesThe vector, where the new position should be added
[in]posThe new position to add to the vector
[in]epsThe distance between two points below which they are assumed to be the same.
Returns
the index of the position in the vector vertices.

This function adds a given position pos to the vector vertices, if it is not already contained.

564 {
565  unsigned int i;
566  for (i=0; i<vertices.size(); i++)
567  {
568  if ((vertices[i]-pos).getLength() < eps)
569  {
570  return i;
571  }
572  }
573 
574  vertices.push_back(pos);
575 
576  return vertices.size() - 1;
577 }
const std::complex< Mdouble > i
Definition: ExtendedMath.h:51

References constants::i.

Referenced by loadFromSTL().

◆ adjustVertexParticleSize()

void Membrane::adjustVertexParticleSize ( )

Calculates an updated radius for every vertex particle in order to account for the correct mass.

Calculates an updated radius for every vertex particle in order to account for the correct mass.

907 {
908  std::vector<unsigned int> numberTriangles(vertexParticle_.size(), 0);
909  std::vector<Mdouble> massParticle(vertexParticle_.size(), 0.0);
910 
911  Mdouble mass;
912  for (auto f: face_)
913  {
914  // 3 Nodes per triangle
915  mass = 1/f->getInvMass()/3;
916  for (auto id: f->getVertexIds())
917  {
918  numberTriangles[id-vertexInitId_] += 1;
919  massParticle[id-vertexInitId_] += mass;
920  }
921  }
922 
923  unsigned int i;
924  Mdouble volume;
925  Mdouble radius;
927  for (i=0; i<vertexParticle_.size(); i++)
928  {
929  mass = massParticle[i]/numberTriangles[i];
930  volume = mass / density;
931  radius = std::cbrt(3.0/4.0 * volume/constants::pi);
932 
933  vertexParticle_[i]->setRadius(radius);
934  }
935 }
double Mdouble
Definition: GeneralDefine.h:34
std::vector< BaseParticle * > vertexParticle_
Definition: Membrane.h:420
std::vector< MeshTriangle * > face_
Definition: Membrane.h:436
Mdouble getDensity() const
Allows density_ to be accessed.
Definition: ParticleSpecies.cc:118
const Mdouble pi
Definition: ExtendedMath.h:45

References face_, ParticleSpecies::getDensity(), constants::i, membraneParticleSpecies_, constants::pi, vertexInitId_, and vertexParticle_.

◆ applyPressure()

void Membrane::applyPressure ( Mdouble  pressure)

Apply a surface pressure to the membrane.

Parameters
[in]pressureValue of the surface pressure acing on the membrane [Pa]

This function applies a given pressure to the surface by calling the function applyPressure of a MeshTriangle.

978  {
979  // Set the pressureforce to 0
980  logger(DEBUG,"Calculating pressure force.");
981 
982  // #pragma omp parallel num_threads(getNumberOfOMPThreads())
983  {
984  // This should also be OMP parralelizable, because addForce pays attention
985  // to OMP
986  // #pragma omp for
987  for (unsigned int i=0; i < face_.size(); i++ ){
988  face_[i]->applyPressure(pressure);
989  }
990  }
991 }

References DEBUG, face_, constants::i, and logger.

◆ buildMesh()

void Membrane::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.

Parameters
[in]p0BaseParticle that will be used for the vertex particles
[in]vertexPositionsThe position of the vertices
[in]edgeVerticesA vector of indices, where two consecutive values define the vertex particles connected by the edge
[in]faceVerticesA vector of indices, where three consecutive values define the vertex particles making up the triangle face

This function creates the actual mesh including the vertex particles and wall elements. The mass of the triangles is set to match the actual weight of the cooresponding triangle. The density of the particle species is set, so that the mass of a vertex particle is equal to 5/3 times the mass of a triangle.

614 {
615  unsigned int triangleCount = faceVertices.size() / 3;
616  unsigned int edgeCount = edgeVertices.size() / 2;
617  // #ifdef MERCURYDPM_USE_MPI
618  // if (PROCESSOR_ID==0)
619  // #endif
620  // logger(INFO, "Creating membrane with % vertices and % faces.", vertexPositions.size(), triangleCount);
621 
622  // Some helper variables
623  unsigned int i;
624  Vec3D v;
625 
626  // Calculate the density of the Membrane particles to approximately get the membrane weight
627  // Note this should be changed. But for the beginning it is a way of setting an approx. okej density for the particles.
628  i = 0;
629  Mdouble membraneDensity = membraneSpecies_->getDensity();
630 
631 
632  // Create particles corresponding to the Vertex Positions
633  createVertexParticles(p0, vertexPositions);
634 
635  Mdouble averageTriangleMass = 0;
636  Mdouble triangleMass = 0;
637  // Create all faces with their initial positions
638  face_.reserve(triangleCount);
639  for (i = 0; i < triangleCount; i++)
640  {
641  // It seems copyAndAddObject does not actually copy in MPI mode
642  MeshTriangle f;
643  f.setVertices(vertexPositions[faceVertices[3*i]],
644  vertexPositions[faceVertices[3*i+1]],
645  vertexPositions[faceVertices[3*i+2]]);
646  f.setVertexIds(vertexParticleId_[faceVertices[3*i]],
647  vertexParticleId_[faceVertices[3*i+1]],
648  vertexParticleId_[faceVertices[3*i+2]]);
650 
651  triangleMass = f.getArea()*thickness_*membraneDensity;
652  averageTriangleMass += triangleMass;
653  f.setMass(triangleMass);
654 
656 
657  face_.push_back(getDPMBase()->wallHandler.copyAndAddObject(f));
658 
659  }
660 
661  averageTriangleMass /= triangleCount;
662  Mdouble particleDensity = 5.0/3.0 * averageTriangleMass/(particleRadius_*particleRadius_*particleRadius_*4/3.0*constants::pi);
663  membraneParticleSpecies_->setDensity(particleDensity);
664 
665  // Create the edges
666  Mdouble invMass1, invMass2;
667 
669  for (i=0; i<edgeCount; i++){
670  Edge e;
671  e.vertexId[0] = vertexParticleId_[edgeVertices[2*i+0]];
672  e.vertexId[1] = vertexParticleId_[edgeVertices[2*i+1]];
673 
674  invMass1 = vertexInvMass;
675  invMass2 = vertexInvMass;
676 
677  e.vertex[0] = getDPMBase()->particleHandler.getObjectById(e.vertexId[0]);
678  e.vertex[1] = getDPMBase()->particleHandler.getObjectById(e.vertexId[1]);
679 
680  // e.vertex[1] = vertexParticle_[edgeVertices[2*i+1]];
681  e.initialState = vertexPositions[edgeVertices[2*i+1]]-vertexPositions[edgeVertices[2*i+0]];
682  e.initialLength = e.initialState.getLength();
683  e.effectiveMass = 1/(invMass1 + invMass2);
684 
685  // e.checkActive();
686  edge_.push_back(e);
687  }
688 
690 
691  for (i=0; i<edge_.size(); i++){
692  edge_[i].checkActive();
693  }
695  for (i=0; i<edge_.size(); i++)
696  {
697  edge_[i].calculateUPre(edge_[i].initialState, edge_[i].initialLength, edge_[i].faceInitialArea);
698  }
699 }
T * getObjectById(const unsigned int id)
Gets a pointer to the Object at the specified index in the BaseHandler.
Definition: BaseHandler.h:565
void setSpecies(const ParticleSpecies *species)
Defines the species of the current wall.
Definition: BaseWall.cc:169
ParticleHandler particleHandler
An object of the class ParticleHandler, contains the pointers to all the particles created.
Definition: DPMBase.h:1437
std::vector< unsigned int > vertexParticleId_
Definition: Membrane.h:425
void updateFaceNeighbors()
Update the faces to have the correct neighbors.
Definition: Membrane.cc:789
void createVertexParticles(BaseParticle &p0, std::vector< Vec3D > vertexPositions)
Handles the actual creation of vertex particles.
Definition: Membrane.cc:705
std::vector< Edge > edge_
Definition: Membrane.h:447
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
MeshTriangle implements a triangle whose vertex positions are defined by three particles.
Definition: MeshTriangle.h:75
void setMass(Mdouble mass)
Definition: MeshTriangle.cc:681
Mdouble getArea() const
Returns the area of the triangle.
Definition: MeshTriangle.h:147
void actionsAfterParticleGhostUpdate() override
actionsPerformed after the position update of (ghost-) particles.
Definition: MeshTriangle.cc:599
void setVertices(Vec3D A, Vec3D B, Vec3D C)
Sets member variables such that the wall represents a triangle with vertices A, B,...
Definition: MeshTriangle.cc:456
void setVertexIds(unsigned int i, unsigned int j, unsigned int k)
sets the ids of the vertex particles. Calls retrieveVertexParticles.
Definition: MeshTriangle.cc:481
void setDensity(Mdouble density)
Definition: ParticleSpecies.cc:108

References MeshTriangle::actionsAfterParticleGhostUpdate(), createVertexParticles(), edge_, Membrane::Edge::effectiveMass, face_, MeshTriangle::getArea(), ParticleSpecies::getDensity(), getDPMBase(), Vec3D::getLength(), BaseHandler< T >::getObjectById(), constants::i, initializeEdgeBendingQuantities(), Membrane::Edge::initialLength, Membrane::Edge::initialState, membraneParticleSpecies_, membraneSpecies_, DPMBase::particleHandler, particleRadius_, constants::pi, ParticleSpecies::setDensity(), MeshTriangle::setMass(), BaseWall::setSpecies(), MeshTriangle::setVertexIds(), MeshTriangle::setVertices(), thickness_, updateFaceNeighbors(), Membrane::Edge::vertex, Membrane::Edge::vertexId, and vertexParticleId_.

Referenced by loadFromSTL().

◆ computeAdditionalForces()

void Membrane::computeAdditionalForces ( )

Compute the forces due to the mass spring system.

This function iterates through all edges to calculate the forces due to the mass spring system. Ideally this function is called from the computeAllForces function of DPMBase.

966  {
967  #pragma omp parallel for schedule(static) num_threads(getDPMBase()->getNumberOfOMPThreads())
968  for (int i=0; i < edge_.size(); i++){
970  }
971 }

References bendingAreaConstant_, critDampCoeff_, edge_, constants::i, Kd_, Ke_, and Kn_.

Referenced by MembraneDemo::computeAdditionalForces(), and MembraneSelfTest::computeAdditionalForces().

◆ copy()

Membrane * Membrane::copy ( ) const

Copy method. It calls the copy constructor of this Object.

Copy assignment operator.

  \param[in] other The Membrane that must be copied.
 &zwj;/

Membrane::Membrane(const Membrane& other) : BaseObject(other) { face_ = other.face_; membraneSpecies_ = other.membraneSpecies_; DPMBase_ = other.DPMBase_;

logger(DEBUG, "Membrane(Membrane&) constructed."); }

Membrane::~Membrane() { logger(DEBUG, "~Membrane() has been called."); }

/*!

Returns
pointer to a Membrane object allocated using new.
81 {
82  return new Membrane(*this);
83 }

References Membrane().

◆ createVertexParticles()

void Membrane::createVertexParticles ( BaseParticle p0,
std::vector< Vec3D vertexPositions 
)
private

Handles the actual creation of vertex particles.

Parameters
[in]p0The template particle used to create the vertex particles

This function creates a particle for every position in the vector vertexPositions

706 {
707 
709 
710  Vec3D pos;
711  for (unsigned int i = 0; i < vertexPositions.size(); i++){
712  pos = vertexPositions[i];
713  // Set the position
714  p0.setPosition(pos);
715 
716  // Add the particle and save a reference
718  vertexParticleId_.push_back(vertexInitId_+i);
719  vertexParticle_.push_back(getDPMBase()->particleHandler.getObjectById(vertexInitId_+i));
720  }
721 
722 }
std::enable_if<!std::is_pointer< U >::value, U * >::type copyAndAddObject(const U &object)
Creates a copy of a Object and adds it to the BaseHandler.
Definition: BaseHandler.h:379
virtual void setPosition(const Vec3D &position)
Sets the position of this BaseInteractable.
Definition: BaseInteractable.h:239
unsigned int getNumberOfRealObjects() const
Returns the number of real objects (on all processors)
Definition: ParticleHandler.cc:1304

References BaseHandler< T >::copyAndAddObject(), getDPMBase(), ParticleHandler::getNumberOfRealObjects(), constants::i, DPMBase::particleHandler, BaseInteractable::setPosition(), vertexInitId_, vertexParticle_, and vertexParticleId_.

Referenced by buildMesh().

◆ getBendingAreaConstant()

Mdouble Membrane::getBendingAreaConstant ( )
inline
286 { return bendingAreaConstant_; }

References bendingAreaConstant_.

◆ getCriticalDampingCoefficient()

Mdouble Membrane::getCriticalDampingCoefficient ( )
inline

◆ getDPMBase()

DPMBase* Membrane::getDPMBase ( )
inline

Get the stored pointer to DPMBase.

306 { return DPMBase_; }
DPMBase * DPMBase_
Definition: Membrane.h:497

References DPMBase_.

Referenced by buildMesh(), createVertexParticles(), initializeEdgeBendingQuantities(), read(), and saveVertexPositions().

◆ getFaces()

std::vector<MeshTriangle*> Membrane::getFaces ( )
inline

Returns a vecter with pointers to the mesh triangles.

316 { return face_; };

References face_.

◆ getKd()

Mdouble Membrane::getKd ( )
inline
284 { return Kd_; }

References Kd_.

Referenced by MembraneDemo::setUpMembrane(), and MembraneSelfTest::setUpMembrane().

◆ getKe()

Mdouble Membrane::getKe ( )
inline
283 { return Ke_; }

References Ke_.

Referenced by MembraneDemo::setUpMembrane(), and MembraneSelfTest::setUpMembrane().

◆ getKn()

Mdouble Membrane::getKn ( )
inline
281 { return Kn_; }

References Kn_.

Referenced by MembraneDemo::setUpMembrane(), and MembraneSelfTest::setUpMembrane().

◆ getName()

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

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

Returns
The string "Membrane".

Implements BaseObject.

282 {
283  return "Membrane";
284 }

◆ getParticleRadius()

Mdouble Membrane::getParticleRadius ( )
inline

Returns the radius of the vertex particles.

296 { return particleRadius_; }

References particleRadius_.

Referenced by MembraneDemo::setUpMembrane(), and MembraneSelfTest::setUpMembrane().

◆ getThickness()

Mdouble Membrane::getThickness ( )
inline
285 { return thickness_; }

References thickness_.

◆ getVertexParticles()

std::vector<BaseParticle*> Membrane::getVertexParticles ( )
inline

Returns a vector with pointers to the vertex particles.

311 { return vertexParticle_; };

References vertexParticle_.

Referenced by MembraneDemo::fixMembraneEdges(), and MembraneSelfTest::fixMembraneEdges().

◆ getVolume()

Mdouble Membrane::getVolume ( )

Calculate the volume of the membrane.

Returns
volume The enclosed volume of the membrane

This function computes the current volume enclosed by the membrane. Note: This function uses the Divergence theorem to integrate along the surface. The volume is therefore only accurate, if the surface is closed.

1067 {
1068  Mdouble volume = 0;
1069 
1070  for (auto f: face_)
1071  {
1072  volume += Vec3D::dot(f->getPosition(), f->getFaceNormal()) * f->getArea();
1073  }
1074 
1075  return volume/3;
1076 }
static Mdouble dot(const Vec3D &a, const Vec3D &b)
Calculates the dot product of two Vec3D: .
Definition: Vector.cc:76

References Vec3D::dot(), and face_.

◆ handleParticleAddition()

void Membrane::handleParticleAddition ( unsigned int  id,
BaseParticle p 
)

Handles the addition of vertex particles to the particle handler.

Parameters
[in]idThe id of the added particle
[in]pPointer to the added particle

This function handles the addition of a particle with the given id to the particle handler. If the particle is contained in an edge, the edge is potentially. reactivated. This function is potentially usable for MPI parallelization

1034 {
1035  for (auto& e: edge_)
1036  {
1037  e.handleParticleAddition(id, p);
1038  }
1039 
1040  // Now find the corresponding vertexParticle_
1041  unsigned i = id - vertexInitId_;
1042  if (vertexParticleId_.size() > i && vertexParticleId_[i] == id)
1043  {
1044  vertexParticle_[i] = p;
1045  } else
1046  {
1047  // Need to search for it
1048  for (i = 0; i < vertexParticleId_.size(); i++)
1049  {
1050  if (vertexParticleId_[i]==id)
1051  {
1052  // now i is the correct index
1053  vertexParticle_[i] = p;
1054  break;
1055  }
1056  }
1057  }
1058 }

References edge_, constants::i, vertexInitId_, vertexParticle_, and vertexParticleId_.

◆ handleParticleRemoval()

void Membrane::handleParticleRemoval ( unsigned int  id)

Handles the removal of vertex particles from the particle handler.

Parameters
[in]idThe id of the removed particle

This function handles the removal of a particle with the given id from the particle handler. If the particle was contained in an edge, the edge is disabled. This function is potentially usable for MPI parallelization

1000 {
1001  for (auto& e: edge_)
1002  {
1003  e.handleParticleRemoval(id);
1004  }
1005 
1006  // Now find the corresponding vertexParticle_
1007  unsigned i = id - vertexInitId_;
1008  if (vertexParticleId_[i] == id)
1009  {
1010  vertexParticle_[i] = nullptr;
1011  } else
1012  {
1013  // Need to search for it
1014  for (i = 0; i < vertexParticleId_.size(); i++)
1015  {
1016  if (vertexParticleId_[i]==id)
1017  {
1018  // now i is the coorect index
1019  vertexParticle_[i] = nullptr;
1020  break;
1021  }
1022  }
1023  }
1024 }

References edge_, constants::i, vertexInitId_, vertexParticle_, and vertexParticleId_.

◆ initializeEdgeBendingQuantities()

void Membrane::initializeEdgeBendingQuantities ( )
private

Compute the forces due to the mass spring system.

This function initializes the quantities and pointer needed on the edges for calculating the forces due to the mass spring system.

728 {
729  //set neighbours
730  unsigned int i, j;
731 
732  for (auto& e: edge_)
733  {
734  e.face[0] = nullptr;
735  e.face[1] = nullptr;
736  e.faceVertexId[0] = UINT_MAX;
737  e.faceVertexId[1] = UINT_MAX;
738  e.faceVertex[0] = nullptr;
739  e.faceVertex[1] = nullptr;
740  e.faceInitialArea[0] = 0;
741  e.faceInitialArea[1] = 0;
742  e.initialSineHalfTheta = 0;
743  i = 0;
744 
745  unsigned int edgeId0 = e.vertexId[0];
746  unsigned int edgeId1 = e.vertexId[1];
747 
748  for (unsigned int f=0; f<face_.size(); f++)
749  {
750  MeshTriangle* face0 = face_[f];
751  std::array<unsigned int, 3> particleIds = face0->getVertexIds();
752  j = UINT_MAX;
753 
754  if ( particleIds[0] == edgeId0)
755  {
756  if (particleIds[1] == edgeId1) j = particleIds[2];
757  if (particleIds[2] == edgeId1) j = particleIds[1];
758  }
759  else if ( particleIds[0] == edgeId1)
760  {
761  if (particleIds[1] == edgeId0) j = particleIds[2];
762  if (particleIds[2] == edgeId0) j = particleIds[1];
763  }
764  else if ( (particleIds[1] == edgeId1 && particleIds[2] == edgeId0)
765  ||(particleIds[1] == edgeId0 && particleIds[2] == edgeId1))
766  {
767  j = particleIds[0];
768  }
769 
770  if (j!=UINT_MAX)
771  {
772  e.faceVertexId[i] = j;
773  e.faceVertex[i] = getDPMBase()->particleHandler.getObjectById(j);
774  e.face[i] = face0;
775  e.faceInitialArea[i] = face0->getArea();
776  if (i==1)
777  {
778  e.initialSineHalfTheta = e.getSineHalfTheta();
779  }
780  i += 1;
781  }
782  }
783  }
784 }
std::array< unsigned int, 3 > getVertexIds() const
Returns an array containing the ids of the vertex particles.
Definition: MeshTriangle.h:160

References edge_, face_, MeshTriangle::getArea(), getDPMBase(), BaseHandler< T >::getObjectById(), MeshTriangle::getVertexIds(), constants::i, and DPMBase::particleHandler.

Referenced by buildMesh().

◆ loadFromSTL() [1/2]

void Membrane::loadFromSTL ( BaseParticle p0,
std::string  fileName 
)

Load the Membrane geometry from a .stl file.

Parameters
[in]p0BaseParticle that will be used for the vertex particles
[in]fileNameThe name of the file to load.

This function loads the geometry from the STL file and calls the function buildMesh.

487 {
488  loadFromSTL(p0, fileName, 1e-8);
489 }
void loadFromSTL(BaseParticle &p0, std::string fileName)
Load the Membrane geometry from a .stl file.
Definition: Membrane.cc:486

Referenced by MembraneDemo::setUpMembrane(), and MembraneSelfTest::setUpMembrane().

◆ loadFromSTL() [2/2]

void Membrane::loadFromSTL ( BaseParticle p0,
std::string  fileName,
Mdouble  eps 
)

Load the Membrane geometry from a .stl file.

Parameters
[in]p0BaseParticle that will be used for the vertex particles
[in]fileNameThe name of the file to load.
[in]epsThe distance between two points below which they are assumed to be the same.

This function loads the geometry from the STL file and calls the function buildMesh.

500 {
501  std::vector<Vec3D> vertexPositions;
502  std::vector<unsigned int> edgeVertices;
503  std::vector<unsigned int> faceVertices;
504 
505 
506  std::ifstream is;
507  is.open(fileName, std::ios_base::in);
508 
509  if (!is.is_open())
510  logger(ERROR, "Membrane::loadFromSTL Cannot open file %.", fileName);
511 
512  std::string dummy;
513  unsigned int i;
514  Vec3D pos;
515 
516  // solid name
517  is >> dummy >> dummy;
518  // endsolid or facet
519  is >> dummy;
520  while ( dummy.compare("endsolid") )
521  {
522  // normal Vec3D
523  is >> dummy >> pos;
524  // outer loop
525  is >> dummy >> dummy;
526  for (i=0; i<3; i++)
527  {
528  // vertex pos
529  is >> dummy >> pos;
530  faceVertices.push_back(addVertex(vertexPositions, pos, eps));
531  }
532  // endloop
533  is >> dummy;
534  // endfacet
535  is >> dummy;
536  // endsolid or facet
537  is >> dummy;
538  }
539  is.close();
540 
541  // Create the neccesary edges;
542  std::map<std::pair<unsigned int, unsigned int>, int> existingEdges;
543  for (i=0; i<faceVertices.size()/3; i++)
544  {
545  addEdge(edgeVertices, existingEdges, faceVertices[3*i+0], faceVertices[3*i+1]);
546  addEdge(edgeVertices, existingEdges, faceVertices[3*i+1], faceVertices[3*i+2]);
547  addEdge(edgeVertices, existingEdges, faceVertices[3*i+2], faceVertices[3*i+0]);
548  }
549 
550  // logger(INFO, "Got % vertices, % edges % faces", vertexPositions.size(), edgeVertices.size()/2, faceVertices.size()/3);
551  buildMesh(p0, vertexPositions, edgeVertices, faceVertices);
552 }
@ ERROR
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 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 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

References addEdge(), addVertex(), buildMesh(), ERROR, constants::i, and logger.

◆ loadVertexPositions()

void Membrane::loadVertexPositions ( std::istream &  is)

Load the positions of the vertex particles from a stream and apply them to existing particles.

384 {
385  unsigned int nVertices, i;
386  std::string dummy;
387  Vec3D pos;
388  is >> dummy >> nVertices;
389  is >> dummy;
390 
391  #ifdef MERCURYDPM_USE_MPI
392  if (PROCESSOR_ID==0)
393  #endif
394  logger.assert_always(nVertices==vertexParticleId_.size(),
395  "Error in Membrane::loadVertexPositions: The number of saved Particles differs from the number of actual particles");
396  for (i=0; i<nVertices; i++)
397  {
398  is >> pos;
399  if (vertexParticle_[i])
400  {
401  vertexParticle_[i]->setPosition(pos);
402  }
403  }
404 
405  // #ifdef MERCURYDPM_USE_MPI
406  // if (PROCESSOR_ID==0)
407  // #endif
408  // logger(INFO, "Loaded % membrane particles", nVertices);
409 }
#define PROCESSOR_ID
Definition: GeneralDefine.h:63

References constants::i, logger, PROCESSOR_ID, vertexParticle_, and vertexParticleId_.

◆ read()

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

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

Parameters
[in]isThe input stream from which the Membrane is read, usually a restart file.

Implements BaseObject.

90 {
91 
92  logger.assert_debug(DPMBase_, "Error in Membrane::read: DPMBase is not knwon.");
93 
94  std::string dummy;
95  unsigned int i, n, Id;
96 
97  // Reading in the Ids of all the vertex particles
98  is >> dummy >> n;
99  vertexParticle_.reserve(n);
100  for (i=0; i<n; i++)
101  {
102  is >> Id;
103  vertexParticleId_.push_back(Id);
104  // this might retrieve a nullptr
105  vertexParticle_.push_back(getDPMBase()->particleHandler.getObjectById(Id));
106  }
107  is >> dummy >> vertexInitId_;
108 
109  // Reading in the Ids of the faces
110  is >> dummy >> n;
111  face_.reserve(n);
112  MeshTriangle* f;
113  for (i=0; i<n; i++)
114  {
115  is >> Id;
116  f = dynamic_cast<MeshTriangle*>(getDPMBase()->wallHandler.getObjectById(Id));
117  face_.push_back(f);
118  }
119 
120  // Reading in the edge information
121  is >> dummy >> n;
122  edge_.reserve(n);
123  for (i=0; i<n; i++)
124  {
125  Edge e;
126  // the vertex ids
127  is >> dummy >> e.vertexId[0];
128  e.vertex[0] = getDPMBase()->particleHandler.getObjectById(e.vertexId[0]);
129  is >> dummy >> e.vertexId[1];
130  e.vertex[1] = getDPMBase()->particleHandler.getObjectById(e.vertexId[1]);
131 
132  is >> dummy >> e.faceVertexId[0];
133  if (e.faceVertexId[0] != UINT_MAX) {
134  e.faceVertex[0] = getDPMBase()->particleHandler.getObjectById(e.faceVertexId[0]);
135  } else {
136  e.faceVertex[0] = nullptr;
137  }
138  is >> dummy >> e.faceVertexId[1];
139  if (e.faceVertexId[1] != UINT_MAX) {
140 
141  e.faceVertex[1] = getDPMBase()->particleHandler.getObjectById(e.faceVertexId[1]);
142  } else {
143  e.faceVertex[1] = nullptr;
144  }
145  is >> dummy >> Id;
146  if (Id != UINT_MAX) {
147  e.face[0] = dynamic_cast<MeshTriangle*>(getDPMBase()->wallHandler.getObjectById(Id));
148  } else {
149  e.face[0] = nullptr;
150  }
151  is >> dummy >> Id;
152  if (Id != UINT_MAX) {
153  e.face[1] = dynamic_cast<MeshTriangle*>(getDPMBase()->wallHandler.getObjectById(Id));
154  } else {
155  e.face[1] = nullptr;
156  }
157 
158  is >> dummy >> e.faceInitialArea[0];
159  is >> dummy >> e.faceInitialArea[1];
160  for (unsigned int i=0; i<8; i++){
161  is >> dummy >> e.uPre[i];
162  }
163 
164  // initial state
165  is >> dummy >> e.initialState;
166  is >> dummy >> e.initialLength;
167 
168  is >> dummy >> e.initialSineHalfTheta;
169 
170  is >> dummy >> e.effectiveMass;
171 
172  e.checkActive();
173  edge_.push_back(e);
174  }
175 
176  // Reading in some general attributes
177  is >> dummy >> particleRadius_; // Mdouble
178 
179  is >> dummy >> Kn_; // Mdouble
180  is >> dummy >> critDampCoeff_; //Mdouble
181 
182  is >> dummy >> Ke_; // Mdouble
183  is >> dummy >> Kd_; // Mdouble
184 
185  is >> dummy >> thickness_;
186 
187  is >> dummy >> bendingAreaConstant_;
188 
189  is >> dummy >> Id;
191  is >> dummy >> Id;
193 
194  for (auto& e: edge_)
195  {
196  e.checkActive();
197  }
198 }
const unsigned n
Definition: CG3DPackingUnitTest.cpp:32
SpeciesHandler speciesHandler
A handler to that stores the species type i.e. LinearViscoelasticSpecies, etc.
Definition: DPMBase.h:1427
WallHandler wallHandler
An object of the class WallHandler. Contains pointers to all the walls created.
Definition: DPMBase.h:1447
void checkActive()
Check if the triangle is considered active.
Definition: MeshTriangle.cc:581

References bendingAreaConstant_, Membrane::Edge::checkActive(), critDampCoeff_, DPMBase_, edge_, Membrane::Edge::effectiveMass, Membrane::Edge::face, face_, Membrane::Edge::faceInitialArea, Membrane::Edge::faceVertex, Membrane::Edge::faceVertexId, getDPMBase(), BaseHandler< T >::getObjectById(), constants::i, Membrane::Edge::initialLength, Membrane::Edge::initialSineHalfTheta, Membrane::Edge::initialState, Kd_, Ke_, Kn_, logger, membraneParticleSpecies_, membraneSpecies_, n, DPMBase::particleHandler, particleRadius_, DPMBase::speciesHandler, thickness_, Membrane::Edge::uPre, Membrane::Edge::vertex, Membrane::Edge::vertexId, vertexInitId_, vertexParticle_, vertexParticleId_, and DPMBase::wallHandler.

◆ saveAsOFF()

void Membrane::saveAsOFF ( unsigned int  d)

Save the Membrane as a .off file.

Parameters
[in]unsignedint d

This function saves the current geometry of the membrane to a file named membrane.$d.off.

417 {
418  std::stringstream lastName("");
419  lastName << "membrane." << d << ".off";
420 
421  std::ofstream membraneOFF;
422  membraneOFF.open(lastName.str(), std::ios_base::out);
423 
424  membraneOFF << "OFF\n";
425  membraneOFF << "# membrane.off\n";
426 
427  // Number of Vertices Faces and Edges (edges are optional)
428  membraneOFF << vertexParticle_.size() << " " << face_.size() << " 0" << std::endl;
429 
430  // Print the coordinates of all vertices
431  for (auto p: vertexParticle_)
432  {
433  membraneOFF << p->getPosition() << std::endl;
434  }
435 
436  // Print all faces
437  for (auto f: face_)
438  {
439  membraneOFF << "3 ";
440  for (auto id: f->getVertexIds())
441  {
442  membraneOFF << id - vertexInitId_ << " ";
443  }
444  membraneOFF << std::endl;
445  }
446 }

References face_, vertexInitId_, and vertexParticle_.

◆ saveAsSTL()

void Membrane::saveAsSTL ( std::string  fileName)

Save the Membrane as a .stl file.

Parameters
[in]fileNameThe name of the file where the data is written to.

This function writes the current geometry of the Membrane to a STL file.

454 {
455  std::ofstream os;
456  os.open(fileName, std::ios_base::out);
457 
458  if (!os.is_open())
459  logger(ERROR, "Membrane::saveAsSTL Cannot open file %.", fileName);
460 
461  os << "solid Membrane" << std::endl;
462 
463  std::array<Vec3D,3> vertexPositions;
464  for (auto face: face_)
465  {
466  os << " facet normal " << face->getFaceNormal() << std::endl;
467  os << " outer loop" << std::endl;
468  vertexPositions = face->getVertices();
469  for (auto pos: vertexPositions)
470  {
471  os << " vertex " << pos << std::endl;
472  }
473  os << " endloop" << std::endl;
474  os << " endfacet" << std::endl;
475  }
476 
477  os << "endsolid Membrane" << std::endl;
478 }

References ERROR, face_, and logger.

◆ saveVertexPositions()

void Membrane::saveVertexPositions ( std::ostream &  os)

save the current positions of the vertex particles to a stream

335 {
336  #ifdef MERCURYDPM_USE_MPI
337  if (PROCESSOR_ID==0)
338  {
339  os << " nMembraneParticles " << vertexParticleId_.size();
340  os << " membraneParticlePositions ";
341  }
342 
343  Vec3D pos, pos1;
344  Mdouble active;
345  BaseParticle* p0;
346  for (unsigned int i=0; i < vertexParticleId_.size(); i++)
347  {
348  pos = Vec3D(0,0,0);
350  if (p0 && !p0->isMPIParticle())
351  {
352  pos = p0->getPosition();
353  active = getMPISum(1.0);
354  }
355  else
356  {
357  active = getMPISum(0.0);
358  }
359 
360  pos1 = getMPISum(pos);
361  if (PROCESSOR_ID==0 && active > 0.1)
362  {
363  os << pos1 << " ";
364  }
365  }
366  if (PROCESSOR_ID==0)
367  {
368  os << "\n";
369  }
370  #else
371  os << " nMembraneParticles " << vertexParticleId_.size();
372  os << " membraneParticlePositions ";
373  BaseParticle* p0;
374  for (unsigned int i=0; i < vertexParticleId_.size(); i++)
375  {
377  os << p0->getPosition() << " ";
378  }
379  #endif
380 }
Vec3D getMPISum(Vec3D &val)
Definition: MpiDataClass.cc:199
const Vec3D & getPosition() const
Returns the position of this BaseInteractable.
Definition: BaseInteractable.h:218
Definition: BaseParticle.h:54
bool isMPIParticle() const
Indicates if this particle is a ghost in the MPI domain.
Definition: BaseParticle.cc:181

References getDPMBase(), getMPISum(), BaseHandler< T >::getObjectById(), BaseInteractable::getPosition(), constants::i, BaseParticle::isMPIParticle(), DPMBase::particleHandler, PROCESSOR_ID, and vertexParticleId_.

◆ setBendingAreaConstant()

void Membrane::setBendingAreaConstant ( bool  areaConstant)

If set to true, the bending penalty is calulated using the initial triangle area instead of the current one.

324 {
325  bendingAreaConstant_ = areaConstant;
326 }

References bendingAreaConstant_.

◆ setCriticalDampingCoefficient()

void Membrane::setCriticalDampingCoefficient ( Mdouble  coeff)

Set damping coefficient for the distance springs.

304 {
305  critDampCoeff_ = coeff;
306 }

References critDampCoeff_.

Referenced by MembraneDemo::setUpMembrane(), and MembraneSelfTest::setUpMembrane().

◆ setDPMBase()

void Membrane::setDPMBase ( DPMBase dpm)
inline

Set a pointer to DPMBase.

301 { DPMBase_ = dpm; }

References DPMBase_.

Referenced by MembraneDemo::read(), MembraneSelfTest::read(), MembraneDemo::setUpMembrane(), and MembraneSelfTest::setUpMembrane().

◆ setElasticModulusAndThickness()

void Membrane::setElasticModulusAndThickness ( Mdouble  E,
Mdouble  thickness 
)

Set the elastic modulus and thickness of the membrane \deltails The supplied values are used to calculate the spring constant of the membrane.

292 {
293  setThickness(thickness);
294  Kn_ = E * thickness_/(sqrt(3)*(1-1/3.0));
295 }
void setThickness(Mdouble thickness)
Set the thickness of the membrane.
Definition: Membrane.cc:317
double E
Elastic modulus.
Definition: TwenteMeshGluing.cpp:64

References Global_Physical_Variables::E, Kn_, setThickness(), and thickness_.

Referenced by MembraneDemo::setUpMembrane(), and MembraneSelfTest::setUpMembrane().

◆ setKeAndKd()

void Membrane::setKeAndKd ( Mdouble  Ke,
Mdouble  Kd 
)

Set the parameters needed for the bending forces.

312 {
313  Ke_ = Ke;
314  Kd_ = Kd;
315 }

References Kd_, and Ke_.

Referenced by MembraneDemo::setUpMembrane(), and MembraneSelfTest::setUpMembrane().

◆ setKnAndCrittDampCoeff()

void Membrane::setKnAndCrittDampCoeff ( Mdouble  Kn,
Mdouble  critDampCoeff 
)

Set the parameters needed for the stretching forces.

286  {
287  Kn_ = Kn;
288  critDampCoeff_ = critDampCoeff;
289 }

References critDampCoeff_, and Kn_.

◆ setParticleRadius()

void Membrane::setParticleRadius ( Mdouble  radius)

Set the radius of the vertex particles.

329 {
330  logger.assert_debug(radius > 0, "Error in Membrane::setParticleRadius: Radius has to be greater than 0.");
331  particleRadius_ = radius;
332 }

References logger, and particleRadius_.

Referenced by MembraneDemo::setUpMembrane(), and MembraneSelfTest::setUpMembrane().

◆ setSpringConstant()

void Membrane::setSpringConstant ( Mdouble  k)

Set the spring constant of the membrane.

298 {
299  logger.assert_debug(k >= 0, "Error in Membrane::setSpringConstant: The constant has to be greater than or equal to 0.");
300  Kn_ = k;
301 }

References Kn_, and logger.

◆ setThickness()

void Membrane::setThickness ( Mdouble  thickness)

Set the thickness of the membrane.

318 {
319  logger.assert_debug(thickness > 0, "Error in Membrane::setThickness: The thickness has to be greater than 0.");
320  thickness_ = thickness;
321 }

References logger, and thickness_.

Referenced by setElasticModulusAndThickness().

◆ updateEdgeMass()

void Membrane::updateEdgeMass ( )

Set the correct edge mass by taking the mass from the conencted vertices.

Set the correct edge mass by taking the mass from the conencted vertices.

941 {
942  Mdouble invMass;
943  for (auto e: edge_)
944  {
945  invMass = vertexParticle_[e.vertexId[0]-vertexInitId_]->getInvMass()
946  + vertexParticle_[e.vertexId[1]-vertexInitId_]->getInvMass();
947 
948  if (invMass > 0)
949  {
950  e.effectiveMass = 1/invMass;
951  }
952  else
953  {
954  e.effectiveMass = 0;
955  }
956 
957  e.checkActive();
958  }
959 }

References edge_, vertexInitId_, and vertexParticle_.

Referenced by MembraneDemo::fixMembraneEdges().

◆ updateFaceNeighbors()

void Membrane::updateFaceNeighbors ( )
private

Update the faces to have the correct neighbors.

This function assigns the coorect neighbors to the wall elements.

790 {
791  //set neighbours
792  unsigned int i, j, k, l;
793 
794  for (i =0 ; i<face_.size(); i++)
795  {
796  face_[i]->vertexNeighbors = std::vector<std::vector<unsigned int>>();
797  for (j=0; j<3; j++)
798  {
799  face_[i]->neighbor[j] = nullptr;
800  std::vector<unsigned int> a;
801  face_[i]->vertexNeighbors.push_back(a);
802  }
803  }
804  for (i =0 ; i<face_.size(); i++)
805  {
806  MeshTriangle* face0 = face_[i];
807  std::array<unsigned int, 3> particleIds0 = face0->getVertexIds();
808  for (j=i+1 ; j<face_.size(); j++)
809  {
810  MeshTriangle* face1 = face_[j];
811  std::array<unsigned int, 3> particleIds1 = face1->getVertexIds();
812  if (particleIds0[0] == particleIds1[0])
813  {
814  if (particleIds0[1] == particleIds1[2])
815  { //edge 0=2
816  face0->neighbor[0] = &*face1;
817  face1->neighbor[2] = &*face0;
818  }
819  else if (particleIds0[2] == particleIds1[1])
820  { //edge 2=0
821  face0->neighbor[2] = &*face1;
822  face1->neighbor[0] = &*face0;
823  }
824  }
825  else if (particleIds0[0] == particleIds1[1])
826  {
827  if (particleIds0[1] == particleIds1[0])
828  { //edge 0=0
829  face0->neighbor[0] = &*face1;
830  face1->neighbor[0] = &*face0;
831  }
832  else if (particleIds0[2] == particleIds1[2])
833  { //edge 2=1
834  face0->neighbor[2] = &*face1;
835  face1->neighbor[1] = &*face0;
836  }
837  }
838  else if (particleIds0[0] == particleIds1[2])
839  {
840  if (particleIds0[1] == particleIds1[1])
841  { //edge 0=1
842  face0->neighbor[0] = &*face1;
843  face1->neighbor[1] = &*face0;
844  }
845  else if (particleIds0[2] == particleIds1[0])
846  { //edge 2=2
847  face0->neighbor[2] = &*face1;
848  face1->neighbor[2] = &*face0;
849  }
850  }
851  else if (particleIds0[1] == particleIds1[0])
852  {
853  if (particleIds0[2] == particleIds1[2])
854  { //edge 1=2
855  face0->neighbor[1] = &*face1;
856  face1->neighbor[2] = &*face0;
857  }
858  }
859  else if (particleIds0[1] == particleIds1[1])
860  {
861  if (particleIds0[2] == particleIds1[0])
862  { //edge 1=0
863  face0->neighbor[1] = &*face1;
864  face1->neighbor[0] = &*face0;
865  }
866  }
867  else if (particleIds0[1] == particleIds1[2])
868  {
869  if (particleIds0[2] == particleIds1[1])
870  { //edge 1=1
871  face0->neighbor[1] = &*face1;
872  face1->neighbor[1] = &*face0;
873  }
874  }
875 
876  }
877  }
878 
879  for (i =0 ; i<face_.size(); i++)
880  {
881  MeshTriangle* face0 = face_[i];
882  std::array<unsigned int, 3> particleIds0 = face0->getVertexIds();
883  for (j=i+1 ; j<face_.size(); j++)
884  {
885  MeshTriangle* face1 = face_[j];
886  std::array<unsigned int, 3> particleIds1 = face1->getVertexIds();
887  for (k=0; k<3; k++)
888  {
889  for (l=0; l<3; l++)
890  {
891  if( particleIds0[k] == particleIds1[l])
892  {
893  face0->vertexNeighbors[k].push_back(face1->getId());
894  face1->vertexNeighbors[l].push_back(face0->getId());
895  }
896  }
897  }
898  }
899  }
900 }
unsigned int getId() const
Returns the unique identifier of any particular object.
Definition: BaseObject.h:125
std::vector< std::vector< unsigned int > > vertexNeighbors
Definition: MeshTriangle.h:262
std::array< MeshTriangle *, 3 > neighbor
Definition: MeshTriangle.h:257

References face_, BaseObject::getId(), MeshTriangle::getVertexIds(), constants::i, MeshTriangle::neighbor, and MeshTriangle::vertexNeighbors.

Referenced by buildMesh().

◆ write()

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

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

Parameters
[in]osThe output stream where the Membrane must be written to, usually a restart file.

Implements BaseObject.

205 {
206 
207  // Write out the vertex-Particle Ids
208  os << " vertexParticle " << vertexParticleId_.size();
209  for (unsigned int i: vertexParticleId_)
210  {
211  os << " " << i;
212  }
213 
214  os << " vertexInitId_ " << vertexInitId_;
215  os << "\n";
216 
217  // Write out the face-wall Ids
218  os << " face " << face_.size();
219  for (auto f: face_)
220  {
221  os << " " << f->getId();
222  }
223  os << "\n";
224 
225  // Write out edge information
226  os << " edge " << edge_.size();
227  for (auto e: edge_)
228  {
229  // the vertex ids
230  os << " v3 " << e.vertexId[0];
231  os << " v4 " << e.vertexId[1];
232 
233 
234  os << " v1 " << e.faceVertexId[0];
235  os << " v2 " << e.faceVertexId[1];
236 
237  if (e.face[1] != nullptr){
238  os << " f1 " << e.face[0]->getId();
239  } else {
240  os << " f1 " << UINT_MAX;
241  }
242  if (e.face[1] != nullptr){
243  os << " f2 " << e.face[1]->getId();
244  } else {
245  os << " f2 " << UINT_MAX;
246  }
247 
248  os << " area1 " << e.faceInitialArea[0];
249  os << " area2 " << e.faceInitialArea[1];
250  for (unsigned int i=0; i<8; i++){
251  os << " uPre " << e.uPre[i];
252  }
253  // initial state
254  os << " initState " << e.initialState;
255  os << " initDist " << e.initialLength;
256  os << " initialSineHalfTheta " << e.initialSineHalfTheta;
257  os << " effectiveMass " << e.effectiveMass;
258  os << "\n";
259  }
260 
261  // Write out general properties
262  os << " particleRadius " << particleRadius_; // Mdouble
263 
264  os << " Kn " << Kn_; // Mdouble
265  os << " critDampCoeff " << critDampCoeff_; //Mdouble
266 
267  os << " Ke " << Ke_; // Mdouble
268  os << " Kd " << Kd_; // Mdouble
269 
270  os << " thickness " << thickness_;
271 
272  os << " bendingAreaConstant " << bendingAreaConstant_;
273 
274  os << " membraneSpecies " << membraneSpecies_->getId();
275  os << " membraneParticleSpecies " << membraneParticleSpecies_->getId();
276 }

References bendingAreaConstant_, critDampCoeff_, edge_, face_, BaseObject::getId(), constants::i, Kd_, Ke_, Kn_, membraneParticleSpecies_, membraneSpecies_, particleRadius_, thickness_, vertexInitId_, and vertexParticleId_.

Member Data Documentation

◆ bendingAreaConstant_

bool Membrane::bendingAreaConstant_
private

◆ critDampCoeff_

Mdouble Membrane::critDampCoeff_
private

◆ DPMBase_

DPMBase* Membrane::DPMBase_ = nullptr
private

Stores a pointer to DPMBase

Referenced by getDPMBase(), read(), and setDPMBase().

◆ edge_

◆ face_

std::vector<MeshTriangle*> Membrane::face_
private

◆ faceVertices_

std::vector<unsigned int> Membrane::faceVertices_
private

Stores the 3 indices per face which can be used to get the ids from vertexParticleId_ // Todo: remove this

◆ Kd_

Mdouble Membrane::Kd_
private

Stores dissipation constant for the bending

Referenced by computeAdditionalForces(), getKd(), Membrane(), read(), setKeAndKd(), and write().

◆ Ke_

Mdouble Membrane::Ke_
private

Stores spring constant for the bending

Referenced by computeAdditionalForces(), getKe(), Membrane(), read(), setKeAndKd(), and write().

◆ Kn_

Mdouble Membrane::Kn_
private

◆ membraneParticleSpecies_

ParticleSpecies* Membrane::membraneParticleSpecies_ = nullptr
private

Stores a pointer to the membrane particle species

Referenced by adjustVertexParticleSize(), buildMesh(), Membrane(), read(), and write().

◆ membraneSpecies_

ParticleSpecies* Membrane::membraneSpecies_ = nullptr
private

Stores a pointer to the membrane species

Referenced by buildMesh(), Membrane(), read(), and write().

◆ particleRadius_

Mdouble Membrane::particleRadius_
private

Stores the radius of membrane particles

Referenced by buildMesh(), getParticleRadius(), Membrane(), read(), setParticleRadius(), and write().

◆ thickness_

Mdouble Membrane::thickness_
private

Stores thickness of the Membrane

Referenced by buildMesh(), getThickness(), Membrane(), read(), setElasticModulusAndThickness(), setThickness(), and write().

◆ vertexInitId_

unsigned int Membrane::vertexInitId_
private

Sores the id of the first particle created for the membrane. This value helps for the eventual MPI parallelization.

Referenced by adjustVertexParticleSize(), createVertexParticles(), handleParticleAddition(), handleParticleRemoval(), Membrane(), read(), saveAsOFF(), updateEdgeMass(), and write().

◆ vertexParticle_

std::vector<BaseParticle*> Membrane::vertexParticle_
private

◆ vertexParticleId_

std::vector<unsigned int> Membrane::vertexParticleId_
private

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