A small application for converting STEP files (demo_CAS_converter.cpp)
Tutorial that teaches how to use the CASCADE module to import, convert and visualize STEP models.
// =============================================================================
// PROJECT CHRONO - http://projectchrono.org
//
// Copyright (c) 2014 projectchrono.org
// All rights reserved.
//
// Use of this source code is governed by a BSD-style license that can be found
// in the LICENSE file at the top level of the distribution and at
// http://projectchrono.org/license-chrono.txt.
//
// =============================================================================
// Authors: Alessandro Tasora, Dario Fusai
// =============================================================================
//
// A small interactive editor to test the convex decomposition settings
// and the STEP conversion features of OpenCASCADE library.
//
// =============================================================================
#include "chrono/core/ChRealtimeStep.h"
#include "chrono/core/ChRandom.h"
#include "chrono/collision/ChConvexDecomposition.h"
#include "chrono/physics/ChSystemNSC.h"
#include "chrono_irrlicht/ChVisualSystemIrrlicht.h"
#include "chrono_cascade/ChCascadeDoc.h"
#include "chrono_cascade/ChCascadeMeshTools.h"
#include "chrono_cascade/ChCascadeIrrMeshTools.h"
// Use the namespace of Chrono
using namespace chrono;
using namespace chrono::irrlicht;
using namespace chrono::cascade;
// Use the main namespace of Irrlicht
using namespace irr;
// Utility global variables
std::shared_ptr<ChTriangleMeshConnected> model_mesh = nullptr;
std::shared_ptr<ChBody> original_body = nullptr;
std::shared_ptr<ChBody> decomposition_body = nullptr;
ChConvexDecompositionHACDv2 decompositionHACDv2;
int hacd_maxhullcount;
int hacd_maxhullmerge;
int hacd_maxhullvertexes;
double hacd_concavity;
double hacd_smallclusterthreshold;
double hacd_fusetolerance;
// Load a triangle mesh using Irrlicht importer
model_mesh = ChTriangleMeshConnected::CreateFromWavefrontFile(filename);
auto model_visshape = chrono_types::make_shared<ChVisualShapeTriangleMesh>(model_mesh, false);
original_body->AddVisualShape(model_visshape);
vis->BindItem(original_body);
}
// Load STEP model using OpenCASCADE
std::cout << "Loading the .STEP model..." << std::endl;
ChCascadeDoc cas_doc;
bool load_ok = cas_doc.LoadSTEP(filename);
if (load_ok) {
// Print hierarchy on screen
cas_doc.Dump(std::cout);
// Find all shapes and get as a single compound
TopoDS_Shape shape;
cas_doc.GetRootShape(shape);
if (!shape.IsNull()) {
model_mesh = chrono_types::make_shared<ChTriangleMeshConnected>();
ChCascadeMeshTools::FillTriangleMeshFromCascade(*model_mesh, shape, ChCascadeTriangulate());
auto model_visshape = chrono_types::make_shared<ChVisualShapeTriangleMesh>(model_mesh, false);
original_body->AddVisualShape(model_visshape);
vis->BindItem(original_body);
std::cout << " ...done" << std::endl;
}
} else {
std::cerr << "Error: unable to load .STEP file\n";
}
}
// Perform convex decomposition using HACDv2
void DecomposeModel(ChVisualSystemIrrlicht* vis) {
if (!model_mesh)
return;
// Perform the convex decomposition using the desired parameters.
decompositionHACDv2.Reset();
decompositionHACDv2.AddTriangleMesh(*model_mesh);
decompositionHACDv2.SetParameters(hacd_maxhullcount, hacd_maxhullmerge, hacd_maxhullvertexes, (float)hacd_concavity,
(float)hacd_smallclusterthreshold, (float)hacd_fusetolerance);
decompositionHACDv2.ComputeConvexDecomposition();
std::cout << "Convex decomposition completed\n";
// Visualize the resulting convex decomposition by creating many colored meshes, each per convex hull.
decomposition_body = chrono_types::make_shared<ChBody>();
auto chmesh_hull = chrono_types::make_shared<ChTriangleMeshSoup>();
decompositionHACDv2.GetConvexHullResult(j, *chmesh_hull);
auto trimesh_connected = chrono_types::make_shared<ChTriangleMeshConnected>();
for (const auto& tri : chmesh_hull->GetTriangles()) {
trimesh_connected->AddTriangle(tri);
}
auto decomposition_visshape = chrono_types::make_shared<ChVisualShapeTriangleMesh>(trimesh_connected, false);
ChColor col(0.1 + 0.5 * ChRandom::Get(), 0.3 + 0.2 * ChRandom::Get(), 0.5 * ChRandom::Get());
decomposition_visshape->SetColor(col);
decomposition_body->AddVisualShape(decomposition_visshape);
}
vis->BindItem(decomposition_body);
vis->UnbindItem(original_body);
}
// Save the convex decomposition as Wavefront .obj to a file
try {
std::ofstream decomposed_objfile(filename);
decompositionHACDv2.WriteConvexHullsAsWavefrontObj(decomposed_objfile);
std::cout << "Saved Wavefront file: " << filename << "\n";
} catch (...) {
vis->GetGUIEnvironment()->addMessageBox(L"Save file error", L"Impossible to write into file.");
}
}
// Save the convex decomposition as chulls list to an .obj file
try {
std::ofstream decomposed_objfile(filename);
decompositionHACDv2.WriteConvexHullsAsChullsFile(decomposed_objfile);
std::cout << "Saved chulls file: " << filename << "\n";
} catch (...) {
vis->GetGUIEnvironment()->addMessageBox(L"Save file error", L"Impossible to write into file.");
}
}
// Define an Irrlicht event receiver class which will be used to manage input from the GUI
class MyEventReceiver : public IEventReceiver {
public:
MyEventReceiver(ChVisualSystemIrrlicht* vsys) : vis(vsys) {
vis->GetDevice()->setEventReceiver(this);
// create menu
menu = vis->GetGUIEnvironment()->addMenu();
menu->addItem(L"File", -1, true, true);
menu->addItem(L"View", -1, true, true);
gui::IGUIContextMenu* submenu;
submenu = menu->getSubMenu(0);
submenu->addItem(L"Load OBJ mesh...", 90);
submenu->addItem(L"Load STEP model...", 95);
submenu->addItem(L"Save convex hulls (.obj)", 91);
submenu->addItem(L"Save convex hulls (.chulls)", 96);
submenu->addSeparator();
submenu->addItem(L"Quit", 92);
submenu = menu->getSubMenu(1);
submenu->addItem(L"View model", 93);
submenu->addItem(L"View decomposition", 94);
text_algo_type =
vis->GetGUIEnvironment()->addStaticText(L"HACDv2 algorithm", core::rect<s32>(510, 35, 650, 50), false);
// ..add a GUI
edit_hacd_maxhullcount = vis->GetGUIEnvironment()->addEditBox(
irr::core::stringw((int)hacd_maxhullcount).c_str(), core::rect<s32>(510, 60, 650, 75), true, 0, 121);
text_hacd_maxhullcount =
vis->GetGUIEnvironment()->addStaticText(L"Max. hull count ", core::rect<s32>(650, 60, 750, 75), false);
// ..add a GUI
edit_hacd_maxhullmerge = vis->GetGUIEnvironment()->addEditBox(
irr::core::stringw((int)hacd_maxhullmerge).c_str(), core::rect<s32>(510, 85, 650, 100), true, 0, 122);
text_hacd_maxhullmerge =
vis->GetGUIEnvironment()->addStaticText(L"Max. hull merge ", core::rect<s32>(650, 85, 750, 100), false);
// ..add a GUI
edit_hacd_maxhullvertexes = vis->GetGUIEnvironment()->addEditBox(
irr::core::stringw((int)hacd_maxhullvertexes).c_str(), core::rect<s32>(510, 110, 650, 125), true, 0, 123);
text_hacd_maxhullvertexes = vis->GetGUIEnvironment()->addStaticText(L"Max. vertices per hull",
core::rect<s32>(650, 110, 750, 125), false);
// ..add a GUI
edit_hacd_concavity = vis->GetGUIEnvironment()->addEditBox(irr::core::stringw(hacd_concavity).c_str(),
core::rect<s32>(510, 135, 650, 150), true, 0, 124);
text_hacd_concavity = vis->GetGUIEnvironment()->addStaticText(L"Max. concavity (0..1)",
core::rect<s32>(650, 135, 750, 150), false);
// ..add a GUI
edit_hacd_smallclusterthreshold = vis->GetGUIEnvironment()->addEditBox(
irr::core::stringw(hacd_smallclusterthreshold).c_str(), core::rect<s32>(510, 160, 650, 175), true, 0, 125);
text_hacd_smallclusterthreshold = vis->GetGUIEnvironment()->addStaticText(
L"Small cluster threshold", core::rect<s32>(650, 160, 750, 175), false);
// ..add a GUI
edit_hacd_fusetolerance = vis->GetGUIEnvironment()->addEditBox(
irr::core::stringw(hacd_fusetolerance).c_str(), core::rect<s32>(510, 185, 650, 200), true, 0, 126);
text_hacd_fusetolerance = vis->GetGUIEnvironment()->addStaticText(L"Vertex fuse tolerance",
core::rect<s32>(650, 185, 750, 200), false);
// .. add buttons..
button_decompose = vis->GetGUIEnvironment()->addButton(core::rect<s32>(510, 210, 650, 225), 0, 106,
L"Decompose", L"Perform convex decomposition");
text_hacd_maxhullcount->setVisible(true);
edit_hacd_maxhullcount->setVisible(true);
text_hacd_maxhullmerge->setVisible(true);
edit_hacd_maxhullmerge->setVisible(true);
text_hacd_maxhullvertexes->setVisible(true);
edit_hacd_maxhullvertexes->setVisible(true);
text_hacd_concavity->setVisible(true);
edit_hacd_concavity->setVisible(true);
text_hacd_smallclusterthreshold->setVisible(true);
edit_hacd_smallclusterthreshold->setVisible(true);
text_hacd_fusetolerance->setVisible(true);
edit_hacd_fusetolerance->setVisible(true);
}
bool OnEvent(const SEvent& event) {
// check if user moved the sliders with mouse..
if (event.EventType == EET_GUI_EVENT) {
s32 id = event.GUIEvent.Caller->getID();
auto env = vis->GetGUIEnvironment();
switch (event.GUIEvent.EventType) {
case gui::EGET_MENU_ITEM_SELECTED: {
// a menu item was clicked
menu = (gui::IGUIContextMenu*)event.GUIEvent.Caller;
id = menu->getItemCommandId(menu->getSelectedItem());
switch (id) {
case 90:
env->addFileOpenDialog(L"Load a mesh file, in .OBJ/.X/.MAX format", true, 0, 80);
break;
case 95:
env->addFileOpenDialog(L"Load a 3D CAD model, in STEP format", true, 0, 85);
break;
case 91:
env->addFileOpenDialog(L"Save decomposed convex hulls in .obj 3D format", true, 0, 81);
break;
case 96:
env->addFileOpenDialog(L"Save decomposed convex hulls in .chulls 3D format", true, 0, 86);
break;
case 92: // File -> Quit
vis->GetDevice()->closeDevice();
break;
case 93: // view model
vis->BindItem(original_body);
vis->UnbindItem(decomposition_body);
break;
case 94: // view decomposition
vis->BindItem(decomposition_body);
vis->UnbindItem(original_body);
break;
}
break;
}
case gui::EGET_FILE_SELECTED: {
// load the model file, selected in the file open dialog
gui::IGUIFileOpenDialog* dialog = (gui::IGUIFileOpenDialog*)event.GUIEvent.Caller;
switch (id) {
case 80:
LoadModel(vis, core::stringc(dialog->getFileName()).c_str());
break;
case 85:
LoadStepModel(vis, core::stringc(dialog->getFileName()).c_str());
break;
case 81:
SaveHullsWavefront(vis, core::stringc(dialog->getFileName()).c_str());
break;
case 86:
SaveHullsChulls(vis, core::stringc(dialog->getFileName()).c_str());
break;
}
} break;
case gui::EGET_EDITBOX_ENTER: {
// load the model file, selected in the file open dialog
gui::IGUIEditBox* medit = (gui::IGUIEditBox*)event.GUIEvent.Caller;
switch (id) {
case 122:
hacd_maxhullmerge = std::atoi(irr::core::stringc(medit->getText()).c_str());
medit->setText(irr::core::stringw(hacd_maxhullmerge).c_str());
break;
case 123:
hacd_maxhullvertexes = std::atoi(irr::core::stringc(medit->getText()).c_str());
medit->setText(irr::core::stringw(hacd_maxhullvertexes).c_str());
break;
case 124:
hacd_concavity = std::atof(irr::core::stringc(medit->getText()).c_str());
medit->setText(irr::core::stringw(hacd_concavity).c_str());
break;
case 125:
hacd_smallclusterthreshold = std::atof(irr::core::stringc(medit->getText()).c_str());
medit->setText(irr::core::stringw(hacd_smallclusterthreshold).c_str());
break;
case 126:
hacd_fusetolerance = std::atof(irr::core::stringc(medit->getText()).c_str());
medit->setText(irr::core::stringw(hacd_fusetolerance).c_str());
break;
}
} break;
case gui::EGET_BUTTON_CLICKED: {
switch (id) {
case 106:
DecomposeModel(vis);
return true;
default:
return false;
}
break;
}
default:
break;
}
}
return false;
}
private:
ChVisualSystemIrrlicht* vis;
gui::IGUIContextMenu* menu;
gui::IGUIButton* button_decompose;
gui::IGUIStaticText* text_algo_type;
gui::IGUIStaticText* text_hacd_maxhullcount;
gui::IGUIEditBox* edit_hacd_maxhullcount;
gui::IGUIStaticText* text_hacd_maxhullmerge;
gui::IGUIEditBox* edit_hacd_maxhullmerge;
gui::IGUIStaticText* text_hacd_maxhullvertexes;
gui::IGUIEditBox* edit_hacd_maxhullvertexes;
gui::IGUIStaticText* text_hacd_concavity;
gui::IGUIEditBox* edit_hacd_concavity;
gui::IGUIStaticText* text_hacd_smallclusterthreshold;
gui::IGUIEditBox* edit_hacd_smallclusterthreshold;
gui::IGUIStaticText* text_hacd_fusetolerance;
gui::IGUIEditBox* edit_hacd_fusetolerance;
};
// Main program
int main(int argc, char* argv[]) {
std::cout << "Copyright (c) 2017 projectchrono.org\nChrono version: " << CHRONO_VERSION << std::endl;
// Create Chrono physical system
ChSystemNSC sys;
sys.SetGravitationalAcceleration({0, 0, 0});
original_body = chrono_types::make_shared<ChBody>();
original_body->SetFixed(true);
sys.Add(original_body);
decomposition_body = chrono_types::make_shared<ChBody>();
decomposition_body->SetFixed(true);
sys.Add(decomposition_body);
// Default settings
hacd_maxhullcount = 512;
hacd_maxhullmerge = 256;
hacd_maxhullvertexes = 64;
hacd_concavity = 0.2;
hacd_smallclusterthreshold = 0.0;
hacd_fusetolerance = 1e-9;
// Create the Irrlicht visualization system
auto vis = chrono_types::make_shared<ChVisualSystemIrrlicht>();
vis->AttachSystem(&sys);
vis->SetWindowSize(800, 600);
vis->SetWindowTitle("Convex decomposition of a mesh");
vis->Initialize();
vis->AddLogo();
vis->AddSkyBox();
vis->AddCamera(ChVector3d(0, 1.5, -2));
// Create a custom event receiver for a GUI
MyEventReceiver receiver(vis.get());
vis->AddUserEventReceiver(&receiver);
// Rendering loop
while (vis->Run()) {
vis->BeginScene();
vis->Render();
vis->EndScene();
}
return 0;
}
static std::shared_ptr< ChTriangleMeshConnected > CreateFromWavefrontFile(const std::string &filename, bool load_normals=true, bool load_uv=false)
Create and return a ChTriangleMeshConnected from a Wavefront OBJ file.
Definition: ChTriangleMeshConnected.cpp:268
void Add(std::shared_ptr< ChPhysicsItem > item)
Attach an arbitrary ChPhysicsItem (e.g.
Definition: ChSystem.cpp:215
bool GetRootShape(TopoDS_Shape &shape, const int num=1) const
Get the root shape.
Definition: ChCascadeDoc.cpp:333
virtual void Render() override
Draw all 3D shapes and GUI elements at the current frame.
Definition: ChVisualSystemIrrlicht.cpp:602
virtual void Initialize() override
Initialize the visualization system.
Definition: ChVisualSystemIrrlicht.cpp:194
void AddSkyBox(const std::string &texture_dir=GetChronoDataFile("skybox/"))
Add a sky box in a 3D scene.
Definition: ChVisualSystemIrrlicht.cpp:389
virtual void Reset() override
Reset the input mesh data.
Definition: ChConvexDecomposition.cpp:282
void Dump(std::ostream &stream) const
Dump the shapes hierarchy on given stream.
Definition: ChCascadeDoc.cpp:308
virtual int ComputeConvexDecomposition() override
Perform the convex decomposition.
Definition: ChConvexDecomposition.cpp:333
virtual bool AddTriangleMesh(const ChTriangleMesh &tm)
Add a triangle mesh soup, by passing an entire ChTriangleMesh object.
Definition: ChConvexDecomposition.cpp:300
Class for wrapping the HACD convex decomposition code revisited by John Ratcliff.
Definition: ChConvexDecomposition.h:147
Class that contains an OCAF document (a tree hierarchy of shapes in the OpenCascade framework).
Definition: ChCascadeDoc.h:38
void SetGravitationalAcceleration(const ChVector3d &gacc)
Set the gravitational acceleration vector (default: [0, 0, 0]).
Definition: ChSystem.h:160
virtual void EndScene() override
End the scene draw at the end of each animation frame.
Definition: ChVisualSystemIrrlicht.cpp:592
virtual void BeginScene() override
Perform any necessary operations at the beginning of each rendering frame.
Definition: ChVisualSystemIrrlicht.cpp:575
virtual bool WriteConvexHullsAsChullsFile(std::ostream &mstream)
Write the convex decomposition to a ".chulls" file, where each hull is a sequence of x y z coords.
Definition: ChConvexDecomposition.cpp:77
virtual int AddCamera(const ChVector3d &pos, ChVector3d targ=VNULL) override
Add a camera in an Irrlicht 3D scene.
Definition: ChVisualSystemIrrlicht.cpp:303
void AddUserEventReceiver(irr::IEventReceiver *receiver)
Attach a custom event receiver to the application.
Definition: ChVisualSystemIrrlicht.cpp:470
virtual bool GetConvexHullResult(unsigned int hullIndex, ChTriangleMesh &convextrimesh) override
Get the n-th computed convex hull, by filling a ChTriangleMesh object that is passed as a parameter.
Definition: ChConvexDecomposition.cpp:394
ChVector3< double > ChVector3d
Alias for double-precision vectors.
Definition: ChVector3.h:287
void SetWindowTitle(const std::string &win_title)
Set the window title (default "").
Definition: ChVisualSystemIrrlicht.cpp:148
void SetParameters(unsigned int mMaxHullCount=256, unsigned int mMaxMergeHullCount=256, unsigned int mMaxHullVertices=64, float mConcavity=0.2f, float mSmallClusterThreshold=0.0f, float mFuseTolerance=1e-9)
Set the parameters for this convex decomposition algorithm.
Definition: ChConvexDecomposition.cpp:308
virtual void AttachSystem(ChSystem *sys) override
Attach another Chrono system to the run-time visualization system.
Definition: ChVisualSystemIrrlicht.cpp:175
virtual unsigned int GetHullCount() override
Get the number of computed hulls after the convex decomposition.
Definition: ChConvexDecomposition.cpp:373
virtual bool Run() override
Run the Irrlicht device.
Definition: ChVisualSystemIrrlicht.cpp:256
irr::scene::ILightSceneNode * AddLight(const ChVector3d &pos, double radius, ChColor color=ChColor(0.7f, 0.7f, 0.7f))
Add a point light to the scene.
Definition: ChVisualSystemIrrlicht.cpp:434
Class for storing settings on OpenCASCADE tessellation of shapes.
Definition: ChCascadeTriangulate.h:25
virtual void BindItem(std::shared_ptr< ChPhysicsItem > item) override
Process the visual assets for the specified physics item.
Definition: ChVisualSystemIrrlicht.cpp:652
bool LoadSTEP(const std::string &filename)
Populate the document with all shapes that are contained in the STEP file.
Definition: ChCascadeDoc.cpp:293
Class for a physical system in which contact is modeled using a non-smooth (complementarity-based) me...
Definition: ChSystemNSC.h:29
void AddLogo(const std::string &logo_filename=GetChronoDataFile("logo_chrono_alpha.png"))
Add a logo in a 3D scene.
Definition: ChVisualSystemIrrlicht.cpp:368
virtual void WriteConvexHullsAsWavefrontObj(std::ostream &mstream) override
Save the computed convex hulls as a Wavefront file using the '.obj' file format, with each hull as a ...
Definition: ChConvexDecomposition.cpp:418
virtual void UnbindItem(std::shared_ptr< ChPhysicsItem > item) override
Remove the visual assets for the specified physics item from this visualization system.
Definition: ChVisualSystemIrrlicht.cpp:669
void SetWindowSize(unsigned int width, unsigned int height)
Set the window size (default 640x480).
Definition: ChVisualSystemIrrlicht.cpp:144
Irrlicht-based Chrono run-time visualization system.
Definition: ChVisualSystemIrrlicht.h:55