Friction¶
Friction Constraints¶
- class ipctk.FrictionConstraints¶
Bases:
pybind11_object
Public Data Attributes:
Public Methods:
__init__
(self)build
(*args, **kwargs)Overloaded function.
compute_potential
(self, mesh, velocity, n]], ...)Compute the friction dissapative potential from the given velocity.
compute_potential_gradient
(self, mesh, ...)Compute the gradient of the friction dissapative potential wrt the velocity.
compute_potential_hessian
(self, mesh, ...)Compute the Hessian of the friction dissapative potential wrt the velocity.
compute_force
(self, mesh, rest_positions, ...)Compute the friction force from the given velocity.
compute_force_jacobian
(self, mesh, ...)Compute the Jacobian of the friction force wrt the velocity.
__len__
(self)Get the number of friction constraints.
empty
(self)Get if the friction constraints are empty.
clear
(self)Clear the friction constraints.
__getitem__
(self, idx)Get a reference to constriant idx.
default_blend_mu
(mu0, mu1)Inherited from
pybind11_object
-
__annotations__ =
{}
¶
-
__module__ =
'ipctk'
¶
- build(*args, **kwargs)¶
Overloaded function.
build(self: ipctk.FrictionConstraints, mesh: ipc::CollisionMesh, vertices: numpy.ndarray[numpy.float64[m, n]], contact_constraints: ipctk.CollisionConstraints, dhat: float, barrier_stiffness: float, mu: float) -> None
build(self: ipctk.FrictionConstraints, mesh: ipc::CollisionMesh, vertices: numpy.ndarray[numpy.float64[m, n]], contact_constraints: ipctk.CollisionConstraints, dhat: float, barrier_stiffness: float, mus: numpy.ndarray[numpy.float64[m, 1]]) -> None
build(self: ipctk.FrictionConstraints, mesh: ipc::CollisionMesh, vertices: numpy.ndarray[numpy.float64[m, n]], contact_constraints: ipctk.CollisionConstraints, dhat: float, barrier_stiffness: float, mus: numpy.ndarray[numpy.float64[m, 1]], blend_mu: Callable[[float, float], float]) -> None
- compute_force(self: ipctk.FrictionConstraints, mesh: ipc::CollisionMesh, rest_positions: numpy.ndarray[numpy.float64[m, n]], lagged_displacements: numpy.ndarray[numpy.float64[m, n]], velocities: numpy.ndarray[numpy.float64[m, n]], dhat: float, barrier_stiffness: float, epsv: float, dmin: float = 0, no_mu: bool = False) numpy.ndarray[numpy.float64[m, 1]] ¶
Compute the friction force from the given velocity.
- Parameters:¶
- mesh
The collision mesh.
- rest_positions
Rest positions of the vertices (rowwise).
- lagged_displacements
Previous displacements of the vertices (rowwise).
- velocities
Current displacements of the vertices (rowwise).
- dhat
Barrier activation distance.
- barrier_stiffness
Barrier stiffness.
- epsv
Mollifier parameter \(\epsilon_v\).
- dmin
Minimum distance to use for the barrier.
- no_mu
whether to not multiply by mu
- Returns:¶
The friction force.
- compute_force_jacobian(self: ipctk.FrictionConstraints, mesh: ipc::CollisionMesh, rest_positions: numpy.ndarray[numpy.float64[m, n]], lagged_displacements: numpy.ndarray[numpy.float64[m, n]], velocities: numpy.ndarray[numpy.float64[m, n]], dhat: float, barrier_stiffness: float, epsv: float, wrt: ipc::FrictionConstraint::DiffWRT, dmin: float = 0) scipy.sparse.csc_matrix[numpy.float64] ¶
Compute the Jacobian of the friction force wrt the velocity.
- Parameters:¶
- mesh
The collision mesh.
- rest_positions
Rest positions of the vertices (rowwise).
- lagged_displacements
Previous displacements of the vertices (rowwise).
- velocities
Current displacements of the vertices (rowwise).
- dhat
Barrier activation distance.
- barrier_stiffness
Barrier stiffness.
- epsv
Mollifier parameter \(\epsilon_v\).
- wrt
The variable to take the derivative with respect to.
- dmin
Minimum distance to use for the barrier.
- Returns:¶
The Jacobian of the friction force wrt the velocity.
- compute_potential(self: ipctk.FrictionConstraints, mesh: ipc::CollisionMesh, velocity: numpy.ndarray[numpy.float64[m, n]], epsv: float) float ¶
Compute the friction dissapative potential from the given velocity.
- compute_potential_gradient(self: ipctk.FrictionConstraints, mesh: ipc::CollisionMesh, velocity: numpy.ndarray[numpy.float64[m, n]], epsv: float) numpy.ndarray[numpy.float64[m, 1]] ¶
Compute the gradient of the friction dissapative potential wrt the velocity.
- compute_potential_hessian(self: ipctk.FrictionConstraints, mesh: ipc::CollisionMesh, velocity: numpy.ndarray[numpy.float64[m, n]], epsv: float, project_hessian_to_psd: bool = False) scipy.sparse.csc_matrix[numpy.float64] ¶
Compute the Hessian of the friction dissapative potential wrt the velocity.
- property ee_constraints : list[ipc::EdgeEdgeFrictionConstraint]¶
- property ev_constraints : list[ipc::EdgeVertexFrictionConstraint]¶
- property fv_constraints : list[ipc::FaceVertexFrictionConstraint]¶
- property vv_constraints : list[ipc::VertexVertexFrictionConstraint]¶
-
__annotations__ =
Friction Constraint¶
- class ipctk.FrictionConstraint¶
Bases:
CollisionStencil
Public Data Attributes:
Contact force magnitude
Coefficient of friction
Weight
Gradient of weight with respect to all DOF
Barycentric coordinates of the closest point(s)
Tangent basis of the contact (max size 3×2)
Public Methods:
__init__
(*args, **kwargs)compute_potential
(self, velocities, edges, ...)Compute the friction dissapative potential.
compute_potential_gradient
(self, velocities, ...)Compute the friction dissapative potential gradient wrt velocities.
compute_potential_hessian
(self, velocities, ...)Compute the friction dissapative potential hessian wrt velocities.
compute_force
(self, rest_positions, ...[, ...])Compute the friction force.
compute_force_jacobian
(self, rest_positions, ...)Compute the friction force Jacobian.
Inherited from
CollisionStencil
__init__
(*args, **kwargs)num_vertices
(self)Get the number of vertices in the contact stencil.
vertex_ids
(self, edges, faces)Get the vertex IDs of the contact stencil.
vertices
(self, vertices, edges, faces)Get the vertex attributes of the contact stencil.
dof
(self, X, edges, faces)Select this stencil's DOF from the full matrix of DOF.
compute_distance
(self, vertices, edges, faces)Compute the distance of the stencil.
compute_distance_gradient
(self, vertices, ...)Compute the distance gradient of the stencil w.r.t.
compute_distance_hessian
(self, vertices, ...)Compute the distance Hessian of the stencil w.r.t.
Inherited from
pybind11_object
- class DiffWRT¶
Bases:
pybind11_object
Members:
REST_POSITIONS
LAGGED_DISPLACEMENTS
VELOCITIES
-
LAGGED_DISPLACEMENTS =
<DiffWRT.LAGGED_DISPLACEMENTS: 1>
¶
-
REST_POSITIONS =
<DiffWRT.REST_POSITIONS: 0>
¶
-
VELOCITIES =
<DiffWRT.VELOCITIES: 2>
¶
-
__annotations__ =
{}
¶
-
__members__ =
{'LAGGED_DISPLACEMENTS': <DiffWRT.LAGGED_DISPLACEMENTS: 1>, 'REST_POSITIONS': <DiffWRT.REST_POSITIONS: 0>, 'VELOCITIES': <DiffWRT.VELOCITIES: 2>}
¶
-
__module__ =
'ipctk'
¶
- __str__()¶
name(self: handle) -> str
- property name : str¶
- property value : int¶
-
LAGGED_DISPLACEMENTS =
-
LAGGED_DISPLACEMENTS =
<DiffWRT.LAGGED_DISPLACEMENTS: 1>
¶
-
REST_POSITIONS =
<DiffWRT.REST_POSITIONS: 0>
¶
-
VELOCITIES =
<DiffWRT.VELOCITIES: 2>
¶
-
__annotations__ =
{}
¶
-
__module__ =
'ipctk'
¶
- property closest_point : numpy.ndarray[numpy.float64[m, 1]]¶
Barycentric coordinates of the closest point(s)
-
compute_force(self, rest_positions: numpy.ndarray[numpy.float64[m, n]], lagged_displacements: numpy.ndarray[numpy.float64[m, n]], velocities: numpy.ndarray[numpy.float64[m, n]], edges: numpy.ndarray[numpy.int32[m, n]], faces: numpy.ndarray[numpy.int32[m, n]], dhat: float, barrier_stiffness: float, epsv: float, dmin: float =
0
, no_mu: bool =False
) numpy.ndarray[numpy.float64[m, 1]] ¶ Compute the friction force.
- Parameters:¶
- rest_positions: numpy.ndarray[numpy.float64[m, n]]¶
Rest positions of the vertices (rowwise).
- lagged_displacements: numpy.ndarray[numpy.float64[m, n]]¶
Previous displacements of the vertices (rowwise).
- velocities: numpy.ndarray[numpy.float64[m, n]]¶
Current displacements of the vertices (rowwise).
- edges: numpy.ndarray[numpy.int32[m, n]]¶
Collision mesh edges
- faces: numpy.ndarray[numpy.int32[m, n]]¶
Collision mesh faces
- dhat: float¶
Barrier activation distance
- barrier_stiffness: float¶
Barrier stiffness
- epsv: float¶
Smooth friction mollifier parameter \(\epsilon_v\).
- dmin: float =
0
¶ Minimum distance
- no_mu: bool =
False
¶ Whether to not multiply by mu
- Returns:¶
Friction force
-
compute_force_jacobian(self, rest_positions: numpy.ndarray[numpy.float64[m, n]], lagged_displacements: numpy.ndarray[numpy.float64[m, n]], velocities: numpy.ndarray[numpy.float64[m, n]], edges: numpy.ndarray[numpy.int32[m, n]], faces: numpy.ndarray[numpy.int32[m, n]], dhat: float, barrier_stiffness: float, epsv: float, wrt: ipctk.FrictionConstraint.DiffWRT, dmin: float =
0
) numpy.ndarray[numpy.float64[m, n]] ¶ Compute the friction force Jacobian.
- Parameters:¶
- rest_positions: numpy.ndarray[numpy.float64[m, n]]¶
Rest positions of the vertices (rowwise).
- lagged_displacements: numpy.ndarray[numpy.float64[m, n]]¶
Previous displacements of the vertices (rowwise).
- velocities: numpy.ndarray[numpy.float64[m, n]]¶
Current displacements of the vertices (rowwise).
- edges: numpy.ndarray[numpy.int32[m, n]]¶
Collision mesh edges
- faces: numpy.ndarray[numpy.int32[m, n]]¶
Collision mesh faces
- dhat: float¶
Barrier activation distance
- barrier_stiffness: float¶
Barrier stiffness
- epsv: float¶
Smooth friction mollifier parameter \(\epsilon_v\).
- wrt: ipctk.FrictionConstraint.DiffWRT¶
Variable to differentiate the friction force with respect to.
- dmin: float =
0
¶ Minimum distance
- Returns:¶
Friction force Jacobian
- compute_potential(self, velocities: numpy.ndarray[numpy.float64[m, n]], edges: numpy.ndarray[numpy.int32[m, n]], faces: numpy.ndarray[numpy.int32[m, n]], epsv: float) float ¶
Compute the friction dissapative potential.
- compute_potential_gradient(self, velocities: numpy.ndarray[numpy.float64[m, n]], edges: numpy.ndarray[numpy.int32[m, n]], faces: numpy.ndarray[numpy.int32[m, n]], epsv: float) numpy.ndarray[numpy.float64[m, 1]] ¶
Compute the friction dissapative potential gradient wrt velocities.
- compute_potential_hessian(self, velocities: numpy.ndarray[numpy.float64[m, n]], edges: numpy.ndarray[numpy.int32[m, n]], faces: numpy.ndarray[numpy.int32[m, n]], epsv: float, project_hessian_to_psd: bool) numpy.ndarray[numpy.float64[m, n]] ¶
Compute the friction dissapative potential hessian wrt velocities.
- Parameters:¶
- velocities: numpy.ndarray[numpy.float64[m, n]]¶
Velocities of the vertices (rowwise)
- edges: numpy.ndarray[numpy.int32[m, n]]¶
Edges of the mesh
- faces: numpy.ndarray[numpy.int32[m, n]]¶
Faces of the mesh
- epsv: float¶
Smooth friction mollifier parameter \(\epsilon_v\).
- project_hessian_to_psd: bool¶
Project the hessian to PSD
- Returns:¶
Hessian of the friction dissapative potential wrt velocities
- property mu : float¶
Coefficient of friction
- property normal_force_magnitude : float¶
Contact force magnitude
- property tangent_basis : numpy.ndarray[numpy.float64[m, n]]¶
Tangent basis of the contact (max size 3×2)
- property weight : float¶
Weight
- property weight_gradient : scipy.sparse.csc_matrix[numpy.float64]¶
Gradient of weight with respect to all DOF
Vertex-Vertex Friction Constraint¶
- class ipctk.VertexVertexFrictionConstraint¶
Bases:
VertexVertexCandidate
,FrictionConstraint
Public Data Attributes:
Inherited from
VertexVertexCandidate
ID of the first vertex
ID of the second vertex
Inherited from
FrictionConstraint
Contact force magnitude
Coefficient of friction
Weight
Gradient of weight with respect to all DOF
Barycentric coordinates of the closest point(s)
Tangent basis of the contact (max size 3×2)
Public Methods:
__init__
(*args, **kwargs)Overloaded function.
Inherited from
VertexVertexCandidate
__init__
(self, vertex0_id, vertex1_id)__str__
(self)__repr__
(self)num_vertices
(self)vertex_ids
(self, edges, faces)Get the indices of the vertices
__eq__
(self, other)__ne__
(self, other)__lt__
(self, other)Compare EdgeVertexCandidates for sorting.
Inherited from
FrictionConstraint
__init__
(*args, **kwargs)compute_potential
(self, velocities, edges, ...)Compute the friction dissapative potential.
compute_potential_gradient
(self, velocities, ...)Compute the friction dissapative potential gradient wrt velocities.
compute_potential_hessian
(self, velocities, ...)Compute the friction dissapative potential hessian wrt velocities.
compute_force
(self, rest_positions, ...[, ...])Compute the friction force.
compute_force_jacobian
(self, rest_positions, ...)Compute the friction force Jacobian.
Inherited from
CollisionStencil
__init__
(*args, **kwargs)num_vertices
(self)Get the number of vertices in the contact stencil.
vertex_ids
(self, edges, faces)Get the vertex IDs of the contact stencil.
vertices
(self, vertices, edges, faces)Get the vertex attributes of the contact stencil.
dof
(self, X, edges, faces)Select this stencil's DOF from the full matrix of DOF.
compute_distance
(self, vertices, edges, faces)Compute the distance of the stencil.
compute_distance_gradient
(self, vertices, ...)Compute the distance gradient of the stencil w.r.t.
compute_distance_hessian
(self, vertices, ...)Compute the distance Hessian of the stencil w.r.t.
Inherited from
pybind11_object
-
__annotations__ =
{}
¶
- __init__(*args, **kwargs)¶
Overloaded function.
__init__(self: ipctk.VertexVertexFrictionConstraint, constraint: ipctk.VertexVertexConstraint) -> None
__init__(self: ipctk.VertexVertexFrictionConstraint, constraint: ipctk.VertexVertexConstraint, vertices: numpy.ndarray[numpy.float64[m, n]], edges: numpy.ndarray[numpy.int32[m, n]], faces: numpy.ndarray[numpy.int32[m, n]], dhat: float, barrier_stiffness: float) -> None
-
__module__ =
'ipctk'
¶
-
__annotations__ =
Edge-Vertex Friction Constraint¶
- class ipctk.EdgeVertexFrictionConstraint¶
Bases:
EdgeVertexCandidate
,FrictionConstraint
Public Data Attributes:
Inherited from
EdgeVertexCandidate
ID of the edge
ID of the vertex
Inherited from
FrictionConstraint
Contact force magnitude
Coefficient of friction
Weight
Gradient of weight with respect to all DOF
Barycentric coordinates of the closest point(s)
Tangent basis of the contact (max size 3×2)
Public Methods:
__init__
(*args, **kwargs)Overloaded function.
Inherited from
EdgeVertexCandidate
__init__
(self, edge_id, vertex_id)__str__
(self)__repr__
(self)num_vertices
(self)vertex_ids
(self, edges, faces)ccd
(self, vertices_t0, vertices_t1, edges, faces)Perform narrow-phase CCD on the candidate.
print_ccd_query
(self, vertices_t0, ...)__eq__
(self, other)__ne__
(self, other)__lt__
(self, other)Compare EdgeVertexCandidates for sorting.
Inherited from
FrictionConstraint
__init__
(*args, **kwargs)compute_potential
(self, velocities, edges, ...)Compute the friction dissapative potential.
compute_potential_gradient
(self, velocities, ...)Compute the friction dissapative potential gradient wrt velocities.
compute_potential_hessian
(self, velocities, ...)Compute the friction dissapative potential hessian wrt velocities.
compute_force
(self, rest_positions, ...[, ...])Compute the friction force.
compute_force_jacobian
(self, rest_positions, ...)Compute the friction force Jacobian.
Inherited from
CollisionStencil
__init__
(*args, **kwargs)num_vertices
(self)Get the number of vertices in the contact stencil.
vertex_ids
(self, edges, faces)Get the vertex IDs of the contact stencil.
vertices
(self, vertices, edges, faces)Get the vertex attributes of the contact stencil.
dof
(self, X, edges, faces)Select this stencil's DOF from the full matrix of DOF.
compute_distance
(self, vertices, edges, faces)Compute the distance of the stencil.
compute_distance_gradient
(self, vertices, ...)Compute the distance gradient of the stencil w.r.t.
compute_distance_hessian
(self, vertices, ...)Compute the distance Hessian of the stencil w.r.t.
Inherited from
ContinuousCollisionCandidate
__init__
(*args, **kwargs)ccd
(self, vertices_t0, vertices_t1, edges, faces)Perform narrow-phase CCD on the candidate.
print_ccd_query
(self, vertices_t0, ...)Inherited from
pybind11_object
-
__annotations__ =
{}
¶
- __init__(*args, **kwargs)¶
Overloaded function.
__init__(self: ipctk.EdgeVertexFrictionConstraint, constraint: ipctk.EdgeVertexConstraint) -> None
__init__(self: ipctk.EdgeVertexFrictionConstraint, constraint: ipctk.EdgeVertexConstraint, vertices: numpy.ndarray[numpy.float64[m, n]], edges: numpy.ndarray[numpy.int32[m, n]], faces: numpy.ndarray[numpy.int32[m, n]], dhat: float, barrier_stiffness: float) -> None
-
__module__ =
'ipctk'
¶
-
__annotations__ =
Edge-Edge Friction Constraint¶
- class ipctk.EdgeEdgeFrictionConstraint¶
Bases:
EdgeEdgeCandidate
,FrictionConstraint
Public Data Attributes:
Inherited from
EdgeEdgeCandidate
ID of the first edge
ID of the second edge
Inherited from
FrictionConstraint
Contact force magnitude
Coefficient of friction
Weight
Gradient of weight with respect to all DOF
Barycentric coordinates of the closest point(s)
Tangent basis of the contact (max size 3×2)
Public Methods:
__init__
(*args, **kwargs)Overloaded function.
Inherited from
EdgeEdgeCandidate
__init__
(self, edge0_id, edge1_id)__str__
(self)__repr__
(self)num_vertices
(self)vertex_ids
(self, edges, faces)ccd
(self, vertices_t0, vertices_t1, edges, faces)Perform narrow-phase CCD on the candidate.
print_ccd_query
(self, vertices_t0, ...)__eq__
(self, other)__ne__
(self, other)__lt__
(self, other)Compare EdgeEdgeCandidates for sorting.
Inherited from
FrictionConstraint
__init__
(*args, **kwargs)compute_potential
(self, velocities, edges, ...)Compute the friction dissapative potential.
compute_potential_gradient
(self, velocities, ...)Compute the friction dissapative potential gradient wrt velocities.
compute_potential_hessian
(self, velocities, ...)Compute the friction dissapative potential hessian wrt velocities.
compute_force
(self, rest_positions, ...[, ...])Compute the friction force.
compute_force_jacobian
(self, rest_positions, ...)Compute the friction force Jacobian.
Inherited from
CollisionStencil
__init__
(*args, **kwargs)num_vertices
(self)Get the number of vertices in the contact stencil.
vertex_ids
(self, edges, faces)Get the vertex IDs of the contact stencil.
vertices
(self, vertices, edges, faces)Get the vertex attributes of the contact stencil.
dof
(self, X, edges, faces)Select this stencil's DOF from the full matrix of DOF.
compute_distance
(self, vertices, edges, faces)Compute the distance of the stencil.
compute_distance_gradient
(self, vertices, ...)Compute the distance gradient of the stencil w.r.t.
compute_distance_hessian
(self, vertices, ...)Compute the distance Hessian of the stencil w.r.t.
Inherited from
ContinuousCollisionCandidate
__init__
(*args, **kwargs)ccd
(self, vertices_t0, vertices_t1, edges, faces)Perform narrow-phase CCD on the candidate.
print_ccd_query
(self, vertices_t0, ...)Inherited from
pybind11_object
-
__annotations__ =
{}
¶
- __init__(*args, **kwargs)¶
Overloaded function.
__init__(self: ipctk.EdgeEdgeFrictionConstraint, constraint: ipctk.EdgeEdgeConstraint) -> None
__init__(self: ipctk.EdgeEdgeFrictionConstraint, constraint: ipctk.EdgeEdgeConstraint, vertices: numpy.ndarray[numpy.float64[m, n]], edges: numpy.ndarray[numpy.int32[m, n]], faces: numpy.ndarray[numpy.int32[m, n]], dhat: float, barrier_stiffness: float) -> None
-
__module__ =
'ipctk'
¶
-
__annotations__ =
Triangle-Vertex Friction Constraint¶
- class ipctk.FaceVertexFrictionConstraint¶
Bases:
FaceVertexCandidate
,FrictionConstraint
Public Data Attributes:
Inherited from
FaceVertexCandidate
ID of the face
ID of the vertex
Inherited from
FrictionConstraint
Contact force magnitude
Coefficient of friction
Weight
Gradient of weight with respect to all DOF
Barycentric coordinates of the closest point(s)
Tangent basis of the contact (max size 3×2)
Public Methods:
__init__
(*args, **kwargs)Overloaded function.
Inherited from
FaceVertexCandidate
__init__
(self, face_id, vertex_id)__str__
(self)__repr__
(self)num_vertices
(self)vertex_ids
(self, edges, faces)ccd
(self, vertices_t0, vertices_t1, edges, faces)print_ccd_query
(self, vertices_t0, ...)__eq__
(self, other)__ne__
(self, other)__lt__
(self, other)Compare FaceVertexCandidate for sorting.
Inherited from
FrictionConstraint
__init__
(*args, **kwargs)compute_potential
(self, velocities, edges, ...)Compute the friction dissapative potential.
compute_potential_gradient
(self, velocities, ...)Compute the friction dissapative potential gradient wrt velocities.
compute_potential_hessian
(self, velocities, ...)Compute the friction dissapative potential hessian wrt velocities.
compute_force
(self, rest_positions, ...[, ...])Compute the friction force.
compute_force_jacobian
(self, rest_positions, ...)Compute the friction force Jacobian.
Inherited from
CollisionStencil
__init__
(*args, **kwargs)num_vertices
(self)Get the number of vertices in the contact stencil.
vertex_ids
(self, edges, faces)Get the vertex IDs of the contact stencil.
vertices
(self, vertices, edges, faces)Get the vertex attributes of the contact stencil.
dof
(self, X, edges, faces)Select this stencil's DOF from the full matrix of DOF.
compute_distance
(self, vertices, edges, faces)Compute the distance of the stencil.
compute_distance_gradient
(self, vertices, ...)Compute the distance gradient of the stencil w.r.t.
compute_distance_hessian
(self, vertices, ...)Compute the distance Hessian of the stencil w.r.t.
Inherited from
ContinuousCollisionCandidate
__init__
(*args, **kwargs)ccd
(self, vertices_t0, vertices_t1, edges, faces)Perform narrow-phase CCD on the candidate.
print_ccd_query
(self, vertices_t0, ...)Inherited from
pybind11_object
-
__annotations__ =
{}
¶
- __init__(*args, **kwargs)¶
Overloaded function.
__init__(self: ipctk.FaceVertexFrictionConstraint, constraint: ipctk.FaceVertexConstraint) -> None
__init__(self: ipctk.FaceVertexFrictionConstraint, constraint: ipctk.FaceVertexConstraint, vertices: numpy.ndarray[numpy.float64[m, n]], edges: numpy.ndarray[numpy.int32[m, n]], faces: numpy.ndarray[numpy.int32[m, n]], dhat: float, barrier_stiffness: float) -> None
-
__module__ =
'ipctk'
¶
-
__annotations__ =
Smooth Mollifier¶
- ipctk.f0_SF(s: float, epsv: float) float ¶
Smooth friction mollifier function.
\[f_0(s)= \begin{cases} -\frac{s^3}{3\epsilon_v^2} + \frac{s^2}{\epsilon_v} + \frac{\epsilon_v}{3}, & |s| < \epsilon_v \newline s, & |s| \geq \epsilon_v \end{cases}\]
- ipctk.f1_SF_over_x(s: float, epsv: float) float ¶
Compute the derivative of f0_SF divided by s (\(\frac{f_0'(s)}{s}\)).
\[f_1(s) = f_0'(s) = \begin{cases} -\frac{s^2}{\epsilon_v^2}+\frac{2 s}{\epsilon_v}, & |s| < \epsilon_v \newline 1, & |s| \geq \epsilon_v \end{cases}\]\[\frac{f_1(s)}{s} = \begin{cases} -\frac{s}{\epsilon_v^2}+\frac{2}{\epsilon_v}, & |s| < \epsilon_v \newline \frac{1}{s}, & |s| \geq \epsilon_v \end{cases}\]
- ipctk.df1_x_minus_f1_over_x3(s: float, epsv: float) float ¶
The derivative of f1 times s minus f1 all divided by s cubed.
\[\frac{f_1'(s) s - f_1(s)}{s^3} = \begin{cases} -\frac{1}{s \epsilon_v^2}, & |s| < \epsilon_v \newline -\frac{1}{s^3}, & |s| \geq \epsilon_v \end{cases}\]
Normal Force Magnitude¶
-
ipctk.compute_normal_force_magnitude(distance_squared: float, dhat: float, barrier_stiffness: float, dmin: float =
0
) float ¶
-
ipctk.compute_normal_force_magnitude_gradient(distance_squared: float, distance_squared_gradient: numpy.ndarray[numpy.float64[m, 1]], dhat: float, barrier_stiffness: float, dmin: float =
0
) numpy.ndarray[numpy.float64[m, 1]] ¶
Tangent Basis¶
- ipctk.point_point_tangent_basis(p0: numpy.ndarray[numpy.float64[m, 1]], p1: numpy.ndarray[numpy.float64[m, 1]]) numpy.ndarray[numpy.float64[m, n]] ¶
Compute a basis for the space tangent to the point-point pair.
- ipctk.point_edge_tangent_basis(p: numpy.ndarray[numpy.float64[m, 1]], e0: numpy.ndarray[numpy.float64[m, 1]], e1: numpy.ndarray[numpy.float64[m, 1]]) numpy.ndarray[numpy.float64[m, n]] ¶
- ipctk.edge_edge_tangent_basis(ea0: numpy.ndarray[numpy.float64[3, 1]], ea1: numpy.ndarray[numpy.float64[3, 1]], eb0: numpy.ndarray[numpy.float64[3, 1]], eb1: numpy.ndarray[numpy.float64[3, 1]]) numpy.ndarray[numpy.float64[3, 2]] ¶
Compute a basis for the space tangent to the edge-edge pair.
- ipctk.point_triangle_tangent_basis(p: numpy.ndarray[numpy.float64[3, 1]], t0: numpy.ndarray[numpy.float64[3, 1]], t1: numpy.ndarray[numpy.float64[3, 1]], t2: numpy.ndarray[numpy.float64[3, 1]]) numpy.ndarray[numpy.float64[3, 2]] ¶
Compute a basis for the space tangent to the point-triangle pair.
\[\begin{bmatrix} \frac{t_1 - t_0}{\|t_1 - t_0\|} & \frac{((t_1 - t_0)\times(t_2 - t_0)) \times(t_1 - t_0)}{\|((t_1 - t_0)\times(t_2 - t_0))\times(t_1 - t_0)\|} \end{bmatrix}\]
Tangent Basis Jacobians¶
- ipctk.point_point_tangent_basis_jacobian(p0: numpy.ndarray[numpy.float64[m, 1]], p1: numpy.ndarray[numpy.float64[m, 1]]) numpy.ndarray[numpy.float64[m, n]] ¶
- ipctk.point_edge_tangent_basis_jacobian(p: numpy.ndarray[numpy.float64[m, 1]], e0: numpy.ndarray[numpy.float64[m, 1]], e1: numpy.ndarray[numpy.float64[m, 1]]) numpy.ndarray[numpy.float64[m, n]] ¶
- ipctk.edge_edge_tangent_basis_jacobian(ea0: numpy.ndarray[numpy.float64[3, 1]], ea1: numpy.ndarray[numpy.float64[3, 1]], eb0: numpy.ndarray[numpy.float64[3, 1]], eb1: numpy.ndarray[numpy.float64[3, 1]]) numpy.ndarray[numpy.float64[36, 2]] ¶
- ipctk.point_triangle_tangent_basis_jacobian(p: numpy.ndarray[numpy.float64[3, 1]], t0: numpy.ndarray[numpy.float64[3, 1]], t1: numpy.ndarray[numpy.float64[3, 1]], t2: numpy.ndarray[numpy.float64[3, 1]]) numpy.ndarray[numpy.float64[36, 2]] ¶
Relative Velocity¶
- ipctk.point_point_relative_velocity(dp0: numpy.ndarray[numpy.float64[m, 1]], dp1: numpy.ndarray[numpy.float64[m, 1]]) numpy.ndarray[numpy.float64[m, 1]] ¶
Compute the relative velocity of two points
- ipctk.point_edge_relative_velocity(dp: numpy.ndarray[numpy.float64[m, 1]], de0: numpy.ndarray[numpy.float64[m, 1]], de1: numpy.ndarray[numpy.float64[m, 1]], alpha: float) numpy.ndarray[numpy.float64[m, 1]] ¶
Compute the relative velocity of a point and an edge
- ipctk.edge_edge_relative_velocity(dea0: numpy.ndarray[numpy.float64[3, 1]], dea1: numpy.ndarray[numpy.float64[3, 1]], deb0: numpy.ndarray[numpy.float64[3, 1]], deb1: numpy.ndarray[numpy.float64[3, 1]], coords: numpy.ndarray[numpy.float64[2, 1]]) numpy.ndarray[numpy.float64[3, 1]] ¶
Compute the relative velocity of the edges.
- Parameters:¶
- dea0: numpy.ndarray[numpy.float64[3, 1]]¶
Velocity of the first endpoint of the first edge
- dea1: numpy.ndarray[numpy.float64[3, 1]]¶
Velocity of the second endpoint of the first edge
- deb0: numpy.ndarray[numpy.float64[3, 1]]¶
Velocity of the first endpoint of the second edge
- deb1: numpy.ndarray[numpy.float64[3, 1]]¶
Velocity of the second endpoint of the second edge
- coords: numpy.ndarray[numpy.float64[2, 1]]¶
Two parametric coordinates of the closest points on the edges
- Returns:¶
The relative velocity of the edges
- ipctk.point_triangle_relative_velocity(dp: numpy.ndarray[numpy.float64[3, 1]], dt0: numpy.ndarray[numpy.float64[3, 1]], dt1: numpy.ndarray[numpy.float64[3, 1]], dt2: numpy.ndarray[numpy.float64[3, 1]], coords: numpy.ndarray[numpy.float64[2, 1]]) numpy.ndarray[numpy.float64[3, 1]] ¶
Compute the relative velocity of the point to the triangle.
- Parameters:¶
- dp: numpy.ndarray[numpy.float64[3, 1]]¶
Velocity of the point
- dt0: numpy.ndarray[numpy.float64[3, 1]]¶
Velocity of the first vertex of the triangle
- dt1: numpy.ndarray[numpy.float64[3, 1]]¶
Velocity of the second vertex of the triangle
- dt2: numpy.ndarray[numpy.float64[3, 1]]¶
Velocity of the third vertex of the triangle
- coords: numpy.ndarray[numpy.float64[2, 1]]¶
Baricentric coordinates of the closest point on the triangle
- Returns:¶
The relative velocity of the point to the triangle
Relative Velocity as Multiplier Matricies¶
- ipctk.point_point_relative_velocity_matrix(dim: int) numpy.ndarray[numpy.float64[m, n]] ¶
Compute the relative velocity premultiplier matrix
- ipctk.point_edge_relative_velocity_matrix(dim: int, alpha: float) numpy.ndarray[numpy.float64[m, n]] ¶
- ipctk.edge_edge_relative_velocity_matrix(dim: int, coords: numpy.ndarray[numpy.float64[2, 1]]) numpy.ndarray[numpy.float64[m, n]] ¶
- ipctk.point_triangle_relative_velocity_matrix(dim: int, coords: numpy.ndarray[numpy.float64[2, 1]]) numpy.ndarray[numpy.float64[m, n]] ¶
Relative Velocity Matrix Jacobians¶
- ipctk.point_point_relative_velocity_matrix_jacobian(dim: int) numpy.ndarray[numpy.float64[m, n]] ¶
Compute the jacobian of the relative velocity premultiplier matrix
- ipctk.point_edge_relative_velocity_matrix_jacobian(dim: int, alpha: float) numpy.ndarray[numpy.float64[m, n]] ¶