NuTo
Numerics Tool
ParticleHandler.h
Go to the documentation of this file.
1 /*
2  * ParticleHandler.h
3  *
4  * Created on: 10 Mar 2014
5  * Author: ttitsche
6  */
7 
8 #pragma once
9 
10 
11 #include <vector>
12 #include <random>
13 #include <Eigen/Core>
14 
15 namespace NuTo
16 {
17 
18 class CollidableParticleBase;
19 class CollidableParticleSphere;
20 class VisualizeUnstructuredGrid;
21 class Specimen;
22 typedef std::vector<CollidableParticleSphere*> ParticleContainer;
23 
26 {
27 public:
29  ParticleHandler(int rNumParticles, Eigen::MatrixXd rParticleBoundingBox, double rVelocityRange, double rGrowthRate,
30  int rSeed = 0);
31 
33  ParticleHandler(Eigen::MatrixXd rSpheres, double rVelocityRange, double rRelativeGrowthRate,
34  double rAbsoluteGrowthRate, int rSeed = 0);
35 
37  ParticleHandler(const std::string& rFileName, double rVelocityRange, double rRelativeGrowthRate,
38  double rAbsoluteGrowthRate, int rSeed = 0);
39 
40 
43 
45  void Sync(const double rTime);
46 
48  double GetKineticEnergy() const;
49 
51  double GetVolume() const;
52 
58  void ExportParticlesToVTU3D(std::string rOutputDirectory, int rTimeStep, double rGlobalTime, bool rFinal) const;
59 
62  void ExportParticlesToVTU2D(std::string rOutputFile, double rZCoord) const;
63 
68  void ExportParticlesToGmsh3D(std::string rOutputFile, Specimen& rSpecimen, double rMeshSize) const;
69 
76  void ExportParticlesToGmsh2D(std::string rOutputFile, Specimen& rSpecimen, double rMeshSize, double rZCoord,
77  double rMinRadius) const;
78 
79 
81  void ResetVelocities();
82 
84  Eigen::MatrixXd GetParticles(bool rInitialRadius = false) const;
85 
90  Eigen::MatrixXd GetParticles2D(double rZCoord, double rMinRadius) const;
91 
93  void ExportParticlesToFile(const std::string& rExportFileName, bool rInitialRadius) const;
94 
96  CollidableParticleSphere* GetParticle(const int rIndex) const;
97 
99  int GetNumParticles() const;
100 
102  double GetAbsoluteMininimalDistance(Specimen& rSpecimen);
103 
105  void SetVisualizationFileName(const std::string& rVisualizationFileName);
106 
108  Eigen::VectorXi GetSubBoxDivisions(Specimen& rSpecimen, const int rParticlesPerBox);
109 
110 private:
111  ParticleContainer mParticles;
112 
116  Eigen::VectorXd GetRandomVector(double rStart, double rEnd);
117 
120  Eigen::VectorXd GetRandomVector(const Eigen::MatrixXd& rBounds);
121 
122  void CreateParticlesFromMatrix(const Eigen::MatrixXd& rSpheres, double rVelocityRange, double rRelativeGrowthRate,
123  double rAbsoluteGrowthRate);
124 
125  std::string mVisualizationFileName;
126 
127  std::mt19937 mRNG;
128 };
129 
130 
131 } /* namespace NuTo */
class for Specimen
Definition: Specimen.h:16
int GetNumParticles() const
getter for the particle list size
Definition: ParticleHandler.cpp:181
void SetVisualizationFileName(const std::string &rVisualizationFileName)
optional: change the file name, default: "spheres_"
Definition: ParticleHandler.cpp:341
Eigen::MatrixXd GetParticles2D(double rZCoord, double rMinRadius) const
cut spheres at a given z-coordinate to create circles (in 2D)
Definition: ParticleHandler.cpp:362
class for spherical collidables
Definition: CollidableParticleSphere.h:28
void ExportParticlesToVTU3D(std::string rOutputDirectory, int rTimeStep, double rGlobalTime, bool rFinal) const
writes a sphere visualization file
Definition: ParticleHandler.cpp:83
handles the particle list
Definition: ParticleHandler.h:25
void Sync(const double rTime)
updates all particles to the same time global time rTime
Definition: ParticleHandler.cpp:154
CollidableParticleSphere * GetParticle(const int rIndex) const
get a single particle from the particle list
Definition: ParticleHandler.cpp:176
void ExportParticlesToGmsh3D(std::string rOutputFile, Specimen &rSpecimen, double rMeshSize) const
writes a gmsh .geo file 3D
Definition: ParticleHandler.cpp:395
void ResetVelocities()
resets all velocities
Definition: ParticleHandler.cpp:323
double GetVolume() const
getter for the volume of all particles
Definition: ParticleHandler.cpp:168
Eigen::VectorXi GetSubBoxDivisions(Specimen &rSpecimen, const int rParticlesPerBox)
calculates approximate sub box length, based on box size and the number of particles per sub box ...
Definition: ParticleHandler.cpp:346
void ExportParticlesToVTU2D(std::string rOutputFile, double rZCoord) const
writes a sphere visualization file
Definition: ParticleHandler.cpp:133
double GetAbsoluteMininimalDistance(Specimen &rSpecimen)
calculates the minimal distance between all particles using sub boxes
Definition: ParticleHandler.cpp:186
Definition: Exception.h:6
void ExportParticlesToGmsh2D(std::string rOutputFile, Specimen &rSpecimen, double rMeshSize, double rZCoord, double rMinRadius) const
writes a gmsh .geo file 3D
Definition: ParticleHandler.cpp:421
~ParticleHandler()
destructor, deletes all particles
Definition: ParticleHandler.cpp:60
void ExportParticlesToFile(const std::string &rExportFileName, bool rInitialRadius) const
exports the particle list to a file
Definition: ParticleHandler.cpp:77
double GetKineticEnergy() const
getter for the kinetic energy of all particles
Definition: ParticleHandler.cpp:160
std::vector< CollidableParticleSphere * > ParticleContainer
Definition: ParticleHandler.h:21
Eigen::MatrixXd GetParticles(bool rInitialRadius=false) const
converts the particle list to a Nx4-matrix
Definition: ParticleHandler.cpp:67
ParticleHandler(int rNumParticles, Eigen::MatrixXd rParticleBoundingBox, double rVelocityRange, double rGrowthRate, int rSeed=0)
constructor, builds rNumParticles equal particles
Definition: ParticleHandler.cpp:24