robotoc
robotoc - efficient ROBOT Optimal Control solvers
Loading...
Searching...
No Matches
contact_status.hxx
Go to the documentation of this file.
1#ifndef ROBOTOC_CONTACT_STATUS_HXX_
2#define ROBOTOC_CONTACT_STATUS_HXX_
3
5
6#include <stdexcept>
7#include <cassert>
8#include <random>
9#include <chrono>
10
11
12namespace robotoc {
13
15 const std::vector<ContactType>& contact_types,
16 const std::vector<std::string>& contact_frame_names,
17 const double default_friction_coefficients)
18 : contact_types_(contact_types),
19 contact_frame_names_(contact_frame_names),
20 is_contact_active_(contact_types.size(), false),
21 contact_placements_(contact_types.size(), SE3::Identity()),
22 contact_positions_(contact_types.size(), Eigen::Vector3d::Zero()),
23 contact_rotations_(contact_types.size(), Eigen::Matrix3d::Identity()),
24 friction_coefficients_(contact_types.size(), 0.7),
25 dimf_(0),
26 max_contacts_(contact_types.size()),
27 max_num_contacts_(contact_types.size()),
28 has_active_contacts_(false) {
29 if (contact_types.size() != contact_frame_names.size()) {
30 throw std::invalid_argument(
31 "[ContactStatus] invalid argument: contact_frame_names.size() must be the same as contact_types.size()!");
32 }
33 if (default_friction_coefficients <= 0.0) {
34 throw std::invalid_argument("[ContactStatus] invalid argument: 'default_friction_coefficients' must be positive!");
35 }
36}
37
38
40 : contact_types_(),
41 contact_frame_names_(),
42 is_contact_active_(),
43 contact_placements_(),
44 contact_positions_(),
45 contact_rotations_(),
46 friction_coefficients_(),
47 dimf_(0),
48 max_contacts_(0),
49 max_num_contacts_(0),
50 has_active_contacts_(false) {
51}
52
53
54inline bool ContactStatus::operator==(const ContactStatus& other) const {
55 assert(other.maxNumContacts() == max_num_contacts_);
56 for (int i=0; i<max_num_contacts_; ++i) {
57 if (other.isContactActive(i) != isContactActive(i)) {
58 return false;
59 }
60 if (!other.contactPlacement(i).isApprox(contactPlacement(i))) {
61 return false;
62 }
63 }
64 return true;
65}
66
67
68inline bool ContactStatus::operator!=(const ContactStatus& other) const {
69 return !(*this == other);
70}
71
72
73inline ContactType ContactStatus::contactType(const int contact_index) const {
74 assert(contact_index >= 0);
75 assert(contact_index < max_num_contacts_);
76 return contact_types_[contact_index];
77}
78
79
80inline const std::vector<ContactType>& ContactStatus::contactTypes() const {
81 return contact_types_;
82}
83
84
85inline const std::string& ContactStatus::contactFrameName(
86 const int contact_index) const {
87 assert(contact_index >= 0);
88 assert(contact_index < max_num_contacts_);
89 if (contact_frame_names_.empty()) {
90 throw std::runtime_error("[ContactStatus] invalid argument: contact_frame_names_ is empty!");
91 }
92 return contact_frame_names_[contact_index];
93}
94
95
96inline const std::vector<std::string>& ContactStatus::contactFrameNames() const {
97 return contact_frame_names_;
98}
99
100
101inline bool ContactStatus::isContactActive(const int contact_index) const {
102 assert(!is_contact_active_.empty());
103 assert(contact_index >= 0);
104 assert(contact_index < is_contact_active_.size());
105 return is_contact_active_[contact_index];
106}
107
108
110 const std::string& contact_frame_name) const {
111 return isContactActive(findContactIndex(contact_frame_name));
112}
113
114
115inline const std::vector<bool>& ContactStatus::isContactActive() const {
116 return is_contact_active_;
117}
118
119
121 return has_active_contacts_;
122}
123
124
125inline int ContactStatus::dimf() const {
126 return dimf_;
127}
128
129
131 return max_num_contacts_;
132}
133
134
135inline void ContactStatus::activateContact(const int contact_index) {
136 assert(contact_index >= 0);
137 assert(contact_index < max_num_contacts_);
138 if (!is_contact_active_[contact_index]) {
139 is_contact_active_[contact_index] = true;
140 switch (contact_types_[contact_index]) {
142 dimf_ += 3;
143 break;
145 dimf_ += 6;
146 break;
147 default:
148 break;
149 }
150 }
151 setHasActiveContacts();
152}
153
154
155inline void ContactStatus::activateContact(const std::string& contact_frame_name) {
156 activateContact(findContactIndex(contact_frame_name));
157}
158
159
160inline void ContactStatus::deactivateContact(const int contact_index) {
161 assert(contact_index >= 0);
162 assert(contact_index < max_num_contacts_);
163 if (is_contact_active_[contact_index]) {
164 is_contact_active_[contact_index] = false;
165 switch (contact_types_[contact_index]) {
167 dimf_ -= 3;
168 break;
170 dimf_ -= 6;
171 break;
172 default:
173 break;
174 }
175 }
176 setHasActiveContacts();
177}
178
179
180inline void ContactStatus::deactivateContact(const std::string& contact_frame_name) {
181 deactivateContact(findContactIndex(contact_frame_name));
182}
183
184
186 const std::vector<int>& contact_indices) {
187 assert(contact_indices.size() <= max_num_contacts_);
188 for (const int e : contact_indices) {
190 }
191}
192
193
195 const std::vector<std::string>& contact_frame_names) {
196 assert(contact_frame_names.size() <= max_num_contacts_);
197 for (const auto& e : contact_frame_names) {
199 }
200}
201
202
204 const std::vector<int>& contact_indices) {
205 assert(contact_indices.size() <= max_num_contacts_);
206 for (const int e : contact_indices) {
208 }
209}
210
211
213 const std::vector<std::string>& contact_frame_names) {
214 assert(contact_frame_names.size() <= max_num_contacts_);
215 for (const auto& e : contact_frame_names) {
217 }
218}
219
220
222 const int contact_index, const Eigen::Vector3d& contact_position) {
223 setContactPlacement(contact_index, contact_position,
224 Eigen::Matrix3d::Identity());
225}
226
227
229 const std::string& contact_frame_name,
230 const Eigen::Vector3d& contact_position) {
231 setContactPlacement(findContactIndex(contact_frame_name), contact_position,
232 Eigen::Matrix3d::Identity());
233}
234
235
237 const int contact_index, const Eigen::Vector3d& contact_position,
238 const Eigen::Matrix3d& contact_rotation) {
239 assert(contact_index >= 0);
240 assert(contact_index < max_num_contacts_);
241 contact_positions_[contact_index] = contact_position;
242 contact_rotations_[contact_index] = contact_rotation;
243 contact_placements_[contact_index] = SE3(contact_rotation, contact_position);
244}
245
246
248 const std::string& contact_frame_name,
249 const Eigen::Vector3d& contact_position,
250 const Eigen::Matrix3d& contact_rotation) {
251 setContactPlacement(findContactIndex(contact_frame_name), contact_position,
252 contact_rotation);
253}
254
255
256inline void ContactStatus::setContactPlacement(const int contact_index,
257 const SE3& contact_placement) {
258 assert(contact_index >= 0);
259 assert(contact_index < max_num_contacts_);
260 contact_positions_[contact_index] = contact_placement.translation();
261 contact_rotations_[contact_index] = contact_placement.rotation();
262 contact_placements_[contact_index] = contact_placement;
263}
264
265
266inline void ContactStatus::setContactPlacement(const std::string& contact_frame_name,
267 const SE3& contact_placement) {
268 setContactPlacement(contact_frame_name, contact_placement);
269}
270
271
273 const std::vector<Eigen::Vector3d>& contact_positions) {
274 assert(contact_positions.size() == max_num_contacts_);
275 for (int i=0; i<max_num_contacts_; ++i) {
276 setContactPlacement(i, contact_positions[i], Eigen::Matrix3d::Identity());
277 }
278}
279
280
282 const std::unordered_map<std::string, Eigen::Vector3d>& contact_positions) {
283 assert(contact_positions.size() == max_num_contacts_);
284 for (int i=0; i<max_num_contacts_; ++i) {
285 setContactPlacement(i, contact_positions.at(contactFrameName(i)),
286 Eigen::Matrix3d::Identity());
287 }
288}
289
290
292 const std::vector<Eigen::Vector3d>& contact_positions,
293 const std::vector<Eigen::Matrix3d>& contact_rotations) {
294 assert(contact_positions.size() == max_num_contacts_);
295 assert(contact_rotations.size() == max_num_contacts_);
296 for (int i=0; i<max_num_contacts_; ++i) {
297 setContactPlacement(i, contact_positions[i], contact_rotations[i]);
298 }
299}
300
301
303 const std::unordered_map<std::string, Eigen::Vector3d>& contact_positions,
304 const std::unordered_map<std::string, Eigen::Matrix3d>& contact_rotations) {
305 assert(contact_positions.size() == max_num_contacts_);
306 assert(contact_rotations.size() == max_num_contacts_);
307 for (int i=0; i<max_num_contacts_; ++i) {
308 setContactPlacement(i, contact_positions.at(contactFrameName(i)),
309 contact_rotations.at(contactFrameName(i)));
310 }
311}
312
313
315 const aligned_vector<SE3>& contact_placements) {
316 assert(contact_placements.size() == max_num_contacts_);
317 for (int i=0; i<max_num_contacts_; ++i) {
318 setContactPlacement(i, contact_placements[i]);
319 }
320}
321
322
324 const aligned_unordered_map<std::string, SE3>& contact_placements) {
325 assert(contact_placements.size() == max_num_contacts_);
326 for (int i=0; i<max_num_contacts_; ++i) {
327 setContactPlacement(i, contact_placements.at(contactFrameName(i)));
328 }
329}
330
331
333 const int contact_index) const {
334 assert(contact_index >= 0);
335 assert(contact_index < max_num_contacts_);
336 return contact_placements_[contact_index];
337}
338
339
341 const std::string& contact_frame_name) const {
342 return contactPlacement(findContactIndex(contact_frame_name));
343}
344
345
346inline const Eigen::Vector3d& ContactStatus::contactPosition(
347 const int contact_index) const {
348 assert(contact_index >= 0);
349 assert(contact_index < max_num_contacts_);
350 return contact_positions_[contact_index];
351}
352
353
354inline const Eigen::Vector3d& ContactStatus::contactPosition(
355 const std::string& contact_frame_name) const {
356 return contactPosition(findContactIndex(contact_frame_name));
357}
358
359
360inline const Eigen::Matrix3d& ContactStatus::contactRotation(
361 const int contact_index) const {
362 assert(contact_index >= 0);
363 assert(contact_index < max_num_contacts_);
364 return contact_rotations_[contact_index];
365}
366
367
368inline const Eigen::Matrix3d& ContactStatus::contactRotation(
369 const std::string& contact_frame_name) const {
370 return contactRotation(findContactIndex(contact_frame_name));
371}
372
373
375 return contact_placements_;
376}
377
378
379inline const std::vector<Eigen::Vector3d>&
381 return contact_positions_;
382}
383
384
385inline const std::vector<Eigen::Matrix3d>&
387 return contact_rotations_;
388}
389
390
391inline void ContactStatus::setFrictionCoefficient(const int contact_index,
392 const double friction_coefficient) {
393 if (friction_coefficient <= 0.0) {
394 throw std::invalid_argument("[ContactStatus] invalid argument: 'friction_coefficient' must be positive!");
395 }
396 assert(contact_index >= 0);
397 assert(contact_index < max_num_contacts_);
398 friction_coefficients_[contact_index] = friction_coefficient;
399}
400
401
402inline void ContactStatus::setFrictionCoefficient(const std::string& contact_frame_name,
403 const double friction_coefficient) {
404 setFrictionCoefficient(findContactIndex(contact_frame_name), friction_coefficient);
405}
406
407
409 const std::vector<double>& friction_coefficients) {
410 assert(friction_coefficients.size() == max_num_contacts_);
411 for (int i=0; i<max_num_contacts_; ++i) {
412 setFrictionCoefficient(i, friction_coefficients[i]);
413 }
414}
415
416
418 const std::unordered_map<std::string, double>& friction_coefficients) {
419 assert(friction_coefficients.size() == max_num_contacts_);
420 for (int i=0; i<max_num_contacts_; ++i) {
421 setFrictionCoefficient(i, friction_coefficients.at(contactFrameName(i)));
422 }
423}
424
425
426inline double ContactStatus::frictionCoefficient(const int contact_index) const {
427 assert(contact_index >= 0);
428 assert(contact_index < max_num_contacts_);
429 return friction_coefficients_[contact_index];
430}
431
432
434 const std::string& contact_frame_name) const {
435 return frictionCoefficient(findContactIndex(contact_frame_name));
436}
437
438
439inline const std::vector<double>& ContactStatus::frictionCoefficients() const {
440 return friction_coefficients_;
441}
442
443
445 const std::string& contact_frame_name) const {
446 if (contact_frame_names_.empty()) {
447 throw std::runtime_error("[ContactStatus] invalid argument: contact_frame_names_ is empty!");
448 }
449 for (int i=0; i<contact_frame_names_.size(); ++i) {
450 if (contact_frame_names_[i] == contact_frame_name) {
451 return i;
452 }
453 }
454 throw std::runtime_error("[ContactStatus] cannot find the input contact_frame_name: '" + contact_frame_name + "'");
455}
456
457
459 std::minstd_rand0 rand(
460 std::chrono::system_clock::now().time_since_epoch().count());
461 for (int i=0; i<max_num_contacts_; ++i) {
462 if (rand()%2 == 0) {
464 }
465 else {
467 }
468 }
469}
470
471
472inline void ContactStatus::setHasActiveContacts() {
473 has_active_contacts_ = (dimf_ > 0);
474}
475
476
477} // namespace robotoc
478
479#endif // ROBOTOC_CONTACT_STATUS_HXX_
Contact status of robot model.
Definition: contact_status.hpp:32
const std::string & contactFrameName(const int contact_index) const
Returns the name of the contact frame.
Definition: contact_status.hxx:85
const std::vector< Eigen::Matrix3d > & contactRotations() const
Gets the contact rotations.
Definition: contact_status.hxx:386
ContactStatus()
Default constructor.
Definition: contact_status.hxx:39
int findContactIndex(const std::string &contact_frame_name) const
Finds the contact index correspoinding to the input contact frame name.
Definition: contact_status.hxx:444
void deactivateContacts(const std::vector< int > &contact_indices)
Deactivates contacts.
Definition: contact_status.hxx:203
const aligned_vector< SE3 > & contactPlacements() const
Gets the contact placements.
Definition: contact_status.hxx:374
bool operator==(const ContactStatus &other) const
Defines a comparison operator.
Definition: contact_status.hxx:54
const Eigen::Matrix3d & contactRotation(const int contact_index) const
Gets the contact rotation.
Definition: contact_status.hxx:360
double frictionCoefficient(const int contact_index) const
Gets the friction coefficint. Default value is 0.7.
Definition: contact_status.hxx:426
bool isContactActive(const int contact_index) const
Returns true if a contact is active and false if not.
Definition: contact_status.hxx:101
ContactType contactType(const int contact_index) const
Returns the type of the contact.
Definition: contact_status.hxx:73
int maxNumContacts() const
Returns the maximum number of the contacts.
Definition: contact_status.hxx:130
int dimf() const
Returns the dimension of the active contact forces.
Definition: contact_status.hxx:125
void setRandom()
Fills contact status randomly.
Definition: contact_status.hxx:458
void deactivateContact(const int contact_index)
Deactivates a contact.
Definition: contact_status.hxx:160
void activateContacts(const std::vector< int > &contact_indices)
Activates contacts.
Definition: contact_status.hxx:185
void setFrictionCoefficient(const int contact_index, const double friction_coefficient)
Gets the friction coefficint.
Definition: contact_status.hxx:391
void setContactPlacement(const int contact_index, const Eigen::Vector3d &contact_position)
Sets a contact placement, that is, the position and rotation of the contact. The contact rotation is ...
Definition: contact_status.hxx:221
const std::vector< double > & frictionCoefficients() const
Gets the friction coefficints. Default value is 0.7.
Definition: contact_status.hxx:439
const std::vector< bool > & isContactActive() const
Returns the activity of the contacts.
Definition: contact_status.hxx:115
void setContactPlacements(const std::vector< Eigen::Vector3d > &contact_positions)
Sets contact placements. The rotation of each contact is set to Eigen::Matrix3d::Identity(),...
Definition: contact_status.hxx:272
const SE3 & contactPlacement(const int contact_index) const
Gets the contact placement.
Definition: contact_status.hxx:332
void setFrictionCoefficients(const std::vector< double > &friction_coefficients)
Sets the friction coefficints.
Definition: contact_status.hxx:408
const Eigen::Vector3d & contactPosition(const int contact_index) const
Gets the contact position.
Definition: contact_status.hxx:346
const std::vector< Eigen::Vector3d > & contactPositions() const
Gets the contact positions.
Definition: contact_status.hxx:380
bool hasActiveContacts() const
Returns true if there are active contacts and false if not.
Definition: contact_status.hxx:120
void activateContact(const int contact_index)
Activates a contact.
Definition: contact_status.hxx:135
const std::vector< std::string > & contactFrameNames() const
Returns the names of the contact frames.
Definition: contact_status.hxx:96
const std::vector< ContactType > & contactTypes() const
Returns the types of the contacts.
Definition: contact_status.hxx:80
bool operator!=(const ContactStatus &other) const
Defines a comparison operator.
Definition: contact_status.hxx:68
Definition: constraint_component_base.hpp:17
ContactType
Types of contacts.
Definition: contact_status.hpp:22
std::vector< T, Eigen::aligned_allocator< T > > aligned_vector
std vector with Eigen::aligned_allocator.
Definition: aligned_vector.hpp:14
std::unordered_map< Key, T, std::hash< Key >, std::equal_to< Key >, Eigen::aligned_allocator< std::pair< Key, T > > > aligned_unordered_map
std unordered_map with Eigen::aligned_allocator.
Definition: aligned_unordered_map.hpp:15
pinocchio::SE3 SE3
Using pinocchio::SE3 without its namespace.
Definition: se3.hpp:15