Load a STEP file (demo_CAS_stepfile.cpp)

Tutorial that teaches how to use the CASCADE module to load CAD models from STEP files.

// =============================================================================
// 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.
//
// =============================================================================
// Show how to use the OpenCASCADE features
// implemented in the unit_CASCADE:
//
// - load a 3D model saved in STEP format from a CAD
// - select some sub assemblies from the STEP model
// - make Chrono::Engine objects out of those parts
// =============================================================================
#include "chrono/core/ChRealtimeStep.h"
#include "chrono/physics/ChSystemNSC.h"
#include "chrono/physics/ChBodyEasy.h"
#include "chrono_cascade/ChBodyEasyCascade.h"
#include "chrono_cascade/ChCascadeDoc.h"
#include "chrono_cascade/ChCascadeShapeAsset.h"
#include "chrono_irrlicht/ChIrrApp.h"
// Use the namespace of Chrono
using namespace chrono;
using namespace chrono::irrlicht;
// Use the main namespaces of Irlicht
using namespace irr;
using namespace core;
using namespace scene;
using namespace video;
using namespace io;
using namespace gui;
// Use the namespace with OpenCascade stuff
using namespace cascade;
//
// This is the program which is executed
//
int main(int argc, char* argv[]) {
// Create a ChronoENGINE physical system: all bodies and constraints
// will be handled by this ChSystemNSC object.
ChSystemNSC my_system;
// Create the Irrlicht visualization (open the Irrlicht device,
// bind a simple user interface, etc. etc.)
ChIrrApp application(&my_system, L"Load a STEP model from file", core::dimension2d<u32>(800, 600), false, true);
// Easy shortcuts to add logo, camera, lights and sky in Irrlicht scene:
ChIrrWizard::add_typical_Logo(application.GetDevice());
ChIrrWizard::add_typical_Sky(application.GetDevice());
ChIrrWizard::add_typical_Lights(application.GetDevice(), core::vector3df(30, 100, 30),
core::vector3df(30, -80, -30), 200, 130);
ChIrrWizard::add_typical_Camera(application.GetDevice(), core::vector3df(0.2f, 0.2f, -0.3f));
//
// Load a STEP file, containing a mechanism. The demo STEP file has been
// created using a 3D CAD (in this case, SolidEdge v.18).
//
// Create the ChCascadeDoc, a container that loads the STEP model
// and manages its subassembles
ChCascadeDoc mydoc;
// load the STEP model using this command:
bool load_ok = mydoc.Load_STEP(GetChronoDataFile("cascade/assembly.stp")
.c_str()); // or specify abs.path: ("C:\\data\\cascade\\assembly.stp");
// print the contained shapes
mydoc.Dump(GetLog());
// In most CADs the Y axis is horizontal, but we want it vertical.
// So define a root transformation for rotating all the imported objects.
ChQuaternion<> rotation1;
rotation1.Q_from_AngAxis(-CH_C_PI / 2, ChVector<>(1, 0, 0)); // 1: rotate 90� on X axis
ChQuaternion<> rotation2;
rotation2.Q_from_AngAxis(CH_C_PI, ChVector<>(0, 1, 0)); // 2: rotate 180� on vertical Y axis
ChQuaternion<> tot_rotation = rotation2 % rotation1; // rotate on 1 then on 2, using quaternion product
ChFrameMoving<> root_frame(ChVector<>(0, 0, 0), tot_rotation);
// Retrieve some sub shapes from the loaded model, using
// the GetNamedShape() function, that can use path/subpath/subsubpath/part
// syntax and * or ? wildcards, etc.
std::shared_ptr<ChBodyEasyCascade> mrigidBody1;
std::shared_ptr<ChBodyEasyCascade> mrigidBody2;
if (load_ok) {
TopoDS_Shape shape1;
if (mydoc.GetNamedShape(shape1, "Assem1/body1")) {
// Create the ChBody using the ChBodyEasyCascade helper:
auto mbody1 = chrono_types::make_shared<ChBodyEasyCascade>(shape1,
1000, // density
true, // add a visualization
false // add a collision model
);
my_system.Add(mbody1);
mbody1->SetBodyFixed(true);
// Move the body as for global displacement/rotation (also mbody1 %= root_frame; )
mbody1->ConcatenatePreTransformation(root_frame);
mrigidBody1 = mbody1;
} else
GetLog() << "Warning. Desired object not found in document \n";
TopoDS_Shape shape2;
if (mydoc.GetNamedShape(shape2, "Assem1/body2")) {
// Create the ChBody using the ChBodyEasyCascade helper (with more detailed info on visualization tesselation):
ChCascadeTriangulateTolerances mvisualtolerances ( 0.02, // chordal deflection for triangulation
false, // chordal deflection is relative
0.5 // angular deflection for triangulation
);
auto mbody2 = chrono_types::make_shared<ChBodyEasyCascade>(shape2,
1000, // density
mvisualtolerances, // add a visualization
false // add a collision model
);
my_system.Add(mbody2);
// Move the body as for global displacement/rotation (also mbody2 %= root_frame; )
mbody2->ConcatenatePreTransformation(root_frame);
mrigidBody2 = mbody2;
} else
GetLog() << "Warning. Desired object not found in document \n";
} else
GetLog() << "Warning. Desired STEP file could not be opened/parsed \n";
// Create a revolute joint between the two parts
// as in a pendulum. We assume we already know in advance
// the aboslute position of the joint (ex. we used measuring tools in the 3D CAD)
ChVector<> measured_joint_pos_mm(0, 48, 120);
double scale = 1. / 1000.; // because we use meters instead of mm
ChVector<> joint_pos =
((ChFrame<>)root_frame) * (measured_joint_pos_mm * scale); // transform because we rotated everything
if (mrigidBody1 && mrigidBody2) {
std::shared_ptr<ChLinkLockRevolute> my_link(new ChLinkLockRevolute);
my_link->Initialize(mrigidBody1, mrigidBody2, ChCoordsys<>(joint_pos));
my_system.AddLink(my_link);
}
// Create a large cube as a floor.
std::shared_ptr<ChBodyEasyBox> mfloor(new ChBodyEasyBox(1, 0.2, 1, 1000));
mfloor->SetPos(ChVector<>(0, -0.3, 0));
mfloor->SetBodyFixed(true);
application.GetSystem()->Add(mfloor);
std::shared_ptr<ChColorAsset> mcolor(new ChColorAsset(0.3f, 0.3f, 0.8f));
mfloor->AddAsset(mcolor);
// Use this function for adding a ChIrrNodeAsset to all items
// Otherwise use application.AssetBind(myitem); on a per-item basis.
application.AssetBindAll();
// Use this function for 'converting' assets into Irrlicht meshes
application.AssetUpdateAll();
//
// THE SOFT-REAL-TIME CYCLE, SHOWING THE SIMULATION
//
application.SetTimestep(0.01);
application.SetTryRealtime(true);
while (application.GetDevice()->run()) {
// Irrlicht must prepare frame to draw
application.BeginScene(true, true, video::SColor(255, 140, 161, 192));
// .. draw solid 3D items (boxes, cylinders, shapes) belonging to Irrlicht scene, if any
application.DrawAll();
application.DoStep();
application.EndScene();
}
return 0;
}
std::string GetChronoDataFile(const std::string &filename)
Obtain the complete path to the specified filename, given relative to the Chrono data directory (thre...
Definition: ChGlobal.cpp:95
void Add(std::shared_ptr< ChPhysicsItem > item)
Attach an arbitrary ChPhysicsItem (e.g.
Definition: ChSystem.cpp:146
static void SetDefaultSuggestedMargin(double mmargin)
Using this function BEFORE you start creating collision shapes, it will make all following collision ...
Definition: ChCollisionModel.cpp:44
Easy-to-use class for quick creation of rigid bodies with a box shape.
Definition: ChBodyEasy.h:102
static void add_typical_Camera(irr::IrrlichtDevice *device, irr::core::vector3df mpos=irr::core::vector3df(0, 0, -8), irr::core::vector3df mtarg=irr::core::vector3df(0, 0, 0))
A very basic and simple function which is just a shortcut to avoid lot of typing when someone wants t...
Definition: ChIrrWizard.cpp:53
COORDSYS:
Definition: ChCoordsys.h:38
ChLog & GetLog()
Global function to get the current ChLog object.
Definition: ChLog.cpp:39
virtual void AddLink(std::shared_ptr< ChLinkBase > link)
Attach a link to the underlying assembly.
Definition: ChSystem.h:268
void SetTimestep(double val)
Set/Get the time step for time integration.
Definition: ChIrrAppInterface.cpp:558
static void add_typical_Logo(irr::IrrlichtDevice *device, const std::string &mlogofilename=GetChronoDataFile("logo_chronoengine_alpha.png"))
A very basic and simple function which is just a shortcut to avoid lot of typing when someone wants t...
Definition: ChIrrWizard.cpp:20
static void add_typical_Sky(irr::IrrlichtDevice *device, const std::string &mtexturedir=GetChronoDataFile("skybox/"))
A very basic and simple function which is just a shortcut to avoid lot of typing when someone wants t...
Definition: ChIrrWizard.cpp:40
static void add_typical_Lights(irr::IrrlichtDevice *device, irr::core::vector3df pos1=irr::core::vector3df(30.f, 100.f, 30.f), irr::core::vector3df pos2=irr::core::vector3df(30.f, 80.f, -30.f), double rad1=290, double rad2=190, irr::video::SColorf col1=irr::video::SColorf(0.7f, 0.7f, 0.7f, 1.0f), irr::video::SColorf col2=irr::video::SColorf(0.7f, 0.8f, 0.8f, 1.0f))
A very basic and simple function which is just a shortcut to avoid lot of typing when someone wants t...
Definition: ChIrrWizard.cpp:25
Base class for assets that carry basic informations about the surface color for visualization assets.
Definition: ChColorAsset.h:28
ChFrame: a class for coordinate systems in 3D space.
Definition: ChFrame.h:42
Class to add some GUI to Irrlicht+ChronoEngine applications.
Definition: ChIrrApp.h:29
virtual void EndScene()
Call this to end the scene draw at the end of each animation frame.
Definition: ChIrrAppInterface.cpp:634
Namespace with classes for the Irrlicht module.
Definition: ChApiIrr.h:48
Definition of general purpose 3d vector variables, such as points in 3D.
Definition: ChVector.h:35
virtual void DoStep()
Call this function inside a loop such as.
Definition: ChIrrAppInterface.cpp:644
ChFrameMoving: a class for coordinate systems in 3D space.
Definition: ChFrameMoving.h:38
Class defining quaternion objects, that is four-dimensional numbers, also known as Euler parameters.
Definition: ChQuaternion.h:45
static void SetDefaultSuggestedEnvelope(double menv)
Using this function BEFORE you start creating collision shapes, it will make all following collision ...
Definition: ChCollisionModel.cpp:40
virtual void BeginScene(bool backBuffer=true, bool zBuffer=true, irr::video::SColor color=irr::video::SColor(255, 0, 0, 0))
Call this to clean the canvas at the beginning of each animation frame.
Definition: ChIrrAppInterface.cpp:617
Main namespace for the Chrono package.
Definition: ChAsset.cpp:18
void Q_from_AngAxis(Real angle, const ChVector< Real > &axis)
Set the quaternion from an angle of rotation and an axis, defined in absolute coords.
Definition: ChQuaternion.h:1129
virtual void DrawAll()
Call this function inside a loop such as.
Definition: ChIrrAppInterface.cpp:757
Class for a physical system in which contact is modeled using a non-smooth (complementarity-based) me...
Definition: ChSystemNSC.h:29
void SetTryRealtime(bool val)
If enabled, the function DoStep() will enforce soft real-time, by spinning in place until simulation ...
Definition: ChIrrAppInterface.h:94