Class PositionBasedDynamics
Defined in File PositionBasedDynamics.h
Class Documentation
-
class PositionBasedDynamics
Public Static Functions
-
static bool solve_DistanceConstraint(const Vector3r &p0, Real invMass0, const Vector3r &p1, Real invMass1, const Real restLength, const Real stiffness, Vector3r &corr0, Vector3r &corr1)
Determine the position corrections for a distance constraint between two particles:\(C(\mathbf{p}_0, \mathbf{p}_1) = \| \mathbf{p}_0 - \mathbf{p}_1\| - l_0 = 0\) More information can be found in the following papers: [12], [2], [5], [6],
- Parameters
p0 – position of first particle
invMass0 – inverse mass of first particle
p1 – position of second particle
invMass1 – inverse mass of second particle
restLength – rest length of distance constraint
stiffness – stiffness coefficient
corr0 – position correction of first particle
corr1 – position correction of second particle
-
static bool solve_DihedralConstraint(const Vector3r &p0, Real invMass0, const Vector3r &p1, Real invMass1, const Vector3r &p2, Real invMass2, const Vector3r &p3, Real invMass3, const Real restAngle, const Real stiffness, Vector3r &corr0, Vector3r &corr1, Vector3r &corr2, Vector3r &corr3)
Determine the position corrections for a dihedral bending constraint. For a pair of adjacent triangles \((\mathbf{p}_1, \mathbf{p}_3, \mathbf{p}_2)\) and \((\mathbf{p}_1, \mathbf{p}_2, \mathbf{p}_4)\) with the common edge \((\mathbf{p}_3, \mathbf{p}_4)\) a bilateral bending constraint is added by the constraint function
\[\begin{equation*} C_{bend}(\mathbf{p}_1, \mathbf{p}_2,\mathbf{p}_3, \mathbf{p}_4) = \text{acos}\left( \frac{\mathbf{p}_{2,1} \times \mathbf{p}_{3,1}}{|\mathbf{p}_{2,1} \times \mathbf{p}_{3,1}|} \cdot \frac{\mathbf{p}_{2,1} \times \mathbf{p}_{4,1}}{|\mathbf{p}_{2,1} \times \mathbf{p}_{4,1}|}\right)-\varphi_0 \end{equation*}\]and stiffness \(k_{bend}\). The scalar \(\varphi_0\) is the initial dihedral angle between the two triangles and \(k_{bend}\)is a global user parameter defining the bending stiffness.
More information can be found in the following papers:
[12], [2], [5], [6],- 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
restAngle – rest angle \(\varphi_0\)
stiffness – stiffness coefficient
corr0 – position correction of first particle
corr1 – position correction of second particle
corr2 – position correction of third particle
corr3 – position correction of fourth particle
-
static bool solve_VolumeConstraint(const Vector3r &p0, Real invMass0, const Vector3r &p1, Real invMass1, const Vector3r &p2, Real invMass2, const Vector3r &p3, Real invMass3, const Real restVolume, const Real stiffness, Vector3r &corr0, Vector3r &corr1, Vector3r &corr2, Vector3r &corr3)
Determine the position corrections for a constraint that conserves the volume of single tetrahedron. Such a constraint has the form
\[\begin{equation*} C(\mathbf{p}_1, \mathbf{p}_2, \mathbf{p}_3, \mathbf{p}_4) = \frac{1}{6} \left(\mathbf{p}_{2,1} \times \mathbf{p}_{3,1}\right) \cdot \mathbf{p}_{4,1} - V_0, \end{equation*}\]where \(\mathbf{p}_1\), \(\mathbf{p}_2\), \(\mathbf{p}_3\) and \(\mathbf{p}_4\) are the four corners of the tetrahedron and \(V_0\)is its rest volume.
More information can be found in the following papers:
[12], [2], [5], [6],- 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
restVolume – rest angle \(V_0\)
stiffness – stiffness coefficient
corr0 – position correction of first particle
corr1 – position correction of second particle
corr2 – position correction of third particle
corr3 – position correction of fourth particle
-
static bool solve_EdgePointDistanceConstraint(const Vector3r &p, Real invMass, const Vector3r &p0, Real invMass0, const Vector3r &p1, Real invMass1, const Real restDist, const Real compressionStiffness, const Real stretchStiffness, Vector3r &corr, Vector3r &corr0, Vector3r &corr1)
Determine the position corrections for a constraint that preserves a rest distance between a point and an edge.
- Parameters
p – position of point particle
invMass – inverse mass of point particle
p0 – position of first edge particle
invMass0 – inverse mass of first edge particle
p1 – position of second edge particle
invMass1 – inverse mass of second edge particle
restDist – rest distance of point and edge
compressionStiffness – stiffness coefficient for compression
stretchStiffness – stiffness coefficient for stretching
corr – position correction of point particle
corr0 – position correction of first edge particle
corr1 – position correction of second edge particle
-
static bool solve_TrianglePointDistanceConstraint(const Vector3r &p, Real invMass, const Vector3r &p0, Real invMass0, const Vector3r &p1, Real invMass1, const Vector3r &p2, Real invMass2, const Real restDist, const Real compressionStiffness, const Real stretchStiffness, Vector3r &corr, Vector3r &corr0, Vector3r &corr1, Vector3r &corr2)
Determine the position corrections for a constraint that preserves a rest distance between a point and a triangle.
- Parameters
p – position of point particle
invMass – inverse mass of point particle
p0 – position of first triangle particle
invMass0 – inverse mass of first triangle particle
p1 – position of second triangle particle
invMass1 – inverse mass of second triangle particle
p2 – position of third triangle particle
invMass2 – inverse mass of third triangle particle
restDist – rest distance of point and triangle
compressionStiffness – stiffness coefficient for compression
stretchStiffness – stiffness coefficient for stretching
corr – position correction of point particle
corr0 – position correction of first triangle particle
corr1 – position correction of second triangle particle
corr2 – position correction of third triangle particle
-
static bool solve_EdgeEdgeDistanceConstraint(const Vector3r &p0, Real invMass0, const Vector3r &p1, Real invMass1, const Vector3r &p2, Real invMass2, const Vector3r &p3, Real invMass3, const Real restDist, const Real compressionStiffness, const Real stretchStiffness, Vector3r &corr0, Vector3r &corr1, Vector3r &corr2, Vector3r &corr3)
Determine the position corrections for a constraint that preserves a rest distance between two edges.
- Parameters
p0 – position of first particle of edge 0
invMass0 – inverse mass of first particle of edge 0
p1 – position of second particle of edge 0
invMass1 – inverse mass of second particle of edge 0
p2 – position of first particle of edge 1
invMass2 – inverse mass of first particle of edge 1
p3 – position of second particle of edge 1
invMass3 – inverse mass of second particle of edge 1
restDist – rest distance between both edges
compressionStiffness – stiffness coefficient for compression
stretchStiffness – stiffness coefficient for stretching
corr0 – position correction of first particle of edge 0
corr1 – position correction of second particle of edge 0
corr2 – position correction of first particle of edge 1
corr3 – position correction of second particle of edge 1
-
static bool init_IsometricBendingConstraint(const Vector3r &p0, const Vector3r &p1, const Vector3r &p2, const Vector3r &p3, Matrix4r &Q)
Initialize the local stiffness matrix Q. The matrix is required by the solver step. It must only be recomputed if the rest shape changes.
Bending is simulated for the angle on (p2, p3) between the triangles (p0, p2, p3) and (p1, p3, p2).
- Parameters
p0 – point 0 of stencil
p1 – point 1 of stencil
p2 – point 2 of stencil
p3 – point 3 of stencil
Q – returns the local stiffness matrix which is required by the solver
-
static bool solve_IsometricBendingConstraint(const Vector3r &p0, Real invMass0, const Vector3r &p1, Real invMass1, const Vector3r &p2, Real invMass2, const Vector3r &p3, Real invMass3, const Matrix4r &Q, const Real stiffness, Vector3r &corr0, Vector3r &corr1, Vector3r &corr2, Vector3r &corr3)
Determine the position corrections for the isometric bending constraint. This constraint can be used for almost inextensible surface models.
More information can be found in:
[6], [4]- 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
Q – local stiffness matrix which must be initialized by calling init_IsometricBendingConstraint()
stiffness – stiffness coefficient for bending
corr0 – position correction of first particle
corr1 – position correction of second particle
corr2 – position correction of third particle
corr3 – position correction of fourth particle
-
static bool init_ShapeMatchingConstraint(const Vector3r x0[], const Real invMasses[], const int numPoints, Vector3r &restCm)
Initialize rest configuration infos for one shape matching cluster which are required by the solver step. It must only be reinitialized if the rest shape changes.
- Parameters
x0 – rest configuration of all particles in the cluster
invMasses – inverse masses of all particles in the cluster
numPoints – number of particles in the cluster
restCm – returns the center of mass of the rest configuration
-
static bool solve_ShapeMatchingConstraint(const Vector3r x0[], const Vector3r x[], const Real invMasses[], const int numPoints, const Vector3r &restCm, const Real stiffness, const bool allowStretch, Vector3r corr[], Matrix3r *rot = NULL)
Determine the position corrections for a shape matching constraint.
More information can be found in:
[6], [2], [5], [11], [3], [8]- Parameters
x0 – rest configuration of all particles in the cluster
x – current configuration of all particles in the cluster
invMasses – invMasses inverse masses of all particles in the cluster
numPoints – number of particles in the cluster
restCm – center of mass of the rest configuration
stiffness – stiffness coefficient
allowStretch – allow stretching
corr – position corrections for all particles in the cluster
rot – returns determined rotation matrix
-
static bool init_StrainTriangleConstraint(const Vector3r &p0, const Vector3r &p1, const Vector3r &p2, Matrix2r &invRestMat)
Initialize rest configuration infos which are required by the solver step. Recomputation is only necessary when rest shape changes.
The triangle is defined in the xy plane.
- Parameters
p0 – point 0 of triangle
p1 – point 1 of triangle
p2 – point 2 of triangle
invRestMat – returns a matrix required by the solver
-
static bool solve_StrainTriangleConstraint(const Vector3r &p0, Real invMass0, const Vector3r &p1, Real invMass1, const Vector3r &p2, Real invMass2, const Matrix2r &invRestMat, const Real xxStiffness, const Real yyStiffness, const Real xyStiffness, const bool normalizeStretch, const bool normalizeShear, Vector3r &corr0, Vector3r &corr1, Vector3r &corr2)
Solve triangle constraint with strain based dynamics and return position corrections.
More information can be found in:
[6], [13]- 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
invRestMat – precomputed matrix determined by init_StrainTriangleConstraint()
xxStiffness – stiffness coefficient for xx stretching
yyStiffness – stiffness coefficient for yy stretching
xyStiffness – stiffness coefficient for xy shearing
normalizeStretch – should stretching be normalized
normalizeShear – should shearing be normalized
corr0 – position correction for point 0
corr1 – position correction for point 1
corr2 – position correction for point 2
-
static bool init_StrainTetraConstraint(const Vector3r &p0, const Vector3r &p1, const Vector3r &p2, const Vector3r &p3, Matrix3r &invRestMat)
Initialize rest configuration infos which are required by the solver step. Recomputation is only necessary when rest shape changes.
- Parameters
p0 – point 0 of tet
p1 – point 1 of tet
p2 – point 2 of tet
p3 – point 3 of tet
invRestMat – returns a matrix required by the solver
-
static bool solve_StrainTetraConstraint(const Vector3r &p0, Real invMass0, const Vector3r &p1, Real invMass1, const Vector3r &p2, Real invMass2, const Vector3r &p3, Real invMass3, const Matrix3r &invRestMat, const Vector3r &stretchStiffness, const Vector3r &shearStiffness, const bool normalizeStretch, const bool normalizeShear, Vector3r &corr0, Vector3r &corr1, Vector3r &corr2, Vector3r &corr3)
-
static void computeGradCGreen(Real restVolume, const Matrix3r &invRestMat, const Matrix3r &sigma, Vector3r *J)
-
static void computeGreenStrainAndPiolaStress(const Vector3r &x1, const Vector3r &x2, const Vector3r &x3, const Vector3r &x4, const Matrix3r &invRestMat, const Real restVolume, const Real mu, const Real lambda, Matrix3r &epsilon, Matrix3r &sigma, Real &energy)
-
static void computeGreenStrainAndPiolaStressInversion(const Vector3r &x1, const Vector3r &x2, const Vector3r &x3, const Vector3r &x4, const Matrix3r &invRestMat, const Real restVolume, const Real mu, const Real lambda, Matrix3r &epsilon, Matrix3r &sigma, Real &energy)
-
static bool init_FEMTriangleConstraint(const Vector3r &p0, const Vector3r &p1, const Vector3r &p2, Real &area, Matrix2r &invRestMat)
Implementation of the finite element method described in
Jan Bender, Dan Koschier, Patrick Charrier and Daniel Weber,
”Position-Based Simulation of Continuous Materials”,
Computers & Graphics 44, 2014
Initialize rest configuration infos which are required by the solver step. Recomputation is only necessary when rest shape changes.
-
static bool solve_FEMTriangleConstraint(const Vector3r &p0, Real invMass0, const Vector3r &p1, Real invMass1, const Vector3r &p2, Real invMass2, const Real &area, const Matrix2r &invRestMat, const Real youngsModulusX, const Real youngsModulusY, const Real youngsModulusShear, const Real poissonRatioXY, const Real poissonRatioYX, Vector3r &corr0, Vector3r &corr1, Vector3r &corr2)
Implementation of the finite element method described in
Jan Bender, Dan Koschier, Patrick Charrier and Daniel Weber,
”Position-Based Simulation of Continuous Materials”,
Computers & Graphics 44, 2014
Solve the continuum mechanical constraint defined for a triangle.
-
static bool init_FEMTetraConstraint(const Vector3r &p0, const Vector3r &p1, const Vector3r &p2, const Vector3r &p3, Real &volume, Matrix3r &invRestMat)
Implementation of the finite element method described in
Jan Bender, Dan Koschier, Patrick Charrier and Daniel Weber,
”Position-Based Simulation of Continuous Materials”,
Computers & Graphics 44, 2014
Initialize rest configuration infos which are required by the solver step. Recomputation is only necessary when rest shape changes.
-
static bool solve_FEMTetraConstraint(const Vector3r &p0, Real invMass0, const Vector3r &p1, Real invMass1, const Vector3r &p2, Real invMass2, const Vector3r &p3, Real invMass3, const Real restVolume, const Matrix3r &invRestMat, const Real youngsModulus, const Real poissonRatio, const bool handleInversion, Vector3r &corr0, Vector3r &corr1, Vector3r &corr2, Vector3r &corr3)
Implementation of the finite element method described in
Jan Bender, Dan Koschier, Patrick Charrier and Daniel Weber,
”Position-Based Simulation of Continuous Materials”,
Computers & Graphics 44, 2014
Solve the continuum mechanical constraint defined for a tetrahedron.
-
static bool init_ParticleTetContactConstraint(const Real invMass0, const Vector3r &x0, const Vector3r &v0, const Real invMass[], const Vector3r x[], const Vector3r v[], const Vector3r &bary, const Vector3r &normal, Eigen::Matrix<Real, 3, 3, Eigen::DontAlign> &constraintInfo)
Initialize contact between a particle and a tetrahedron and return info which is required by the solver step.
- Parameters
invMass0 – inverse mass of particle which collides with tet
x0 – particle position
v0 – particle velocity
invMass – inverse masses of tet particles
x – positions of tet particles
v – velocities of tet particles
bary – barycentric coordinates of contact point in tet
normal – contact normal in body 1
constraintInfo –
Stores the local and global position of the contact points and the contact normal.
The joint info contains the following columns:
0: contact normal in body 1 (global)
1: contact tangent (global)
0,2: 1.0 / normal^T * K * normal
1,2: maximal impulse in tangent direction
-
static bool solve_ParticleTetContactConstraint(const Real invMass0, const Vector3r &x0, const Real invMass[], const Vector3r x[], const Vector3r &bary, Eigen::Matrix<Real, 3, 3, Eigen::DontAlign> &constraintInfo, Real &lambda, Vector3r &corr0, Vector3r corr[])
Perform a solver step for a contact constraint between a particle and a tetrahedron. A contact constraint handles collisions and resting contacts between the bodies. The contact info must be generated in each time step.
- Parameters
invMass0 – inverse mass of particle which collides with tet
x0 – particle position
invMass – inverse masses of tet particles
x – positions of tet particles
bary – barycentric coordinates of contact point in tet
constraintInfo – information which is required by the solver. This information must be generated in the beginning by calling init_RigidBodyContactConstraint().
corr0 – position correction of particle
corr – position corrections of tet particles
-
static bool velocitySolve_ParticleTetContactConstraint(const Real invMass0, const Vector3r &x0, const Vector3r &v0, const Real invMass[], const Vector3r x[], const Vector3r v[], const Vector3r &bary, const Real lambda, const Real frictionCoeff, Eigen::Matrix<Real, 3, 3, Eigen::DontAlign> &constraintInfo, Vector3r &corr_v0, Vector3r corr_v[])
Perform a solver step for a contact constraint between a particle and a tetrahedron. A contact constraint handles collisions and resting contacts between the bodies. The contact info must be generated in each time step.
- Parameters
invMass0 – inverse mass of particle which collides with tet
x0 – particle position
v0 – particle velocity
invMass – inverse masses of tet particles
x – positions of tet particles
v – velocities of tet particles
bary – barycentric coordinates of contact point in tet
frictionCoeff – friction coefficient
constraintInfo – information which is required by the solver. This information must be generated in the beginning by calling init_RigidBodyContactConstraint().
corr_v0 – velocity correction of particle
corr_v – velocity corrections of tet particles
-
static bool solve_DistanceConstraint(const Vector3r &p0, Real invMass0, const Vector3r &p1, Real invMass1, const Real restLength, const Real stiffness, Vector3r &corr0, Vector3r &corr1)