NetDEM v1.0
Loading...
Searching...
No Matches
90_packing_demo.cpp

This is an example of how to use the netdem library.

#include "data_dumper.hpp"
#include "gen_pack.hpp"
#include "gravity.hpp"
#include "shape_sphere.hpp"
#include "simulation.hpp"
#include <filesystem>
#include <iostream>
#include <unordered_map>
using namespace netdem;
using namespace std;
int PackingDemo(int argc, char **argv) {
Simulation *sim = new Simulation();
sim->domain_manager.SetBound(-0.6, -0.6, -0.6, 0.6, 0.6, 0.6);
sim->domain_manager.SetCellSpacing(0.3, 0.3, 0.3);
ContactSolverSettings::SolverType::sdf;
LinearSpring contact_model = LinearSpring(2.0e3, 1.0e4, 0.3, 0.0);
sim->scene.InsertContactModel(&contact_model);
sim->scene.SetCollisionModel(0, 0, contact_model.label);
int shape_case;
if (argc == 2) {
shape_case = -1;
} else {
shape_case = atof(argv[2]);
}
string root_dir;
switch (shape_case) {
case 0: {
// convex
PolySuperEllipsoid ellipsoid(0.5, 1, 1.5, 0.5, 1, 1.5, 1.5, 1.5);
ellipsoid.SetSize(0.1);
sim->scene.InsertShape(&ellipsoid);
root_dir = "tmp/examples/sdf_dem/packing/ellipsoid/";
filesystem::create_directories(root_dir);
} break;
case 1: {
// convex
PolySuperQuadrics quadrics(0.5, 1.0, 2.5, 0.5, 1.7, 0.5, 0.5, 1.0, 1.5, 1.0,
1.2, 0.8);
quadrics.SetSize(0.1);
sim->scene.InsertShape(&quadrics);
root_dir = "tmp/examples/sdf_dem/packing/quadrics/";
filesystem::create_directories(root_dir);
} break;
case 10: {
// non-convex
PolySuperEllipsoid ellipsoid(0.5, 1, 1.5, 0.5, 1, 1.5, 2.5, 2.5);
ellipsoid.SetSize(0.1);
sim->scene.InsertShape(&ellipsoid);
root_dir = "tmp/examples/sdf_dem/packing/ellipsoid_non_convex/";
filesystem::create_directories(root_dir);
} break;
case 11: {
// non-convex
PolySuperQuadrics quadrics(0.5, 1.0, 2.5, 0.5, 1.7, 0.5, 0.5, 1.0, 1.5, 2.0,
2.2, 1.8);
quadrics.SetSize(0.1);
sim->scene.InsertShape(&quadrics);
root_dir = "tmp/examples/sdf_dem/packing/quadrics_non_convex/";
filesystem::create_directories(root_dir);
} break;
case 2: {
sh.InitFromSTL("data/particle_template.stl");
sh.SetSize(0.1);
sim->scene.InsertShape(&sh);
root_dir = "tmp/examples/sdf_dem/packing/spherical_harmonics/";
filesystem::create_directories(root_dir);
} break;
case 3: {
TriMesh tri_mesh;
tri_mesh.InitFromSTL("data/copyleft/bolt_nut/bolt.stl");
tri_mesh.AlignAxes();
tri_mesh.SetSize(0.1);
sim->scene.InsertShape(&tri_mesh);
tri_mesh.InitFromSTL("data/copyleft/bolt_nut/nut.stl");
tri_mesh.AlignAxes();
tri_mesh.SetSize(0.1);
sim->scene.InsertShape(&tri_mesh);
root_dir = "tmp/examples/sdf_dem/packing/nut_bolt/";
filesystem::create_directories(root_dir);
} break;
case 4: {
LevelSet level_set;
level_set.InitFromSTL("data/particle_template.stl");
level_set.AlignAxes();
level_set.SetSize(0.1);
sim->scene.InsertShape(&level_set);
root_dir = "tmp/examples/sdf_dem/packing/level_set/";
filesystem::create_directories(root_dir);
} break;
case 5: {
PolySuperEllipsoid ellipsoid(0.5, 1, 1.5, 0.5, 1, 1.5, 1.5, 1.5);
ellipsoid.SetSize(0.1);
sim->scene.InsertShape(&ellipsoid);
PolySuperQuadrics quadrics(0.5, 1.0, 2.5, 0.5, 1.7, 0.5, 0.5, 1.0, 1.5, 1.0,
1.2, 0.8);
quadrics.SetSize(0.1);
sim->scene.InsertShape(&quadrics);
sh.InitFromSTL("data/particle_template.stl");
sh.SetSize(0.1);
sim->scene.InsertShape(&sh);
TriMesh tri_mesh;
tri_mesh.InitFromSTL("data/copyleft/bolt_nut/bolt.stl");
tri_mesh.AlignAxes();
tri_mesh.SetSize(0.1);
sim->scene.InsertShape(&tri_mesh);
tri_mesh.InitFromSTL("data/copyleft/bolt_nut/nut.stl");
tri_mesh.AlignAxes();
tri_mesh.SetSize(0.1);
sim->scene.InsertShape(&tri_mesh);
LevelSet level_set;
level_set.InitFromSTL("data/particle_template.stl");
level_set.AlignAxes();
level_set.SetSize(0.1);
sim->scene.InsertShape(&level_set);
root_dir = "tmp/examples/sdf_dem/packing/all_shapes/";
filesystem::create_directories(root_dir);
} break;
case 15: {
PolySuperEllipsoid ellipsoid(0.5, 1, 1.5, 0.5, 1, 1.5, 2.5, 2.5);
ellipsoid.SetSize(0.1);
sim->scene.InsertShape(&ellipsoid);
PolySuperQuadrics quadrics(0.5, 1.0, 2.5, 0.5, 1.7, 0.5, 0.5, 1.0, 1.5, 2.0,
2.2, 1.8);
quadrics.SetSize(0.1);
sim->scene.InsertShape(&quadrics);
sh.InitFromSTL("data/particle_template.stl");
sh.SetSize(0.1);
sim->scene.InsertShape(&sh);
TriMesh tri_mesh;
tri_mesh.InitFromSTL("data/copyleft/bolt_nut/bolt.stl");
tri_mesh.AlignAxes();
tri_mesh.SetSize(0.1);
sim->scene.InsertShape(&tri_mesh);
tri_mesh.InitFromSTL("data/copyleft/bolt_nut/nut.stl");
tri_mesh.AlignAxes();
tri_mesh.SetSize(0.1);
sim->scene.InsertShape(&tri_mesh);
LevelSet level_set;
level_set.InitFromSTL("data/particle_template.stl");
level_set.AlignAxes();
level_set.SetSize(0.1);
sim->scene.InsertShape(&level_set);
root_dir = "tmp/examples/sdf_dem/packing/all_shapes/";
filesystem::create_directories(root_dir);
} break;
default: {
TriMesh tri_mesh;
tri_mesh.InitFromSTL("data/particle_template.stl");
tri_mesh.Decimate(200);
tri_mesh.AlignAxes();
tri_mesh.SetSize(0.1);
sim->scene.InsertShape(&tri_mesh);
root_dir = "tmp/examples/sdf_dem/packing/trimesh/";
filesystem::create_directories(root_dir);
} break;
}
VecXT<Particle> particle_list = PackGenerator::GetGridPack(
1, 1, 0.2, 0, 0, 0.4, 5, 5, 1, sim->scene.GetShapes());
for (auto &p : particle_list) {
p.SetVelocity(0, 0, -2.0);
p.damp_numerical = 0.7;
}
WallBoxPlane wall_box(1, 1, 1, 0, 0, 0);
wall_box.ImportToScene(&(sim->scene));
Gravity grav;
grav.Init(sim);
sim->modifier_manager.Insert(&grav);
DataDumper data_dumper;
data_dumper.Init(sim);
data_dumper.SetRootPath(root_dir);
data_dumper.SetSaveByCycles(100);
data_dumper.SaveShapeInfoAsSTL();
data_dumper.dump_wall_info = true;
data_dumper.dump_contact_info = true;
sim->modifier_manager.Insert(&data_dumper);
sim->modifier_manager.Enable(data_dumper.label);
for (int i = 0; i < 20; i++) {
sim->scene.InsertParticle(particle_list);
sim->Run(0.1);
}
sim->Run(3.0);
delete sim;
return 0;
}
std::string label
Definition contact_model.hpp:34
ContactSolverSettings settings
Definition contact_solver_factory.hpp:64
SolverType solver_type
Definition contact_solver_factory.hpp:27
int sdf_potential_type
Definition contact_solver_factory.hpp:43
ContactSolverFactory contact_solver_factory
Definition dem_solver.hpp:50
A class used to dump particle data into vtk files. This is a post-modifier, which will be executed at...
Definition data_dumper.hpp:12
void Init(Simulation *sim) override
Initializes the DataDumper instance.
Definition data_dumper.cpp:24
bool dump_wall_info
A flag that determines whether to dump wall information.
Definition data_dumper.hpp:21
bool dump_contact_info
A flag that determines whether to dump contact information.
Definition data_dumper.hpp:25
void SaveShapeInfoAsSTL()
Saves shape information as an STL file.
Definition data_dumper.cpp:2245
void SetSaveByCycles(double interval)
Sets the interval for saving data by cycles.
Definition data_dumper.cpp:39
void SetRootPath(std::string const &root_path)
Sets the root directory path for the output file.
Definition data_dumper.cpp:26
void SetBound(double bmin_x, double bmin_y, double bmin_z, double bmax_x, double bmax_y, double bmax_z)
Sets the lower and upper bounds of the domain.
Definition domain_manager.cpp:83
void SetCellSpacing(double s_x, double s_y, double s_z)
Sets the spacing between cells in each dimension.
Definition domain_manager.cpp:105
A class used to apply gravity to particles in a DEM simulation.
Definition gravity.hpp:10
void Init(Simulation *sim) override
Initializes the Gravity instance.
Definition gravity.cpp:14
A class for representing a level set function as a Shape object.
Definition shape_level_set.hpp:25
void InitFromSTL(std::string const &file, int mesh_res=40)
Initialize the LevelSet object from an STL file.
void SetSize(double d) override
Set the size of the LevelSet object.
Definition shape_level_set.cpp:218
void AlignAxes()
Align the axes of the LevelSet object.
Definition shape_level_set.cpp:135
Contact model that uses linear spring elements to evaluate contact forces and moments.
Definition model_linear_spring.hpp:16
std::string label
Definition modifier.hpp:20
Modifier * Insert(Modifier *e)
Inserts new modifier into the simulation.
Definition modifier_manager.cpp:17
void Enable(std::string const &label)
Enables a modifier in the simulation.
A class representing a poly superellipsoid with two different orders and three different axes.
Definition shape_poly_super_ellipsoid.hpp:27
void SetSize(double d) override
Set the size of the PolySuperEllipsoid object.
Definition shape_poly_super_ellipsoid.cpp:66
A class representing a poly superquadric with three different axes and orders.
Definition shape_poly_super_quadrics.hpp:30
void SetSize(double d) override
Set the size of the PolySuperQuadrics object.
Definition shape_poly_super_quadrics.cpp:72
ContactModel * InsertContactModel(const ContactModel *const cm_ptr)
Insert a contact model into this scene.
Definition scene.cpp:290
void SetNumberOfMaterials(int num)
Set the number of materials in this scene and initialize the contact lookup table accordingly.
Definition scene.cpp:346
Shape * InsertShape(const Shape *const s_ptr)
Insert a single shape into this scene.
Definition scene.cpp:19
void SetCollisionModel(int mat_type_1, int mat_type_2, ContactModel *const cnt_model)
Set the collision model between two materials.
Definition scene.cpp:397
Particle * InsertParticle(const Particle *const p_ptr)
Insert a single particle into this scene.
Definition scene.cpp:76
VecXT< Shape * > GetShapes()
Return a vector of pointers to all shapes in this scene.
Definition scene.cpp:316
Class for managing a DEM simulation.
Definition simulation.hpp:21
DEMSolver dem_solver
Implements DEM algorithms to solve the scene.
Definition simulation.hpp:47
void Run(double time)
Runs the simulation for a specified period of time.
Definition simulation.cpp:28
DomainManager domain_manager
Manager for domain and sub-domain calculations.
Definition simulation.hpp:31
ModifierManager modifier_manager
Manages add-on features (i.e., customized evaluations not hard-coded in the DEM calculation cycle).
Definition simulation.hpp:53
Scene scene
Contains and manages basic DEM objects (e.g., shapes, particles, walls) for a DEM simulation.
Definition simulation.hpp:42
A class representing a spherical harmonics object.
Definition shape_spherical_harmonics.hpp:24
void InitFromSTL(std::string const &file)
Initialize the SphericalHarmonics object from an STL file.
void SetSize(double d) override
Set the size of the SphericalHarmonics object.
Definition shape_spherical_harmonics.cpp:175
A class representing a triangular mesh in 3D space.
Definition shape_trimesh.hpp:23
void InitFromSTL(std::string const &file)
Initialize the TriMesh object from an STL file.
void Decimate(int num_nodes)
Decimate the TriMesh object.
Definition shape_trimesh.cpp:121
void AlignAxes()
Align the axes of the TriMesh object.
Definition shape_trimesh.cpp:107
void SetSize(double d) override
Set the size of the TriMesh object.
Definition shape_trimesh.cpp:207
A class for generating a box of six walls.
Definition gen_wall_box_plane.hpp:23
void ImportToScene(Scene *scene)
Imports the shapes and walls generated by the wall box into a scene.
Definition gen_wall_box_plane.hpp:102
Definition bond_entry.hpp:7
std::vector< T > VecXT
Definition utils_macros.hpp:31