world.h
Go to the documentation of this file.
1 #ifndef WORLDH
2 #define WORLDH
3 
4 #include "vector.h"
5 #include "environment.h"
6 #include "mechanism.h"
7 #include "cuiksystem.h"
8 #include "equations.h"
9 #include "plot3d.h"
10 #include "box.h"
11 #include "csmapping.h"
12 #include "parameters.h"
13 #include "filename.h"
14 #include "statistics.h"
15 #include "list.h"
16 #include "jacobian.h"
17 #include "cd.h"
18 
19 #include <stdlib.h>
20 
21 
32 /*****************************************************************************/
33 
39 #define INITIAL_FRAME_DELAY 5.0
40 
45 #define FRAME_RATE (1.0/10)
46 
53 #define FINAL_FRAME_DELAY 500.0
54 /*****************************************************************************/
55 /*****************************************************************************/
56 
61 #define WORLD_IN_DEFINITION 0
62 
67 #define WORLD_DEFINED 1
68 
69 /*****************************************************************************/
70 
80 typedef struct {
82  double sign;
86 } TBranchStep;
87 
108 typedef struct {
110 } Tbranch;
111 
124 typedef struct {
125  unsigned int stage;
133  unsigned int nl;
134  unsigned int nb;
135  unsigned int nj;
136  unsigned int no;
137  unsigned int np;
140  unsigned int endEffector;
142  double maxCoord;
145  boolean **checkCollisionsLL;
154  boolean **checkCollisionsLO;
167  unsigned int nvars;
168  unsigned int nsvars;
169  boolean *systemVars;
174  unsigned int *sortedLinks;
176  unsigned int *jointTo;
203  boolean *openLink;
210  unsigned int nwcd;
216  unsigned int ndof;
217  unsigned int *joint2dof;
218  unsigned int *dof2joint;
219  unsigned int *dof2range;
221 } Tworld;
222 
223 
238 void InitWorld(Tworld *w);
239 
257  Tworld *w);
258 
274 unsigned int AddLink2World(Tlink *l,boolean object,Tworld *w);
275 
276 
287 unsigned int AddJoint2World(Tjoint *j,Tworld *w);
288 
298 void AddObstacle2World(char *name,Tpolyhedron *o,Tworld *w);
299 
313 unsigned int GetWorldMobility(Tworld *w);
314 
325 
337 unsigned int GetWorldLinkID(char *linkName,Tworld *w);
338 
350 unsigned int GetWorldObstacleID(char *obsName,Tworld *w);
351 
363 Tlink *GetWorldLink(unsigned int linkID,Tworld *w);
364 
365 
377 Tjoint *GetWorldJoint(unsigned int jointID,Tworld *w);
378 
388 unsigned int GetWorldNLinks(Tworld *w);
389 
399 unsigned int GetWorldNJoints(Tworld *w);
400 
412 unsigned int GetWorldNObstacles(Tworld *w);
413 
426 unsigned int GetWorldNConvexBodiesInLinks(Tworld *w);
427 
442 unsigned int GetWorldNConvexBodies(Tworld *w);
443 
444 
458 unsigned int GetWorldNDOF(Tworld *w);
459 
460 
470 void GetWorldRangeDOF(unsigned int ndof,Tinterval *r,Tworld *w);
471 
472 
490 void GetWorldJointLabel(unsigned int ndof,char **string,Tworld *w);
491 
507 void GetWorldDOFLabel(unsigned int ndof,char **string,Tworld *w);
508 
509 
520 boolean IsWorldPolynomial(Tworld *w);
521 
534 void CheckAllCollisions(unsigned int fl,unsigned int fo,Tworld *w);
535 
549 void NoCheckAllCollisions(unsigned int fl,unsigned int fo,Tworld *w);
550 
560 boolean AnyCollision(Tworld *w);
561 
571 void PrintCollisions(FILE *f,Tworld *w);
572 
583 void CheckSelfCollisions(unsigned int fl,Tworld *w);
584 
595 void NoCheckSelfCollisions(unsigned int fl,Tworld *w);
596 
613 void CheckLinkLinkCollision(unsigned int a,unsigned int b,Tworld *w);
614 
627 void NoCheckLinkLinkCollision(unsigned int a,unsigned int b,Tworld *w);
628 
646 void CheckLinkObstacleCollision(unsigned int a,unsigned int b,Tworld *w);
647 
661 void NoCheckLinkObstacleCollision(unsigned int a,unsigned int b,Tworld *w);
662 
681 void InitWorldCD(Tparameters *pr,unsigned int mt,Tworld *w);
682 
695 boolean WorldCanCollide(Tworld *w);
696 
713 void GetLinkTransformsFromSolutionPoint(Tparameters *p,boolean simp,double *sol,
714  THTransform **tl,Tworld *w);
715 
734  double **sol,Tworld *w);
735 
746 
759 boolean WorldContinuousCD(Tworld *w);
760 
782 boolean MoveAndCheckCDFromTransforms(boolean all,unsigned int tID,
783  THTransform *tl,THTransform *tlPrev,Tworld *w);
784 
813 boolean MoveAndCheckCD(Tparameters *p,boolean all,unsigned int tID,
814  double *sol,double *solPrev,Tworld *w);
815 
816 
844 boolean NewtonInWorld(Tparameters *p,double *v,Tbox *b_sol,Tworld *w);
845 
856 unsigned int GetWorldSystemVars(boolean **sv,Tworld *w);
857 
868 unsigned int GetWorldVarTopology(unsigned int vID,Tworld *w);
869 
888 unsigned int GetWorldSimpVariableMask(Tparameters *p,boolean **sv,Tworld *w);
889 
890 
903 void GetWorldVarNames(char **vn,Tworld *w);
904 
920 
950  Tworld *w);
951 
970  Tworld *w);
971 
986 void GetWorldJacobian(TJacobian *J,Tworld *w);
987 
1004 
1021 unsigned int GetWorldSimpTopology(Tparameters *p,
1022  unsigned int **t,Tworld *w);
1023 
1043 void GetWorldKinJacobian(Tparameters *p,unsigned int *nr,unsigned int *nc,
1044  Tequation ***J,Tworld *w);
1060 unsigned int WorldSimpCuikNewton(Tparameters *p,double *pt,Tworld *w);
1061 
1078 void EvaluateWorldJacobian(double *p,double ***J,Tworld *w);
1079 
1092 void WorldEvaluateEquations(double *pt,double *r,Tworld *w);
1093 
1105 boolean WorldInequalitiesHold(double *pt,Tworld *w);
1106 
1119 boolean WorldSimpInequalitiesHold(Tparameters *p,double *pt,Tworld *w);
1120 
1132 double WorldErrorInSimpInequalities(Tparameters *p,double *pt,Tworld *w);
1133 
1147 void WorldEvaluateSimpEquations(Tparameters *p,double *pt,double *r,
1148  Tworld *w);
1149 
1164 void WorldEvaluateSubSetSimpEquations(Tparameters *p,boolean *se,double *pt,double *r,
1165  Tworld *w);
1166 
1181 double WorldErrorInEquations(double *pt,Tworld *w);
1182 
1198 double WorldErrorInSimpEquations(Tparameters *p,double *pt,Tworld *w);
1199 
1213 double WorldErrorFromDOFs(Tparameters *p,double *dof,Tworld *w);
1214 
1226 double EvaluateWorldCost(Tparameters *p,boolean simp,double *pt,void *w);
1227 
1238 unsigned int GetWorldNumVariables(Tworld *w);
1239 
1250 unsigned int GetWorldNumSystemVariables(Tworld *w);
1251 
1269 unsigned int RegenerateWorldOriginalPoint(Tparameters *p,double *s,
1270  double **o,Tworld *w);
1271 
1289 unsigned int WorldGenerateSimplifiedPoint(Tparameters *p,double *o,
1290  double **s,Tworld *w);
1291 
1303 void GetWorldInitialBox(Tbox *b,Tworld *w);
1304 
1315 
1316 
1334 unsigned int MaxWorldReduction(Tparameters *p,Tbox *b,
1335  double *reduction,Tworld *w);
1336 
1351 void PlotWorld(Tparameters *pr,Tplot3d *pt,double axesLength,Tworld *w);
1352 
1367 unsigned int RegenerateWorldSolutionPoint(Tparameters *pr,double *p,
1368  double **v,Tworld *w);
1369 
1400 void WorldPrintAtoms(Tparameters *pr,FILE *f,double *pt,Tworld *w);
1401 
1419 void WorldStoreRigidGroups(Tparameters *pr,FILE *f,double *pt,Tworld *w);
1420 
1441 void WorldAtomJacobian(Tparameters *pr,double *sol,unsigned int *nr,unsigned int *nc,double ***J,Tworld *w);
1442 
1465 void MoveWorld(Tparameters *pr,Tplot3d *pt,Tbox *b,Tworld *w);
1466 
1480 void PrintWorldAxes(Tparameters *pr,FILE *f,Tbox *b,Tworld *w);
1481 
1500 void PrintWorldCollisionInfo(FILE *f,char *fname,Tworld *w);
1501 
1521 void MoveWorldDOF(Tparameters *pr,Tplot3d *pt,double *dof,Tworld *w);
1522 
1541 unsigned int WorldDOF2Sol(Tparameters *p,double *dof,double **sol,Tbox *b,Tworld *w);
1542 
1559 void WorldSample2DOF(Tparameters *p,double *sample,double *dof,Tworld *w);
1560 
1585 void AnimateWorld(Tparameters *pr,char *pname,double axesLength,
1586  double frameDelay,Tlist *p,Tworld *w);
1587 
1598 void PrintWorldCS(Tparameters *p,Tfilename *fname,Tworld *w);
1599 
1614 void PrintWorld(char *fname,int argc,char **arg,Tworld *w);
1615 
1623 void DeleteWorld(Tworld *w);
1624 
1625 #endif
unsigned int RegenerateWorldSolutionPoint(Tparameters *pr, double *p, double **v, Tworld *w)
Computes the missing values in a kinematic solution.
Definition: world.c:2276
boolean IsWorldPolynomial(Tworld *w)
Checks if the system of equations is polynomial.
Definition: world.c:1607
double WorldErrorInSimpInequalities(Tparameters *p, double *pt, Tworld *w)
Determines the maximum error in the inequalites for the simplified system.
Definition: world.c:2089
double WorldErrorInEquations(double *pt, Tworld *w)
Evalates the norm of the error in the kinematic equations for a given point.
Definition: world.c:2113
unsigned int nl
Definition: world.h:133
void MoveWorld(Tparameters *pr, Tplot3d *pt, Tbox *b, Tworld *w)
Moves the mechanisms defined in a world information to a given configuration.
Definition: world.c:2387
void MoveWorldDOF(Tparameters *pr, Tplot3d *pt, double *dof, Tworld *w)
Moves a mechanisms to a configuration given by the degrees of freedom.
Definition: world.c:2514
boolean WorldCanCollide(Tworld *w)
Determines if any collision is potentially possible.
Definition: world.c:1046
void WorldStoreRigidGroups(Tparameters *pr, FILE *f, double *pt, Tworld *w)
Generates a file with the atoms grouped in rigid clusters.
Definition: world.c:2318
void GetWorldRangeDOF(unsigned int ndof, Tinterval *r, Tworld *w)
Gets the range for a given degree of freedom.
Definition: world.c:1564
void PrintWorld(char *fname, int argc, char **arg, Tworld *w)
Prints the world.
Definition: world.c:2742
Tvector steps
Definition: world.h:109
Relation between two links.
Definition: joint.h:177
boolean ** checkCollisionsLL
Definition: world.h:145
void PrintWorldCollisionInfo(FILE *f, char *fname, Tworld *w)
Prints information collected durint last collision check.
Definition: world.c:2508
TCuikSystem kinCS
Definition: world.h:165
Data structure to hold the information about the name of a file.
Definition: filename.h:248
Auxiliary data for the collision detection.
Definition: cd.h:257
void WorldEvaluateEquations(double *pt, double *r, Tworld *w)
Evaluates the kinematic equations.
Definition: world.c:2065
boolean * systemVars
Definition: world.h:169
A step in a kinematic branch.
Definition: world.h:80
unsigned int WorldDOF2Sol(Tparameters *p, double *dof, double **sol, Tbox *b, Tworld *w)
Transforms from degrees of freedom to cuik variables.
Definition: world.c:2543
A homgeneous transform in R^3.
Definition of the TJacobian type and the associated functions.
boolean MoveAndCheckCD(Tparameters *p, boolean all, unsigned int tID, double *sol, double *solPrev, Tworld *w)
Checks a point for collision.
Definition: world.c:1076
Tmechanism m
Definition: world.h:131
boolean WorldInequalitiesHold(double *pt, Tworld *w)
Check if the inequalities hold.
Definition: world.c:2073
boolean WorldContinuousCD(Tworld *w)
Determines the type of collision library used.
Definition: world.c:1051
unsigned int GetWorldVarTopology(unsigned int vID, Tworld *w)
Get the topology of a given variable.
Definition: world.c:1818
void InitWorldFromFile(Tparameters *p, Tfilename *f, Tworld *w)
Constructor.
boolean MoveAndCheckCDFromTransforms(boolean all, unsigned int tID, THTransform *tl, THTransform *tlPrev, Tworld *w)
Checks a point for collision.
Definition: world.c:1060
void EvaluateWorldJacobian(double *p, double ***J, Tworld *w)
Evaluates the kinematic Jacobian.
Definition: world.c:2057
void GetWorldKinJacobian(Tparameters *p, unsigned int *nr, unsigned int *nc, Tequation ***J, Tworld *w)
Gets the simplified kinematic Jacobian.
void GetWorldJacobian(TJacobian *J, Tworld *w)
Gets the kinematic Jacobian.
Definition: world.c:2027
void WorldSample2DOF(Tparameters *p, double *sample, double *dof, Tworld *w)
Transforms a sample degrees of freedom.
Definition: world.c:2625
Tenvironment e
Definition: world.h:130
void GetWorldVarNames(char **vn, Tworld *w)
Return the variable names.
Definition: world.c:1842
double WorldErrorFromDOFs(Tparameters *p, double *dof, Tworld *w)
Error in equations from DOFs.
Definition: world.c:2129
void WorldPrintAtoms(Tparameters *pr, FILE *f, double *pt, Tworld *w)
Generates a file with the atom centers in global coordiantes.
Definition: world.c:2304
void GetWorldDOFLabel(unsigned int ndof, char **string, Tworld *w)
Returns a label identifying each degree of freedom.
Definition: world.c:1585
Definition of the Tfilename type and the associated functions.
Definition of the Tplot3d type and the associated functions.
boolean ** checkCollisionsLO
Definition: world.h:154
void WorldEvaluateSubSetSimpEquations(Tparameters *p, boolean *se, double *pt, double *r, Tworld *w)
Evaluates a subset of the simplified kinematic equations.
Definition: world.c:2105
unsigned int WorldSimpCuikNewton(Tparameters *p, double *pt, Tworld *w)
Tries to reach the kinematic manifold.
Definition: world.c:2052
void PrintWorldAxes(Tparameters *pr, FILE *f, Tbox *b, Tworld *w)
Prints the axes of the mechanism.
Definition: world.c:2455
unsigned int np
Definition: world.h:137
unsigned int GetWorldNumSystemVariables(Tworld *w)
Number of system variables in the kinematic subsystem.
Definition: world.c:2158
void NoCheckLinkObstacleCollision(unsigned int a, unsigned int b, Tworld *w)
Desactivates the possible collision between a particular link and an object in the environment...
Definition: world.c:1712
All the necessary information to generate equations for mechanisms.
Definition: world.h:124
double WorldErrorInSimpEquations(Tparameters *p, double *pt, Tworld *w)
Evalates the norm of the error in the kinematic equations for a given point.
Definition: world.c:2121
A mechanism description.
Definition: mechanism.h:51
unsigned int GetWorldNObstacles(Tworld *w)
Gets the number of obstacles in the environment included in the world.
Definition: world.c:1467
void NoCheckAllCollisions(unsigned int fl, unsigned int fo, Tworld *w)
Desactivates all the possible collision between links and links and obstacles.
Definition: world.c:1637
unsigned int * dof2range
Definition: world.h:219
void WorldAtomJacobian(Tparameters *pr, double *sol, unsigned int *nr, unsigned int *nc, double ***J, Tworld *w)
Jacobian of the atom position w.r.t. the variables.
Definition: world.c:2332
Definition of the Tstatistics type and the associated functions.
void GetWorldInitialBox(Tbox *b, Tworld *w)
Gets the kinematic search space for a given problem.
Definition: world.c:2204
unsigned int no
Definition: world.h:136
void DeleteWorld(Tworld *w)
Destructor.
Definition: world.c:2792
A polyhedron.
Definition: polyhedron.h:124
unsigned int * joint2dof
Definition: world.h:217
unsigned int MaxWorldReduction(Tparameters *p, Tbox *b, double *reduction, Tworld *w)
Reduces the system variables as much as possible using the kinematic constraints. ...
Definition: world.c:2220
A 3D plot.
Definition: plot3d.h:54
double EvaluateWorldCost(Tparameters *p, boolean simp, double *pt, void *w)
Evaluates the functions cost defined in a world.
Definition: world.c:2142
boolean IsMechanismInWorldAllSpheres(Tworld *w)
TRUE if the mechanisms in the world is based on spheres.
Definition: world.c:1432
Tjoint * joint
Definition: world.h:81
void CheckAllCollisions(unsigned int fl, unsigned int fo, Tworld *w)
Activates all the possible collision between links and links and obstacles.
Definition: world.c:1615
boolean NewtonInWorld(Tparameters *p, double *v, Tbox *b_sol, Tworld *w)
Generates a sample from a the kinematic cuiksystem in the world using the Newton algorithm.
Definition: world.c:1762
unsigned int nj
Definition: world.h:135
unsigned int GetWorldObstacleID(char *obsName, Tworld *w)
Gets the identifier of an obstacle from its name.
Definition: world.c:1442
void GenerateWorldEquations(Tparameters *p, Tworld *w)
Generates all the cuiksystems derived from the world information.
Definition: world.c:1850
void CheckSelfCollisions(unsigned int fl, Tworld *w)
Activates all the possible collision between links.
Definition: world.c:1656
void AnimateWorld(Tparameters *pr, char *pname, double axesLength, double frameDelay, Tlist *p, Tworld *w)
Produces an animation along a path.
Definition: world.c:2677
Set of equations.
Definition: equations.h:81
An equation.
Definition: equation.h:236
A collection of obstacles (convex polyhedrons) with their names.
Definition: environment.h:39
void AddObstacle2World(char *name, Tpolyhedron *o, Tworld *w)
Adds an obstacle to the environment in the world.
Definition: world.c:1393
Definition of the Tbox type and the associated functions.
unsigned int AddLink2World(Tlink *l, boolean object, Tworld *w)
Adds a link to the mechanism in the world.
Definition: world.c:1333
A generic list.
Definition: list.h:46
TworldCD * wcd
Definition: world.h:214
unsigned int nsvars
Definition: world.h:168
Definition of the Tlist type and the associated functions.
boolean WorldSimpInequalitiesHold(Tparameters *p, double *pt, Tworld *w)
Check if the inequalities hold for the simplified system.
Definition: world.c:2081
void GenerateWorldSingularityEquations(Tparameters *p, char *ln, TCuikSystem *cs, Tworld *w)
Adds variables/equations to detect singularities.
Definition: world.c:1929
void PrintWorldCS(Tparameters *p, Tfilename *fname, Tworld *w)
Prints the cuiksystems derived from a world.
Definition: world.c:2721
unsigned int * dof2joint
Definition: world.h:218
Definition of the Tmechanism type and the associated functions.
unsigned int nvars
Definition: world.h:167
void GetWorldJointLabel(unsigned int ndof, char **string, Tworld *w)
Returns a label identifying each joint.
Definition: world.c:1575
unsigned int GetWorldNLinks(Tworld *w)
Gets the number of links in the mechanism included in the world.
Definition: world.c:1457
unsigned int nb
Definition: world.h:134
unsigned int GetWorldSimpTopology(Tparameters *p, unsigned int **t, Tworld *w)
Gets the topology of the variables.
Definition: world.c:2043
A table of parameters.
Definition of the TCuikSystem type and the associated functions.
unsigned int AddJoint2World(Tjoint *j, Tworld *w)
Adds a joint to the mechanism in the world.
Definition: world.c:1383
boolean * openLink
Definition: world.h:203
boolean AnyCollision(Tworld *w)
Determines if we want to avoid any collision.
Definition: world.c:1721
Tequations refEqs
Definition: world.h:186
void WorldEvaluateSimpEquations(Tparameters *p, double *pt, double *r, Tworld *w)
Evaluates the simplified kinematic equations.
Definition: world.c:2097
A generic vector.
Definition: vector.h:227
unsigned int GetWorldNConvexBodies(Tworld *w)
Gets the number of convex polyhedrons that define the mechanism and the environment included in the w...
Definition: world.c:1477
A kinematic branch.
Definition: world.h:108
TJacobian refEqsJ
Definition: world.h:196
void CheckLinkLinkCollision(unsigned int a, unsigned int b, Tworld *w)
Activates the possible collision between a particular pair of links.
Definition: world.c:1685
unsigned int stage
Definition: world.h:125
A box.
Definition: box.h:83
void GenerateWorldTWSEquations(Tparameters *p, char *ln, TCuikSystem *cs, Tworld *w)
Adds variables/equations to detect translational workspace boundaries.
Definition: world.c:1994
A cuiksystem, i.e., a set of variables and equations defining a position analysis problem...
Definition: cuiksystem.h:181
unsigned int GetWorldLinkID(char *linkName, Tworld *w)
Gets the identifier of a link from its name.
Definition: world.c:1437
void DeleteLinkTransforms(THTransform *tl, Tworld *w)
Deletes transforms for each link.
Definition: world.c:1290
unsigned int RegenerateWorldOriginalPoint(Tparameters *p, double *s, double **o, Tworld *w)
Reconstruct a point in original kinematic system from a simplified point.
Definition: world.c:2166
Definition of the Tvector type and the associated functions.
void GetLinkTransformsFromSolutionPoint(Tparameters *p, boolean simp, double *sol, THTransform **tl, Tworld *w)
Define transforms for the links from the a solution point.
Definition: world.c:1165
void PlotWorld(Tparameters *pr, Tplot3d *pt, double axesLength, Tworld *w)
Adds a world (environment plus mechanism) in a 3D scene.
Definition: world.c:2260
unsigned int GetWorldSystemVars(boolean **sv, Tworld *w)
Gets the system vars of the kinematic cuiksystem.
Definition: world.c:1810
Tlink * GetWorldLink(unsigned int linkID, Tworld *w)
Gets a link from its identifier.
Definition: world.c:1447
void NoCheckLinkLinkCollision(unsigned int a, unsigned int b, Tworld *w)
Desactivates the possible collision between a particular pair of links.
Definition: world.c:1694
Tjoint * GetWorldJoint(unsigned int jointID, Tworld *w)
Gets a joint from its identifier.
Definition: world.c:1452
unsigned int GetWorldNJoints(Tworld *w)
Gets the number of joints in the mechanism included in the world.
Definition: world.c:1462
Definition of the Tenvironment type and the associated functions.
The Jacobian of a set of equations.
Definition: jacobian.h:23
unsigned int ndof
Definition: world.h:216
Tbranch * branch2Link
Definition: world.h:178
void InitWorld(Tworld *w)
Constructor.
Definition: world.c:1308
void NoCheckSelfCollisions(unsigned int fl, Tworld *w)
Desactivates all the possible collision between links.
Definition: world.c:1671
Headers of the interface with the collision detection engine.
unsigned int GetWorldSimpVariableMask(Tparameters *p, boolean **sv, Tworld *w)
Identifies pose related variable that survied in the simplified system.
Definition: world.c:1823
unsigned int * sortedLinks
Definition: world.h:174
void GetWorldSimpJacobian(Tparameters *p, TJacobian *J, Tworld *w)
Gets the simplified kinematic Jacobian.
Definition: world.c:2035
unsigned int endEffector
Definition: world.h:140
unsigned int WorldGenerateSimplifiedPoint(Tparameters *p, double *o, double **s, Tworld *w)
Reconstruct a point in simplified kinematic system from an original point.
Definition: world.c:2196
Definition of the Tmapping type and the associated functions.
unsigned int GetWorldNDOF(Tworld *w)
Gets the number of degrees of freedom in the world.
Definition: world.c:1559
void InitWorldCD(Tparameters *pr, unsigned int mt, Tworld *w)
Initializes the collision detector.
Definition: world.c:1013
void PrintCollisions(FILE *f, Tworld *w)
Stores the collision information into a file.
Definition: world.c:1733
unsigned int GetWorldNumVariables(Tworld *w)
Number of variables in the kinematic subsystem.
Definition: world.c:2150
Defines a interval.
Definition: interval.h:33
unsigned int GetWorldNConvexBodiesInLinks(Tworld *w)
Gets the number of convex polyhedrons that define the mechanism included in the world.
Definition: world.c:1472
Definition of the Tparameters type and the associated functions.
double sign
Definition: world.h:82
Definition of the Tequations type and the associated functions.
void GetWorldSimpInitialBox(Tparameters *p, Tbox *b, Tworld *w)
Gets the kinematic simplified search space for a given problem.
Definition: world.c:2212
unsigned int nwcd
Definition: world.h:210
unsigned int GetWorldMobility(Tworld *w)
Returns the number of degrees of freedom of the mechanism in the world.
Definition: world.c:1427
void CheckLinkObstacleCollision(unsigned int a, unsigned int b, Tworld *w)
Activates the possible collision between a particular link and an object in the environment.
Definition: world.c:1703
unsigned int * jointTo
Definition: world.h:176
double maxCoord
Definition: world.h:142
unsigned int GetSolutionPointFromLinkTransforms(Tparameters *p, THTransform *tl, double **sol, Tworld *w)
Determines the mechanisms configuration from the pose of all links.
Definition: world.c:1215