Skip to content

Commit 2172383

Browse files
[#1392] Split qpos/qvel into translation/attitude StateData
1 parent c81eb27 commit 2172383

6 files changed

Lines changed: 596 additions & 70 deletions

File tree

src/simulation/mujocoDynamics/_GeneralModuleFiles/MJJoint.cpp

Lines changed: 33 additions & 22 deletions
Original file line numberDiff line numberDiff line change
@@ -102,31 +102,33 @@ void MJScalarJoint::writeJointStateMessage(uint64_t CurrentSimNanos)
102102
{
103103
checkInitialized();
104104

105+
auto& scene = body.getSpec().getScene();
106+
105107
ScalarJointStateMsgPayload stateOutMsgPayload;
106-
auto state = body.getSpec().getScene().getQposState()->getState();
107-
stateOutMsgPayload.state = state(this->qposAdr.value());
108-
this->stateOutMsg.write(&stateOutMsgPayload, body.getSpec().getScene().moduleID, 0);
108+
stateOutMsgPayload.state = scene.qposGet(this->qposAdr.value());
109+
this->stateOutMsg.write(&stateOutMsgPayload, scene.moduleID, 0);
109110

110111
ScalarJointStateMsgPayload stateDotOutMsgPayload;
111-
auto stateDot = body.getSpec().getScene().getQvelState()->getState();
112-
stateDotOutMsgPayload.state = stateDot(this->qvelAdr.value());
113-
this->stateDotOutMsg.write(&stateDotOutMsgPayload, body.getSpec().getScene().moduleID, 0);
112+
stateDotOutMsgPayload.state = scene.qvelGet(this->qvelAdr.value());
113+
this->stateDotOutMsg.write(&stateDotOutMsgPayload, scene.moduleID, 0);
114114
}
115115

116116
void MJScalarJoint::setPosition(double value)
117117
{
118118
checkInitialized();
119119

120-
body.getSpec().getScene().getQposState()->state(this->qposAdr.value()) = value;
121-
body.getSpec().getScene().markKinematicsAsStale();
120+
auto& scene = body.getSpec().getScene();
121+
scene.qposRef(this->qposAdr.value()) = value;
122+
scene.markKinematicsAsStale();
122123
}
123124

124125
void MJScalarJoint::setVelocity(double value)
125126
{
126127
checkInitialized();
127128

128-
body.getSpec().getScene().getQvelState()->state(this->qvelAdr.value()) = value;
129-
body.getSpec().getScene().markKinematicsAsStale();
129+
auto& scene = body.getSpec().getScene();
130+
scene.qvelRef(this->qvelAdr.value()) = value;
131+
scene.markKinematicsAsStale();
130132
}
131133

132134
MJSingleJointEquality
@@ -139,45 +141,54 @@ void MJFreeJoint::setPosition(const Eigen::Vector3d& position)
139141
{
140142
checkInitialized();
141143

142-
auto& qposState = body.getSpec().getScene().getQposState()->state;
144+
auto& scene = body.getSpec().getScene();
143145
auto i = this->qposAdr.value();
144-
qposState.middleRows(i, 3) = position;
146+
scene.qposRef(i + 0) = position[0];
147+
scene.qposRef(i + 1) = position[1];
148+
scene.qposRef(i + 2) = position[2];
145149

146-
body.getSpec().getScene().markKinematicsAsStale();
150+
scene.markKinematicsAsStale();
147151
}
148152

149153
void MJFreeJoint::setVelocity(const Eigen::Vector3d& velocity)
150154
{
151155
checkInitialized();
152156

153-
auto& qvelState = body.getSpec().getScene().getQvelState()->state;
157+
auto& scene = body.getSpec().getScene();
154158
auto i = this->qvelAdr.value();
155-
qvelState.middleRows(i, 3) = velocity;
159+
scene.qvelRef(i + 0) = velocity[0];
160+
scene.qvelRef(i + 1) = velocity[1];
161+
scene.qvelRef(i + 2) = velocity[2];
156162

157-
body.getSpec().getScene().markKinematicsAsStale();
163+
scene.markKinematicsAsStale();
158164
}
159165

160166
void MJFreeJoint::setAttitude(const Eigen::MRPd& attitude)
161167
{
162168
checkInitialized();
163169

164-
auto& qposState = body.getSpec().getScene().getQposState()->state;
170+
auto& scene = body.getSpec().getScene();
165171
auto i = this->qposAdr.value();
166172

167173
auto mat = attitude.toRotationMatrix();
168174
auto quat = Eigen::Quaterniond(mat);
169175

170-
qposState.middleRows(i + 3, 4) = Eigen::Vector4d{quat.w(), quat.x(), quat.y(), quat.z()};
171-
body.getSpec().getScene().markKinematicsAsStale();
176+
scene.qposRef(i + 3) = quat.w();
177+
scene.qposRef(i + 4) = quat.x();
178+
scene.qposRef(i + 5) = quat.y();
179+
scene.qposRef(i + 6) = quat.z();
180+
scene.markKinematicsAsStale();
172181
}
173182

174183
void MJFreeJoint::setAttitudeRate(const Eigen::Vector3d& attitudeRate)
175184
{
176185
checkInitialized();
177186

178-
auto& qvelState = body.getSpec().getScene().getQvelState()->state;
187+
auto& scene = body.getSpec().getScene();
179188
auto i = this->qvelAdr.value();
180-
qvelState.middleRows(i + 3, 3) = attitudeRate;
189+
scene.qvelRef(i + 3) = attitudeRate[0];
190+
scene.qvelRef(i + 4) = attitudeRate[1];
191+
scene.qvelRef(i + 5) = attitudeRate[2];
181192

182-
body.getSpec().getScene().markKinematicsAsStale();
193+
scene.markKinematicsAsStale();
183194
}
Lines changed: 68 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,68 @@
1+
/*
2+
ISC License
3+
4+
Copyright (c) 2026, Autonomous Vehicle Systems Lab, University of Colorado at Boulder
5+
6+
Permission to use, copy, modify, and/or distribute this software for any
7+
purpose with or without fee is hereby granted, provided that the above
8+
copyright notice and this permission notice appear in all copies.
9+
10+
THE SOFTWARE IS PROVIDED "AS IS" AND THE AUTHOR DISCLAIMS ALL WARRANTIES
11+
WITH REGARD TO THIS SOFTWARE INCLUDING ALL IMPLIED WARRANTIES OF
12+
MERCHANTABILITY AND FITNESS. IN NO EVENT SHALL THE AUTHOR BE LIABLE FOR
13+
ANY SPECIAL, DIRECT, INDIRECT, OR CONSEQUENTIAL DAMAGES OR ANY DAMAGES
14+
WHATSOEVER RESULTING FROM LOSS OF USE, DATA OR PROFITS, WHETHER IN AN
15+
ACTION OF CONTRACT, NEGLIGENCE OR OTHER TORTIOUS ACTION, ARISING OUT OF
16+
OR IN CONNECTION WITH THE USE OR PERFORMANCE OF THIS SOFTWARE.
17+
18+
*/
19+
20+
#include "MJQPosAttitudeStateData.h"
21+
22+
#include <mujoco/mujoco.h>
23+
24+
#include <stdexcept>
25+
26+
std::unique_ptr<StateData> MJQPosAttitudeStateData::clone() const
27+
{
28+
return std::unique_ptr<StateData>();
29+
}
30+
31+
void MJQPosAttitudeStateData::configure(int numQuaternions, int numLinearDofs)
32+
{
33+
this->numQuaternions = numQuaternions;
34+
this->numLinearDofs = numLinearDofs;
35+
36+
int stateSize = 4 * numQuaternions + numLinearDofs;
37+
int derivSize = 3 * numQuaternions + numLinearDofs;
38+
39+
this->state.resize(stateSize, 1);
40+
this->stateDeriv.resize(derivSize, 1);
41+
}
42+
43+
void MJQPosAttitudeStateData::propagateState(double dt, std::vector<double> pseudoStep)
44+
{
45+
// Quaternion blocks: state has 4 elements per quaternion, stateDeriv
46+
// has 3 elements per angular velocity.
47+
for (int i = 0; i < numQuaternions; ++i) {
48+
mju_quatIntegrate(this->state.data() + 4 * i,
49+
this->stateDeriv.data() + 3 * i,
50+
dt);
51+
}
52+
53+
// Linear DOFs follow the quaternion block in both vectors.
54+
int linearStart_state = 4 * numQuaternions;
55+
int linearStart_deriv = 3 * numQuaternions;
56+
for (int i = 0; i < numLinearDofs; ++i) {
57+
this->state(linearStart_state + i) +=
58+
this->stateDeriv(linearStart_deriv + i) * dt;
59+
}
60+
61+
if (getNumNoiseSources() > 0 && pseudoStep.size() == 0) {
62+
auto errorMsg = "State " + this->getName() + " has stochastic dynamics, but "
63+
+ "the integrator tried to propagate it without pseudoSteps. Are you sure "
64+
+ "you are using a stochastic integrator?";
65+
bskLogger.bskLog(BSK_ERROR, errorMsg.c_str());
66+
throw std::invalid_argument(errorMsg);
67+
}
68+
}
Lines changed: 84 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,84 @@
1+
/*
2+
ISC License
3+
4+
Copyright (c) 2026, Autonomous Vehicle Systems Lab, University of Colorado at Boulder
5+
6+
Permission to use, copy, modify, and/or distribute this software for any
7+
purpose with or without fee is hereby granted, provided that the above
8+
copyright notice and this permission notice appear in all copies.
9+
10+
THE SOFTWARE IS PROVIDED "AS IS" AND THE AUTHOR DISCLAIMS ALL WARRANTIES
11+
WITH REGARD TO THIS SOFTWARE INCLUDING ALL IMPLIED WARRANTIES OF
12+
MERCHANTABILITY AND FITNESS. IN NO EVENT SHALL THE AUTHOR BE LIABLE FOR
13+
ANY SPECIAL, DIRECT, INDIRECT, OR CONSEQUENTIAL DAMAGES OR ANY DAMAGES
14+
WHATSOEVER RESULTING FROM LOSS OF USE, DATA OR PROFITS, WHETHER IN AN
15+
ACTION OF CONTRACT, NEGLIGENCE OR OTHER TORTIOUS ACTION, ARISING OUT OF
16+
OR IN CONNECTION WITH THE USE OR PERFORMANCE OF THIS SOFTWARE.
17+
18+
*/
19+
20+
#ifndef MJQPOS_ATTITUDE_STATE_DATA_H
21+
#define MJQPOS_ATTITUDE_STATE_DATA_H
22+
23+
#include "simulation/dynamics/_GeneralModuleFiles/stateData.h"
24+
25+
/**
26+
* @brief Position state data for the attitude/hinges portion of an `MJScene`.
27+
*
28+
* Holds the rotational and scalar-joint position DOFs of an MJScene's qpos
29+
* vector, segregated from the free-joint translational DOFs (which live in a
30+
* separate plain `StateData`). The split lets the adaptive integrator scale
31+
* tolerances independently for translation (orbital scales) and attitude
32+
* (radian-scale): without it, the orbital state norm would dominate the
33+
* tolerance computation and either force absurdly small steps everywhere or
34+
* (with looser relTol) allow steps too large for the stiff panel-joint ODE.
35+
*
36+
* Layout convention (assembled by the owning `MJScene`):
37+
*
38+
* * `state` [4 * numQuaternions + numLinearDofs]
39+
* - First `4 * numQuaternions` entries: 4-element unit quaternions, one
40+
* per free-joint rotation or ball joint.
41+
* - Remaining `numLinearDofs` entries: scalar joint positions (hinges
42+
* and slides), one entry each.
43+
* * `stateDeriv` [3 * numQuaternions + numLinearDofs]
44+
* - First `3 * numQuaternions` entries: 3-element angular velocities
45+
* (one per quaternion block, in the same order).
46+
* - Remaining `numLinearDofs` entries: scalar joint velocities.
47+
*
48+
* `propagateState` integrates the quaternion blocks with `mju_quatIntegrate`
49+
* (proper SO(3) update with renormalization) and the linear DOFs with the
50+
* default Euler step `state += stateDeriv * dt`.
51+
*/
52+
class MJQPosAttitudeStateData : public StateData
53+
{
54+
public:
55+
/**
56+
* @brief Constructs an `MJQPosAttitudeStateData` with the given name and
57+
* initial state matrix (typically empty until `configure` runs).
58+
*/
59+
MJQPosAttitudeStateData(std::string inName, const Eigen::MatrixXd& newState)
60+
: StateData(std::move(inName), newState) {}
61+
62+
/** Clone constructor (returns null to match the existing MuJoCo state convention). */
63+
std::unique_ptr<StateData> clone() const override;
64+
65+
/**
66+
* @brief Sizes the `state` and `stateDeriv` vectors and records the layout.
67+
*
68+
* @param numQuaternions Number of 4-element quaternion blocks (free + ball joints).
69+
* @param numLinearDofs Number of 1-element scalar-joint DOFs (hinges + slides).
70+
*/
71+
void configure(int numQuaternions, int numLinearDofs);
72+
73+
/**
74+
* @brief Propagates the state by `dt`. Quaternions are updated with
75+
* `mju_quatIntegrate` and linear DOFs with `state += stateDeriv * dt`.
76+
*/
77+
void propagateState(double dt, std::vector<double> pseudoStep = {}) override;
78+
79+
protected:
80+
int numQuaternions = 0;
81+
int numLinearDofs = 0;
82+
};
83+
84+
#endif

0 commit comments

Comments
 (0)