Class PositionBasedElasticRods

Class Documentation

class PositionBasedElasticRods

Public Static Functions

static bool solve_PerpendiculaBisectorConstraint(const Vector3r &p0, Real invMass0, const Vector3r &p1, Real invMass1, const Vector3r &p2, Real invMass2, const Real stiffness, Vector3r &corr0, Vector3r &corr1, Vector3r &corr2)

Determine the position corrections for a perpendicular bisector constraint:\(C(\mathbf{p}_0, \mathbf{p}_1, \mathbf{p}_2) = ( \mathbf{p}_1 - 0.5 (\mathbf{p}_0 + \mathbf{p}_1))^T (\mathbf{p}_1 - \mathbf{p}_0) = 0\)

Parameters
  • p0 – position of first particle

  • invMass0 – inverse mass of first particle

  • p1 – position of second particle

  • invMass1 – inverse mass of second particle

  • p2 – position of third particle

  • invMass2 – inverse mass of third particle

  • stiffness – stiffness coefficient

  • corr0 – position correction of first particle

  • corr1 – position correction of second particle

  • corr2 – position correction of third particle

static bool solve_GhostPointEdgeDistanceConstraint(const Vector3r &p0, Real invMass0, const Vector3r &p1, Real invMass1, const Vector3r &p2, Real invMass2, const Real stiffness, const Real ghostEdgeRestLength, Vector3r &corr0, Vector3r &corr1, Vector3r &corr2)

Determine the position corrections for a constraint that enforces a rest length between an edge and a ghost point:\(C(\mathbf{p}_0, \mathbf{p}_1, \mathbf{p}_2) = \| ( 0.5 (\mathbf{p}_0 + \mathbf{p}_1) - \mathbf{p}_2 \| - L_0 = 0\)

Parameters
  • p0 – position of first particle

  • invMass0 – inverse mass of first particle

  • p1 – position of second particle

  • invMass1 – inverse mass of second particle

  • p2 – position of third particle

  • invMass2 – inverse mass of third particle

  • stiffness – stiffness coefficient

  • ghostEdgeRestLength – rest length

  • corr0 – position correction of first particle

  • corr1 – position correction of second particle

  • corr2 – position correction of third particle

static bool solve_DarbouxVectorConstraint(const Vector3r &p0, Real invMass0, const Vector3r &p1, Real invMass1, const Vector3r &p2, Real invMass2, const Vector3r &p3, Real invMass3, const Vector3r &p4, Real invMass4, const Vector3r &bendingAndTwistingKs, const Real midEdgeLength, const Vector3r &restDarbouxVector, Vector3r &oa, Vector3r &ob, Vector3r &oc, Vector3r &od, Vector3r &oe)

Determine the position corrections for the Darboux vector constraint (eq. 21 in the paper). See the paper appendix for derivation details

Parameters
  • p0 – position of first particle

  • invMass0 – inverse mass of first particle

  • p1 – position of second particle

  • invMass1 – inverse mass of second particle

  • p2 – position of third particle

  • invMass2 – inverse mass of third particle

  • p3 – position of fourth particle

  • invMass3 – inverse mass of fourth particle

  • p4 – position of fifth particle

  • invMass4 – inverse mass of fifth particle

  • bendingAndTwistingKs – stiffness coefficients for bending and twisting

  • midEdgeLength – average edge length

  • restDarbouxVector – Darboux vector in rest state

  • corr0 – position correction of first particle

  • corr1 – position correction of second particle

  • corr2 – position correction of third particle

  • corr3 – position correction of fourth particle

  • corr4 – position correction of fifth particle

static bool computeMaterialFrame(const Vector3r &p0, const Vector3r &p1, const Vector3r &p2, Matrix3r &frame)

Computes the material frame (eq. 3 in the paper)

static bool computeDarbouxVector(const Matrix3r &dA, const Matrix3r &dB, const Real mid_edge_length, Vector3r &darboux_vector)

Computes the Darboux Vector (eq. 10 in the paper)

static bool computeMaterialFrameDerivative(const Vector3r &p0, const Vector3r &p1, const Vector3r &p2, const Matrix3r &d, Matrix3r &d1p0, Matrix3r &d1p1, Matrix3r &d1p2, Matrix3r &d2p0, Matrix3r &d2p1, Matrix3r &d2p2, Matrix3r &d3p0, Matrix3r &d3p1, Matrix3r &d3p2)

Computes the material frame derivatives (eq. 43, 44 and 45 in the appendix)

static bool computeDarbouxGradient(const Vector3r &darboux_vector, const Real length, const Matrix3r &dA, const Matrix3r &dB, const Matrix3r[3][3], const Matrix3r[3][3], Matrix3r &omega_pa, Matrix3r &omega_pb, Matrix3r &omega_pc, Matrix3r &omega_pd, Matrix3r &omega_pe)

Compute the Darboux gradient in respect to each point (eq. 49-53 in the appendix)