occ
Loading...
Searching...
No Matches
occ::mults::RigidBodyState Struct Reference

State of a rigid body for molecular dynamics simulations. More...

#include <rigid_body.h>

Collaboration diagram for occ::mults::RigidBodyState:
[legend]

Classes

struct  BodySite
 Body-frame site: multipole expansion at an offset from the center of mass. More...
 

Public Member Functions

 RigidBodyState ()
 Default constructor.
 
 RigidBodyState (const std::string &mol_name, int mol_id, const Vec3 &pos, double mol_mass, int max_rank=4)
 Construct with basic parameters.
 
Mat3 rotation_matrix () const
 Get rotation matrix from current quaternion.
 
void set_euler_angles (double alpha, double beta, double gamma)
 Set orientation from Euler angles (ZYZ convention)
 
Vec3 get_euler_angles ()
 Get Euler angles from current quaternion (cached if possible)
 
void set_quaternion (const Eigen::Quaterniond &q)
 Set orientation from quaternion.
 
void set_rotation_matrix (const Mat3 &R)
 Set orientation from rotation matrix.
 
Vec3 get_angle_axis () const
 Get orientation as angle-axis vector.
 
void set_angle_axis (const Vec3 &p)
 Set orientation from angle-axis vector.
 
Mat3 angle_axis_jacobian () const
 Compute Jacobian relating angle-axis to Euler angle derivatives.
 
std::array< Mat3, 3 > rotation_matrix_derivatives () const
 Compute rotation matrix derivatives with respect to angle-axis components.
 
std::array< Mat3, 9 > rotation_matrix_second_derivatives () const
 Compute second derivatives of rotation matrix w.r.t.
 
void update_lab_multipoles ()
 Update lab-frame multipoles from body-frame multipoles.
 
const occ::dma::Multget_lab_multipoles ()
 Get multipoles in lab frame (updating if necessary)
 
void invalidate_lab_multipoles ()
 Mark lab multipoles as invalid (call after orientation change)
 
bool is_multi_site () const
 Check if this molecule uses multi-site representation.
 
int num_sites () const
 Get number of multipole sites.
 
void set_single_site (const occ::dma::Mult &mult)
 Set single-site multipole (clears any multi-site data)
 
void set_sites (const std::vector< BodySite > &sites)
 Set multi-site data (replaces single-site multipole_body)
 
void set_inertia_tensor (const Mat3 &I)
 Set inertia tensor (and compute inverse)
 
void set_spherical_inertia (double I_moment)
 Set isotropic (spherical) inertia.
 
void set_diagonal_inertia (double Ixx, double Iyy, double Izz)
 Set diagonal inertia tensor.
 
Mat3 inertia_space () const
 Get inertia tensor in space frame (rotated from body frame)
 
double translational_kinetic_energy () const
 Compute translational kinetic energy.
 
double rotational_kinetic_energy () const
 Compute rotational kinetic energy.
 
double kinetic_energy () const
 Compute total kinetic energy.
 
Vec3 angular_momentum_body () const
 Compute angular momentum in body frame.
 
Vec3 angular_momentum_space () const
 Compute angular momentum in space frame.
 

Public Attributes

std::string name
 
int id
 
Vec3 position
 
Vec3 velocity
 
Vec3 force
 
double mass
 
Eigen::Quaterniond quaternion
 
Vec3 angular_velocity_body
 
Vec3 torque_body
 
Vec3 euler_angles
 
Vec3 torque_euler
 
bool euler_valid
 
Mat3 inertia_body
 
Mat3 inertia_inv_body
 
occ::dma::Mult multipole_body
 
occ::dma::Mult multipole_lab
 
bool multipole_lab_valid
 
std::vector< BodySitesites_body
 

Detailed Description

State of a rigid body for molecular dynamics simulations.

This structure contains all kinematic and dynamic state information for a rigid molecule, including position, orientation, linear/angular velocities, forces, torques, and multipole moments.

Orientation Representations

Multiple representations are supported:

  • Quaternion (primary storage): Used for MD integration, numerically stable
  • Angle-axis (optimization): Minimal parametrization, no singularities
  • Euler angles (I/O): User-friendly input/output (ZYZ convention)

When to Use Each:

  • Quaternions: MD time integration, internal storage
  • Angle-axis: Geometry optimization, gradient-based minimization
  • Euler angles: Reading/writing input files, user interface

Conversion Example:

mol.set_euler_angles(0.5, 0.3, 0.2); // User input
Vec3 p = mol.get_angle_axis(); // For optimization
// ... optimize ...
mol.set_angle_axis(p_optimized); // Update state
Vec3 euler = mol.get_euler_angles(); // For output
Eigen::Vector3d Vec3
Definition linear_algebra.h:64
State of a rigid body for molecular dynamics simulations.
Definition rigid_body.h:39
void set_euler_angles(double alpha, double beta, double gamma)
Set orientation from Euler angles (ZYZ convention)
void set_angle_axis(const Vec3 &p)
Set orientation from angle-axis vector.
Vec3 get_angle_axis() const
Get orientation as angle-axis vector.
Vec3 get_euler_angles()
Get Euler angles from current quaternion (cached if possible)

Constructor & Destructor Documentation

◆ RigidBodyState() [1/2]

occ::mults::RigidBodyState::RigidBodyState ( )
inline

Default constructor.

◆ RigidBodyState() [2/2]

occ::mults::RigidBodyState::RigidBodyState ( const std::string &  mol_name,
int  mol_id,
const Vec3 pos,
double  mol_mass,
int  max_rank = 4 
)
inline

Construct with basic parameters.

Member Function Documentation

◆ angle_axis_jacobian()

Mat3 occ::mults::RigidBodyState::angle_axis_jacobian ( ) const

Compute Jacobian relating angle-axis to Euler angle derivatives.

For gradient transformation in optimization: grad_aa = J^T * grad_euler where grad_euler = [dE/dalpha, dE/dbeta, dE/dgamma]

Implements Orient's aaderiv formula from interact.f90

Returns
3x3 Jacobian matrix

◆ angular_momentum_body()

Vec3 occ::mults::RigidBodyState::angular_momentum_body ( ) const

Compute angular momentum in body frame.

◆ angular_momentum_space()

Vec3 occ::mults::RigidBodyState::angular_momentum_space ( ) const

Compute angular momentum in space frame.

◆ get_angle_axis()

Vec3 occ::mults::RigidBodyState::get_angle_axis ( ) const

Get orientation as angle-axis vector.

Returns p = psi * n where:

  • psi = rotation angle in [0, pi]
  • n = unit rotation axis

This is the parametrization Orient uses for optimization.

Returns
Angle-axis vector (3 components)

◆ get_euler_angles()

Vec3 occ::mults::RigidBodyState::get_euler_angles ( )

Get Euler angles from current quaternion (cached if possible)

◆ get_lab_multipoles()

const occ::dma::Mult & occ::mults::RigidBodyState::get_lab_multipoles ( )

Get multipoles in lab frame (updating if necessary)

◆ inertia_space()

Mat3 occ::mults::RigidBodyState::inertia_space ( ) const

Get inertia tensor in space frame (rotated from body frame)

◆ invalidate_lab_multipoles()

void occ::mults::RigidBodyState::invalidate_lab_multipoles ( )
inline

Mark lab multipoles as invalid (call after orientation change)

◆ is_multi_site()

bool occ::mults::RigidBodyState::is_multi_site ( ) const
inline

Check if this molecule uses multi-site representation.

◆ kinetic_energy()

double occ::mults::RigidBodyState::kinetic_energy ( ) const

Compute total kinetic energy.

◆ num_sites()

int occ::mults::RigidBodyState::num_sites ( ) const
inline

Get number of multipole sites.

◆ rotation_matrix()

Mat3 occ::mults::RigidBodyState::rotation_matrix ( ) const

Get rotation matrix from current quaternion.

◆ rotation_matrix_derivatives()

std::array< Mat3, 3 > occ::mults::RigidBodyState::rotation_matrix_derivatives ( ) const

Compute rotation matrix derivatives with respect to angle-axis components.

Computes ∂M/∂p_i for i=1,2,3 where M is the rotation matrix and p is the angle-axis vector. This implements Orient's aaderiv routine.

The derivatives are computed analytically using the Rodrigues formula: M = I + sin(ψ)·[n]× + (1-cos(ψ))·[n]ײ

where ψ = |p|, n = p/ψ, and [n]× is the skew-symmetric (cross-product) matrix.

For small rotations (ψ < 1e-8), uses first-order approximation. For general rotations, uses exact derivatives via chain rule.

Returns
Array of 3 matrices, where result[i] = ∂M/∂p_i

◆ rotation_matrix_second_derivatives()

std::array< Mat3, 9 > occ::mults::RigidBodyState::rotation_matrix_second_derivatives ( ) const

Compute second derivatives of rotation matrix w.r.t.

angle-axis components

Computes ∂²M/∂p_k∂p_l for k,l = 0,1,2 (9 matrices total). The result is stored in row-major order: result[3*k + l] = ∂²M/∂p_k∂p_l

For small rotations (ψ < 1e-8), returns zeros (second-order term). For general rotations, computed analytically from Rodrigues formula.

Returns
Array of 9 matrices, where result[3*k + l] = ∂²M/∂p_k∂p_l

◆ rotational_kinetic_energy()

double occ::mults::RigidBodyState::rotational_kinetic_energy ( ) const

Compute rotational kinetic energy.

KE_rot = (1/2) * omega^T * I_body * omega

◆ set_angle_axis()

void occ::mults::RigidBodyState::set_angle_axis ( const Vec3 p)

Set orientation from angle-axis vector.

Parameters
pAngle-axis vector where |p| = rotation angle

◆ set_diagonal_inertia()

void occ::mults::RigidBodyState::set_diagonal_inertia ( double  Ixx,
double  Iyy,
double  Izz 
)

Set diagonal inertia tensor.

◆ set_euler_angles()

void occ::mults::RigidBodyState::set_euler_angles ( double  alpha,
double  beta,
double  gamma 
)

Set orientation from Euler angles (ZYZ convention)

◆ set_inertia_tensor()

void occ::mults::RigidBodyState::set_inertia_tensor ( const Mat3 I)

Set inertia tensor (and compute inverse)

◆ set_quaternion()

void occ::mults::RigidBodyState::set_quaternion ( const Eigen::Quaterniond &  q)

Set orientation from quaternion.

◆ set_rotation_matrix()

void occ::mults::RigidBodyState::set_rotation_matrix ( const Mat3 R)

Set orientation from rotation matrix.

◆ set_single_site()

void occ::mults::RigidBodyState::set_single_site ( const occ::dma::Mult mult)
inline

Set single-site multipole (clears any multi-site data)

◆ set_sites()

void occ::mults::RigidBodyState::set_sites ( const std::vector< BodySite > &  sites)
inline

Set multi-site data (replaces single-site multipole_body)

◆ set_spherical_inertia()

void occ::mults::RigidBodyState::set_spherical_inertia ( double  I_moment)

Set isotropic (spherical) inertia.

◆ translational_kinetic_energy()

double occ::mults::RigidBodyState::translational_kinetic_energy ( ) const

Compute translational kinetic energy.

◆ update_lab_multipoles()

void occ::mults::RigidBodyState::update_lab_multipoles ( )

Update lab-frame multipoles from body-frame multipoles.

This rotates the body-frame multipoles using the current orientation and caches the result in multipole_lab.

Member Data Documentation

◆ angular_velocity_body

Vec3 occ::mults::RigidBodyState::angular_velocity_body

◆ euler_angles

Vec3 occ::mults::RigidBodyState::euler_angles

◆ euler_valid

bool occ::mults::RigidBodyState::euler_valid

◆ force

Vec3 occ::mults::RigidBodyState::force

◆ id

int occ::mults::RigidBodyState::id

◆ inertia_body

Mat3 occ::mults::RigidBodyState::inertia_body

◆ inertia_inv_body

Mat3 occ::mults::RigidBodyState::inertia_inv_body

◆ mass

double occ::mults::RigidBodyState::mass

◆ multipole_body

occ::dma::Mult occ::mults::RigidBodyState::multipole_body

◆ multipole_lab

occ::dma::Mult occ::mults::RigidBodyState::multipole_lab

◆ multipole_lab_valid

bool occ::mults::RigidBodyState::multipole_lab_valid

◆ name

std::string occ::mults::RigidBodyState::name

◆ position

Vec3 occ::mults::RigidBodyState::position

◆ quaternion

Eigen::Quaterniond occ::mults::RigidBodyState::quaternion

◆ sites_body

std::vector<BodySite> occ::mults::RigidBodyState::sites_body

◆ torque_body

Vec3 occ::mults::RigidBodyState::torque_body

◆ torque_euler

Vec3 occ::mults::RigidBodyState::torque_euler

◆ velocity

Vec3 occ::mults::RigidBodyState::velocity

The documentation for this struct was generated from the following file: