Candidates¶
Candidates¶
- class Candidates;¶
A class for storing and managing collision candidates.
Public Functions¶
- Candidates() = default;¶
-
void build(const CollisionMesh& mesh,
Eigen::ConstRef<Eigen::MatrixXd> vertices,
const double inflation_radius = 0,
BroadPhase* broad_phase = nullptr);¶ Initialize the set of discrete collision detection candidates.
- Parameters:¶
- const CollisionMesh &mesh¶
The surface of the collision mesh.
- Eigen::ConstRef<Eigen::MatrixXd> vertices¶
Surface vertex positions (rowwise).
- const double inflation_radius = 0¶
Amount to inflate the bounding boxes.
- BroadPhase *broad_phase = nullptr¶
Broad phase method to use.
-
void build(const CollisionMesh& mesh,
Eigen::ConstRef<Eigen::MatrixXd> vertices_t0,
Eigen::ConstRef<Eigen::MatrixXd> vertices_t1,
const double inflation_radius = 0,
BroadPhase* broad_phase = nullptr);¶ Initialize the set of continuous collision detection candidates.
Note
Assumes the trajectory is linear.
- Parameters:¶
- const CollisionMesh &mesh¶
The surface of the collision mesh.
- Eigen::ConstRef<Eigen::MatrixXd> vertices_t0¶
Surface vertex starting positions (rowwise).
- Eigen::ConstRef<Eigen::MatrixXd> vertices_t1¶
Surface vertex ending positions (rowwise).
- const double inflation_radius = 0¶
Amount to inflate the bounding boxes.
- BroadPhase *broad_phase = nullptr¶
Broad phase method to use.
- size_t size() const;¶
Get the number of collision candidates.
- Returns:¶
The number of collision candidates.
- bool empty() const;¶
Check if there are no collision candidates.
- Returns:¶
True if there are no collision candidates, false otherwise.
- void clear();¶
Clear all collision candidates.
- CollisionStencil& operator[](size_t i);¶
Get a collision stencil by index.
- const CollisionStencil& operator[](size_t i) const;¶
Get a collision stencil by index.
-
bool is_step_collision_free(const CollisionMesh& mesh,
Eigen::ConstRef<Eigen::MatrixXd> vertices_t0,
Eigen::ConstRef<Eigen::MatrixXd> vertices_t1,
const double min_distance = 0.0,
const NarrowPhaseCCD& narrow_phase_ccd
= DEFAULT_NARROW_PHASE_CCD) const;¶ Determine if the step is collision free from the set of candidates.
Note
Assumes the trajectory is linear.
- Parameters:¶
- const CollisionMesh &mesh¶
The collision mesh.
- Eigen::ConstRef<Eigen::MatrixXd> vertices_t0¶
Surface vertex starting positions (rowwise).
- Eigen::ConstRef<Eigen::MatrixXd> vertices_t1¶
Surface vertex ending positions (rowwise).
- const double min_distance = 0.0¶
The minimum distance allowable between any two elements.
- const NarrowPhaseCCD &narrow_phase_ccd = DEFAULT_NARROW_PHASE_CCD¶
The narrow phase CCD algorithm to use.
- Returns:¶
True if any collisions occur.
-
double compute_collision_free_stepsize(const CollisionMesh& mesh,
Eigen::ConstRef<Eigen::MatrixXd> vertices_t0,
Eigen::ConstRef<Eigen::MatrixXd> vertices_t1,
const double min_distance = 0.0,
const NarrowPhaseCCD& narrow_phase_ccd
= DEFAULT_NARROW_PHASE_CCD) const;¶ Computes a maximal step size that is collision free using the set of collision candidates.
Note
Assumes the trajectory is linear.
- Parameters:¶
- const CollisionMesh &mesh¶
The collision mesh.
- Eigen::ConstRef<Eigen::MatrixXd> vertices_t0¶
Surface vertex starting positions (rowwise). Assumed to be intersection free.
- Eigen::ConstRef<Eigen::MatrixXd> vertices_t1¶
Surface vertex ending positions (rowwise).
- const double min_distance = 0.0¶
The minimum distance allowable between any two elements.
- const NarrowPhaseCCD &narrow_phase_ccd = DEFAULT_NARROW_PHASE_CCD¶
The narrow phase CCD algorithm to use.
- Returns:¶
A step-size \(\in [0, 1]\) that is collision free. A value of 1.0 if a full step and 0.0 is no step.
-
double compute_noncandidate_conservative_stepsize(
const CollisionMesh& mesh,
Eigen::ConstRef<Eigen::MatrixXd> displacements,
const double dhat) const;¶ Computes a conservative bound on the largest-feasible step size for surface primitives not in collision.
-
double compute_cfl_stepsize(const CollisionMesh& mesh,
Eigen::ConstRef<Eigen::MatrixXd> vertices_t0,
Eigen::ConstRef<Eigen::MatrixXd> vertices_t1, const double dhat,
const double min_distance = 0.0,
BroadPhase* broad_phase = nullptr,
const NarrowPhaseCCD& narrow_phase_ccd
= DEFAULT_NARROW_PHASE_CCD) const;¶ Computes a CFL-inspired CCD maximum step step size.
- Parameters:¶
- const CollisionMesh &mesh¶
The collision mesh.
- Eigen::ConstRef<Eigen::MatrixXd> vertices_t0¶
Surface vertex starting positions (rowwise).
- Eigen::ConstRef<Eigen::MatrixXd> vertices_t1¶
Surface vertex ending positions (rowwise).
- const double dhat¶
Barrier activation distance.
- const double min_distance = 0.0¶
The minimum distance allowable between any two elements.
- BroadPhase *broad_phase = nullptr¶
The broad phase algorithm to use.
- const NarrowPhaseCCD &narrow_phase_ccd = DEFAULT_NARROW_PHASE_CCD¶
The narrow phase CCD algorithm to use.
- Returns:¶
A step-size \(\in [0, 1]\) that is collision free.
-
Eigen::VectorXd compute_per_vertex_safe_distances(
const CollisionMesh& mesh,
Eigen::ConstRef<Eigen::MatrixXd> vertices,
const double inflation_radius,
const double min_distance = 0.0) const;¶ Compute the maximum distance every vertex can move (independently) without colliding with any other element.
Note
Cap the value at the inflation radius used to build the candidates.
- Parameters:¶
- const CollisionMesh &mesh¶
The collision mesh.
- Eigen::ConstRef<Eigen::MatrixXd> vertices¶
Collision mesh vertex positions (rowwise).
- const double inflation_radius¶
The inflation radius used to build the candidates.
- const double min_distance = 0.0¶
The minimum distance allowable between any two elements.
- Returns:¶
A vector of minimum distances, one for each vertex.
-
std::vector<VertexVertexCandidate> edge_vertex_to_vertex_vertex(
const CollisionMesh& mesh,
Eigen::ConstRef<Eigen::MatrixXd> vertices,
const std::function<bool(double)>& is_active
= default_is_active) const;¶ Convert edge-vertex candidates to vertex-vertex candidates.
- Parameters:¶
- const CollisionMesh &mesh¶
The collision mesh.
- Eigen::ConstRef<Eigen::MatrixXd> vertices¶
Collision mesh vertex positions (rowwise).
- const std::function<bool(double)> &is_active = default_is_active¶
(Optional) Function to determine if a candidate is active.
- Returns:¶
Vertex-vertex candidates derived from edge-vertex candidates.
-
std::vector<VertexVertexCandidate> face_vertex_to_vertex_vertex(
const CollisionMesh& mesh,
Eigen::ConstRef<Eigen::MatrixXd> vertices,
const std::function<bool(double)>& is_active
= default_is_active) const;¶ Convert face-vertex candidates to vertex-vertex candidates.
- Parameters:¶
- const CollisionMesh &mesh¶
The collision mesh.
- Eigen::ConstRef<Eigen::MatrixXd> vertices¶
Collision mesh vertex positions (rowwise).
- const std::function<bool(double)> &is_active = default_is_active¶
(Optional) Function to determine if a candidate is active.
- Returns:¶
Vertex-vertex candidates derived from face-vertex candidates.
-
std::vector<EdgeVertexCandidate> face_vertex_to_edge_vertex(
const CollisionMesh& mesh,
Eigen::ConstRef<Eigen::MatrixXd> vertices,
const std::function<bool(double)>& is_active
= default_is_active) const;¶ Convert face-vertex candidates to edge-vertex candidates.
- Parameters:¶
- const CollisionMesh &mesh¶
The collision mesh.
- Eigen::ConstRef<Eigen::MatrixXd> vertices¶
Collision mesh vertex positions (rowwise).
- const std::function<bool(double)> &is_active = default_is_active¶
(Optional) Function to determine if a candidate is active.
- Returns:¶
Edge-vertex candidates derived from face-vertex candidates.
-
std::vector<EdgeVertexCandidate> edge_edge_to_edge_vertex(
const CollisionMesh& mesh,
Eigen::ConstRef<Eigen::MatrixXd> vertices,
const std::function<bool(double)>& is_active
= default_is_active) const;¶ Convert edge-edge candidates to edge-vertex candidates.
- Parameters:¶
- const CollisionMesh &mesh¶
The collision mesh.
- Eigen::ConstRef<Eigen::MatrixXd> vertices¶
Collision mesh vertex positions (rowwise).
- const std::function<bool(double)> &is_active = default_is_active¶
(Optional) Function to determine if a candidate is active.
- Returns:¶
Edge-vertex candidates derived from edge-edge candidates.
-
bool save_obj(const std::string& filename,
Eigen::ConstRef<Eigen::MatrixXd> vertices,
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const;¶ Save the collision candidates to an OBJ file.
- Parameters:¶
- const std::string &filename¶
The name of the file to save the candidates to.
- Eigen::ConstRef<Eigen::MatrixXd> vertices¶
Collision mesh vertex positions (rowwise).
- Eigen::ConstRef<Eigen::MatrixXi> edges¶
Collision mesh edge indices (rowwise).
- Eigen::ConstRef<Eigen::MatrixXi> faces¶
Collision mesh face indices (rowwise).
- Returns:¶
True if the file was saved successfully.
Public Members¶
- std::vector<VertexVertexCandidate> vv_candidates;¶
- std::vector<EdgeVertexCandidate> ev_candidates;¶
- std::vector<EdgeEdgeCandidate> ee_candidates;¶
- std::vector<FaceVertexCandidate> fv_candidates;¶
Collision Stencil¶
- class CollisionStencil;¶
Inheritance diagram for ipc::CollisionStencil:
A stencil representing a collision between at most four vertices.
Subclassed by ipc::EdgeEdgeCandidate, ipc::EdgeVertexCandidate, ipc::FaceVertexCandidate, ipc::NormalCollision, ipc::TangentialCollision, ipc::VertexVertexCandidate
Public Functions¶
- virtual ~CollisionStencil() = default;¶
- virtual int num_vertices() const = 0;¶
Get the number of vertices in the collision stencil.
-
virtual std::array<index_t, STENCIL_SIZE> vertex_ids(
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const
= 0;¶ Get the vertex IDs of the collision stencil.
- Parameters:¶
- Returns:¶
The vertex IDs of the collision stencil. Elements i > num_vertices() are -1.
-
inline std::array<VectorMax3d, STENCIL_SIZE> vertices(
Eigen::ConstRef<Eigen::MatrixXd> vertices,
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const;¶ Get the vertex attributes of the collision stencil.
-
inline VectorMax12d dof(Eigen::ConstRef<Eigen::MatrixXd> X,
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const;¶ Select this stencil’s DOF from the full matrix of DOF.
-
inline double compute_distance(
Eigen::ConstRef<Eigen::MatrixXd> vertices,
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const;¶ Compute the distance of the stencil.
-
inline VectorMax12d compute_distance_gradient(
Eigen::ConstRef<Eigen::MatrixXd> vertices,
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const;¶ Compute the distance gradient of the stencil w.r.t.
the stencil’s vertex positions.
-
inline MatrixMax12d compute_distance_hessian(
Eigen::ConstRef<Eigen::MatrixXd> vertices,
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const;¶ Compute the distance Hessian of the stencil w.r.t.
the stencil’s vertex positions.
-
inline VectorMax4d compute_coefficients(
Eigen::ConstRef<Eigen::MatrixXd> vertices,
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const;¶ Compute the coefficients of the stencil s.t.
\(d(x) = \|\sum c_i \mathbf{x}_i\|^2\).
-
inline VectorMax3d compute_normal(
Eigen::ConstRef<Eigen::MatrixXd> vertices,
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces,
const bool flip_if_negative = true,
double* sign = nullptr) const;¶ Compute the normal of the stencil.
- Parameters:¶
- Eigen::ConstRef<Eigen::MatrixXd> vertices¶
Collision mesh vertices
- Eigen::ConstRef<Eigen::MatrixXi> edges¶
Collision mesh edges
- Eigen::ConstRef<Eigen::MatrixXi> faces¶
Collision mesh faces
- const bool flip_if_negative = true¶
If true, flip the normal if the point is on the negative side.
- double *sign = nullptr¶
If not nullptr, set to the sign of the normal before any flipping.
- Returns:¶
Normal of the stencil.
-
inline MatrixMax<double, 3, 12> compute_normal_jacobian(
Eigen::ConstRef<Eigen::MatrixXd> vertices,
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces,
const bool flip_if_negative = true) const;¶ Compute the Jacobian of the normal of the stencil.
-
virtual double compute_distance(
Eigen::ConstRef<VectorMax12d> positions) const
= 0;¶ Compute the distance of the stencil.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Distance of the stencil.
-
virtual VectorMax12d compute_distance_gradient(
Eigen::ConstRef<VectorMax12d> positions) const
= 0;¶ Compute the distance gradient of the stencil w.r.t.
the stencil’s vertex positions.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Distance gradient of the stencil w.r.t. the stencil’s vertex positions.
-
virtual MatrixMax12d compute_distance_hessian(
Eigen::ConstRef<VectorMax12d> positions) const
= 0;¶ Compute the distance Hessian of the stencil w.r.t.
the stencil’s vertex positions.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Distance Hessian of the stencil w.r.t. the stencil’s vertex positions.
-
virtual VectorMax4d compute_coefficients(
Eigen::ConstRef<VectorMax12d> positions) const
= 0;¶ Compute the coefficients of the stencil s.t.
d(x) = ‖∑ cᵢ xᵢ‖².
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Coefficients of the stencil.
-
VectorMax3d compute_normal(Eigen::ConstRef<VectorMax12d> positions,
bool flip_if_negative = true, double* sign = nullptr) const;¶ Compute the normal of the stencil.
-
MatrixMax<double, 3, 12> compute_normal_jacobian(
Eigen::ConstRef<VectorMax12d> positions,
bool flip_if_negative = true) const;¶ Compute the Jacobian of the normal of the stencil.
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- bool flip_if_negative = true¶
If true, flip the normal if the point is on the negative side.
- Returns:¶
Jacobian of the normal of the stencil.
-
virtual bool ccd(Eigen::ConstRef<VectorMax12d> vertices_t0,
Eigen::ConstRef<VectorMax12d> vertices_t1, double& toi,
const double min_distance = 0.0, const double tmax = 1.0,
const NarrowPhaseCCD& narrow_phase_ccd
= DEFAULT_NARROW_PHASE_CCD) const
= 0;¶ Perform narrow-phase CCD on the candidate.
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> vertices_t0¶
[in] Stencil vertices at the start of the time step.
- Eigen::ConstRef<VectorMax12d> vertices_t1¶
[in] Stencil vertices at the end of the time step.
- double &toi¶
[out] Computed time of impact (normalized).
- const double min_distance = 0.0¶
[in] Minimum separation distance between primitives.
- const double tmax = 1.0¶
[in] Maximum time (normalized) to look for collisions.
- const NarrowPhaseCCD &narrow_phase_ccd = DEFAULT_NARROW_PHASE_CCD¶
[in] The narrow phase CCD algorithm to use.
- Returns:¶
If the candidate had a collision over the time interval.
-
std::ostream& write_ccd_query(std::ostream& out,
Eigen::ConstRef<VectorMax12d> vertices_t0,
Eigen::ConstRef<VectorMax12d> vertices_t1) const;¶ Write the CCD query to a stream.
- Parameters:¶
- std::ostream &out¶
Stream to write to.
- Eigen::ConstRef<VectorMax12d> vertices_t0¶
Stencil vertices at the start of the time step.
- Eigen::ConstRef<VectorMax12d> vertices_t1¶
Stencil vertices at the end of the time step.
- Returns:¶
The stream.
Public Static Attributes¶
- static constexpr int STENCIL_SIZE = 4;¶
The maximum number of vertices in a collision stencil.
Protected Functions¶
-
virtual VectorMax3d compute_unnormalized_normal(
Eigen::ConstRef<VectorMax12d> positions) const
= 0;¶ Compute the unnormalized normal of the stencil.
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Unnormalized normal of the stencil.
-
virtual MatrixMax<double, 3, 12>
compute_unnormalized_normal_jacobian(
Eigen::ConstRef<VectorMax12d> positions) const
= 0;¶ Compute the Jacobian of the unnormalized normal of the stencil.
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Jacobian of the unnormalized normal of the stencil.
Vertex-Vertex Candidate¶
- class VertexVertexCandidate : public virtual ipc::CollisionStencil;¶
Inheritance diagram for ipc::VertexVertexCandidate:
Collaboration diagram for ipc::VertexVertexCandidate:
A candidate for vertex-vertex collision detection.
Subclassed by ipc::VertexVertexNormalCollision, ipc::VertexVertexTangentialCollision
Public Functions¶
- VertexVertexCandidate(index_t vertex0_id, index_t vertex1_id);¶
- inline virtual int num_vertices() const override;¶
Get the number of vertices in the collision stencil.
-
inline virtual std::array<index_t, 4> vertex_ids(
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const override;¶ Get the indices of the vertices.
-
virtual double compute_distance(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the distance of the stencil.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Distance of the stencil.
-
virtual VectorMax12d compute_distance_gradient(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the distance gradient of the stencil w.r.t.
the stencil’s vertex positions.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Distance gradient of the stencil w.r.t. the stencil’s vertex positions.
-
virtual MatrixMax12d compute_distance_hessian(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the distance Hessian of the stencil w.r.t.
the stencil’s vertex positions.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Distance Hessian of the stencil w.r.t. the stencil’s vertex positions.
-
virtual VectorMax4d compute_coefficients(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the coefficients of the stencil s.t.
d(x) = ‖∑ cᵢ xᵢ‖².
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Coefficients of the stencil.
-
virtual bool ccd(Eigen::ConstRef<VectorMax12d> vertices_t0,
Eigen::ConstRef<VectorMax12d> vertices_t1, double& toi,
const double min_distance = 0.0, const double tmax = 1.0,
const NarrowPhaseCCD& narrow_phase_ccd
= DEFAULT_NARROW_PHASE_CCD) const override;¶ Perform narrow-phase CCD on the candidate.
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> vertices_t0¶
[in] Stencil vertices at the start of the time step.
- Eigen::ConstRef<VectorMax12d> vertices_t1¶
[in] Stencil vertices at the end of the time step.
- double &toi¶
[out] Computed time of impact (normalized).
- const double min_distance = 0.0¶
[in] Minimum separation distance between primitives.
- const double tmax = 1.0¶
[in] Maximum time (normalized) to look for collisions.
- const NarrowPhaseCCD &narrow_phase_ccd = DEFAULT_NARROW_PHASE_CCD¶
[in] The narrow phase CCD algorithm to use.
- Returns:¶
If the candidate had a collision over the time interval.
- bool operator==(const VertexVertexCandidate& other) const;¶
- bool operator!=(const VertexVertexCandidate& other) const;¶
- bool operator<(const VertexVertexCandidate& other) const;¶
Compare EdgeVertexCandidates for sorting.
-
inline VectorMax4d compute_coefficients(
Eigen::ConstRef<Eigen::MatrixXd> vertices,
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const;¶ Compute the coefficients of the stencil s.t.
\(d(x) = \|\sum c_i \mathbf{x}_i\|^2\).
-
virtual VectorMax4d compute_coefficients(
Eigen::ConstRef<VectorMax12d> positions) const
= 0; Compute the coefficients of the stencil s.t.
d(x) = ‖∑ cᵢ xᵢ‖².
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Coefficients of the stencil.
-
inline double compute_distance(
Eigen::ConstRef<Eigen::MatrixXd> vertices,
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const;¶ Compute the distance of the stencil.
-
virtual double compute_distance(
Eigen::ConstRef<VectorMax12d> positions) const
= 0; Compute the distance of the stencil.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Distance of the stencil.
-
inline VectorMax12d compute_distance_gradient(
Eigen::ConstRef<Eigen::MatrixXd> vertices,
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const;¶ Compute the distance gradient of the stencil w.r.t.
the stencil’s vertex positions.
-
virtual VectorMax12d compute_distance_gradient(
Eigen::ConstRef<VectorMax12d> positions) const
= 0; Compute the distance gradient of the stencil w.r.t.
the stencil’s vertex positions.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Distance gradient of the stencil w.r.t. the stencil’s vertex positions.
-
inline MatrixMax12d compute_distance_hessian(
Eigen::ConstRef<Eigen::MatrixXd> vertices,
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const;¶ Compute the distance Hessian of the stencil w.r.t.
the stencil’s vertex positions.
-
virtual MatrixMax12d compute_distance_hessian(
Eigen::ConstRef<VectorMax12d> positions) const
= 0; Compute the distance Hessian of the stencil w.r.t.
the stencil’s vertex positions.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Distance Hessian of the stencil w.r.t. the stencil’s vertex positions.
Public Members¶
- index_t vertex0_id;¶
ID of the first vertex.
- index_t vertex1_id;¶
ID of the second vertex.
Protected Functions¶
-
virtual VectorMax3d compute_unnormalized_normal(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the unnormalized normal of the stencil.
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Unnormalized normal of the stencil.
-
virtual MatrixMax<double, 3, 12>
compute_unnormalized_normal_jacobian(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the Jacobian of the unnormalized normal of the stencil.
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Jacobian of the unnormalized normal of the stencil.
Edge-Vertex Candidate¶
- class EdgeVertexCandidate : public virtual ipc::CollisionStencil;¶
Inheritance diagram for ipc::EdgeVertexCandidate:
Collaboration diagram for ipc::EdgeVertexCandidate:
A candidate for edge-vertex collision detection.
Subclassed by ipc::EdgeVertexNormalCollision, ipc::EdgeVertexTangentialCollision
Public Functions¶
- inline virtual int num_vertices() const override;¶
Get the number of vertices in the collision stencil.
-
inline virtual std::array<index_t, 4> vertex_ids(
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const override;¶ Get the vertex IDs for the edge-vertex pair.
-
virtual double compute_distance(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the distance of the stencil.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Distance of the stencil.
-
virtual VectorMax12d compute_distance_gradient(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the distance gradient of the stencil w.r.t.
the stencil’s vertex positions.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Distance gradient of the stencil w.r.t. the stencil’s vertex positions.
-
virtual MatrixMax12d compute_distance_hessian(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the distance Hessian of the stencil w.r.t.
the stencil’s vertex positions.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Distance Hessian of the stencil w.r.t. the stencil’s vertex positions.
-
virtual VectorMax4d compute_coefficients(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the coefficients of the stencil s.t.
d(x) = ‖∑ cᵢ xᵢ‖².
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Coefficients of the stencil.
-
virtual bool ccd(Eigen::ConstRef<VectorMax12d> vertices_t0,
Eigen::ConstRef<VectorMax12d> vertices_t1, double& toi,
const double min_distance = 0.0, const double tmax = 1.0,
const NarrowPhaseCCD& narrow_phase_ccd
= DEFAULT_NARROW_PHASE_CCD) const override;¶ Perform narrow-phase CCD on the candidate.
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> vertices_t0¶
[in] Stencil vertices at the start of the time step.
- Eigen::ConstRef<VectorMax12d> vertices_t1¶
[in] Stencil vertices at the end of the time step.
- double &toi¶
[out] Computed time of impact (normalized).
- const double min_distance = 0.0¶
[in] Minimum separation distance between primitives.
- const double tmax = 1.0¶
[in] Maximum time (normalized) to look for collisions.
- const NarrowPhaseCCD &narrow_phase_ccd = DEFAULT_NARROW_PHASE_CCD¶
[in] The narrow phase CCD algorithm to use.
- Returns:¶
If the candidate had a collision over the time interval.
- inline virtual PointEdgeDistanceType known_dtype() const;¶
- bool operator==(const EdgeVertexCandidate& other) const;¶
- bool operator!=(const EdgeVertexCandidate& other) const;¶
- bool operator<(const EdgeVertexCandidate& other) const;¶
Compare EdgeVertexCandidates for sorting.
-
inline VectorMax4d compute_coefficients(
Eigen::ConstRef<Eigen::MatrixXd> vertices,
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const;¶ Compute the coefficients of the stencil s.t.
\(d(x) = \|\sum c_i \mathbf{x}_i\|^2\).
-
virtual VectorMax4d compute_coefficients(
Eigen::ConstRef<VectorMax12d> positions) const
= 0; Compute the coefficients of the stencil s.t.
d(x) = ‖∑ cᵢ xᵢ‖².
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Coefficients of the stencil.
-
inline double compute_distance(
Eigen::ConstRef<Eigen::MatrixXd> vertices,
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const;¶ Compute the distance of the stencil.
-
virtual double compute_distance(
Eigen::ConstRef<VectorMax12d> positions) const
= 0; Compute the distance of the stencil.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Distance of the stencil.
-
inline VectorMax12d compute_distance_gradient(
Eigen::ConstRef<Eigen::MatrixXd> vertices,
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const;¶ Compute the distance gradient of the stencil w.r.t.
the stencil’s vertex positions.
-
virtual VectorMax12d compute_distance_gradient(
Eigen::ConstRef<VectorMax12d> positions) const
= 0; Compute the distance gradient of the stencil w.r.t.
the stencil’s vertex positions.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Distance gradient of the stencil w.r.t. the stencil’s vertex positions.
-
inline MatrixMax12d compute_distance_hessian(
Eigen::ConstRef<Eigen::MatrixXd> vertices,
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const;¶ Compute the distance Hessian of the stencil w.r.t.
the stencil’s vertex positions.
-
virtual MatrixMax12d compute_distance_hessian(
Eigen::ConstRef<VectorMax12d> positions) const
= 0; Compute the distance Hessian of the stencil w.r.t.
the stencil’s vertex positions.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Distance Hessian of the stencil w.r.t. the stencil’s vertex positions.
Protected Functions¶
-
virtual VectorMax3d compute_unnormalized_normal(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the unnormalized normal of the stencil.
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Unnormalized normal of the stencil.
-
virtual MatrixMax<double, 3, 12>
compute_unnormalized_normal_jacobian(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the Jacobian of the unnormalized normal of the stencil.
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Jacobian of the unnormalized normal of the stencil.
Edge-Edge Candidate¶
- class EdgeEdgeCandidate : public virtual ipc::CollisionStencil;¶
Inheritance diagram for ipc::EdgeEdgeCandidate:
Collaboration diagram for ipc::EdgeEdgeCandidate:
A candidate for edge-edge collision detection.
Subclassed by ipc::EdgeEdgeNormalCollision, ipc::EdgeEdgeTangentialCollision
Public Functions¶
- inline virtual int num_vertices() const override;¶
Get the number of vertices in the collision stencil.
-
inline virtual std::array<index_t, 4> vertex_ids(
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const override;¶ Get the vertex IDs for the edge-edge pair.
-
virtual double compute_distance(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the distance of the stencil.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Distance of the stencil.
-
virtual VectorMax12d compute_distance_gradient(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the distance gradient of the stencil w.r.t.
the stencil’s vertex positions.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Distance gradient of the stencil w.r.t. the stencil’s vertex positions.
-
virtual MatrixMax12d compute_distance_hessian(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the distance Hessian of the stencil w.r.t.
the stencil’s vertex positions.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Distance Hessian of the stencil w.r.t. the stencil’s vertex positions.
-
virtual VectorMax4d compute_coefficients(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the coefficients of the stencil s.t.
d(x) = ‖∑ cᵢ xᵢ‖².
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Coefficients of the stencil.
-
virtual bool ccd(Eigen::ConstRef<VectorMax12d> vertices_t0,
Eigen::ConstRef<VectorMax12d> vertices_t1, double& toi,
const double min_distance = 0.0, const double tmax = 1.0,
const NarrowPhaseCCD& narrow_phase_ccd
= DEFAULT_NARROW_PHASE_CCD) const override;¶ Perform narrow-phase CCD on the candidate.
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> vertices_t0¶
[in] Stencil vertices at the start of the time step.
- Eigen::ConstRef<VectorMax12d> vertices_t1¶
[in] Stencil vertices at the end of the time step.
- double &toi¶
[out] Computed time of impact (normalized).
- const double min_distance = 0.0¶
[in] Minimum separation distance between primitives.
- const double tmax = 1.0¶
[in] Maximum time (normalized) to look for collisions.
- const NarrowPhaseCCD &narrow_phase_ccd = DEFAULT_NARROW_PHASE_CCD¶
[in] The narrow phase CCD algorithm to use.
- Returns:¶
If the candidate had a collision over the time interval.
- inline virtual EdgeEdgeDistanceType known_dtype() const;¶
- bool operator==(const EdgeEdgeCandidate& other) const;¶
- bool operator!=(const EdgeEdgeCandidate& other) const;¶
- bool operator<(const EdgeEdgeCandidate& other) const;¶
Compare EdgeEdgeCandidates for sorting.
-
inline VectorMax4d compute_coefficients(
Eigen::ConstRef<Eigen::MatrixXd> vertices,
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const;¶ Compute the coefficients of the stencil s.t.
\(d(x) = \|\sum c_i \mathbf{x}_i\|^2\).
-
virtual VectorMax4d compute_coefficients(
Eigen::ConstRef<VectorMax12d> positions) const
= 0; Compute the coefficients of the stencil s.t.
d(x) = ‖∑ cᵢ xᵢ‖².
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Coefficients of the stencil.
-
inline double compute_distance(
Eigen::ConstRef<Eigen::MatrixXd> vertices,
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const;¶ Compute the distance of the stencil.
-
virtual double compute_distance(
Eigen::ConstRef<VectorMax12d> positions) const
= 0; Compute the distance of the stencil.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Distance of the stencil.
-
inline VectorMax12d compute_distance_gradient(
Eigen::ConstRef<Eigen::MatrixXd> vertices,
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const;¶ Compute the distance gradient of the stencil w.r.t.
the stencil’s vertex positions.
-
virtual VectorMax12d compute_distance_gradient(
Eigen::ConstRef<VectorMax12d> positions) const
= 0; Compute the distance gradient of the stencil w.r.t.
the stencil’s vertex positions.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Distance gradient of the stencil w.r.t. the stencil’s vertex positions.
-
inline MatrixMax12d compute_distance_hessian(
Eigen::ConstRef<Eigen::MatrixXd> vertices,
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const;¶ Compute the distance Hessian of the stencil w.r.t.
the stencil’s vertex positions.
-
virtual MatrixMax12d compute_distance_hessian(
Eigen::ConstRef<VectorMax12d> positions) const
= 0; Compute the distance Hessian of the stencil w.r.t.
the stencil’s vertex positions.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Distance Hessian of the stencil w.r.t. the stencil’s vertex positions.
Protected Functions¶
-
virtual VectorMax3d compute_unnormalized_normal(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the unnormalized normal of the stencil.
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Unnormalized normal of the stencil.
-
virtual MatrixMax<double, 3, 12>
compute_unnormalized_normal_jacobian(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the Jacobian of the unnormalized normal of the stencil.
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Jacobian of the unnormalized normal of the stencil.
Edge-Face Candidate¶
- class EdgeFaceCandidate;¶
Candidate for intersection between edge and face.
Not included in Candidates because it is not a collision candidate.
Public Functions¶
- bool operator==(const EdgeFaceCandidate& other) const;¶
- bool operator!=(const EdgeFaceCandidate& other) const;¶
- bool operator<(const EdgeFaceCandidate& other) const;¶
Compare EdgeFaceCandidate for sorting.
Face-Vertex Candidate¶
- class FaceVertexCandidate : public virtual ipc::CollisionStencil;¶
Inheritance diagram for ipc::FaceVertexCandidate:
Collaboration diagram for ipc::FaceVertexCandidate:
A candidate for face-vertex collision detection.
Subclassed by ipc::FaceVertexNormalCollision, ipc::FaceVertexTangentialCollision
Public Functions¶
- inline virtual int num_vertices() const override;¶
Get the number of vertices in the collision stencil.
-
inline virtual std::array<index_t, 4> vertex_ids(
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const override;¶ Get the vertex IDs for the face-vertex pair.
-
virtual double compute_distance(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the distance of the stencil.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Distance of the stencil.
-
virtual VectorMax12d compute_distance_gradient(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the distance gradient of the stencil w.r.t.
the stencil’s vertex positions.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Distance gradient of the stencil w.r.t. the stencil’s vertex positions.
-
virtual MatrixMax12d compute_distance_hessian(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the distance Hessian of the stencil w.r.t.
the stencil’s vertex positions.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Distance Hessian of the stencil w.r.t. the stencil’s vertex positions.
-
virtual VectorMax4d compute_coefficients(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the coefficients of the stencil s.t.
d(x) = ‖∑ cᵢ xᵢ‖².
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Coefficients of the stencil.
-
virtual bool ccd(Eigen::ConstRef<VectorMax12d> vertices_t0,
Eigen::ConstRef<VectorMax12d> vertices_t1, double& toi,
const double min_distance = 0.0, const double tmax = 1.0,
const NarrowPhaseCCD& narrow_phase_ccd
= DEFAULT_NARROW_PHASE_CCD) const override;¶ Perform narrow-phase CCD on the candidate.
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> vertices_t0¶
[in] Stencil vertices at the start of the time step.
- Eigen::ConstRef<VectorMax12d> vertices_t1¶
[in] Stencil vertices at the end of the time step.
- double &toi¶
[out] Computed time of impact (normalized).
- const double min_distance = 0.0¶
[in] Minimum separation distance between primitives.
- const double tmax = 1.0¶
[in] Maximum time (normalized) to look for collisions.
- const NarrowPhaseCCD &narrow_phase_ccd = DEFAULT_NARROW_PHASE_CCD¶
[in] The narrow phase CCD algorithm to use.
- Returns:¶
If the candidate had a collision over the time interval.
- inline virtual PointTriangleDistanceType known_dtype() const;¶
- bool operator==(const FaceVertexCandidate& other) const;¶
- bool operator!=(const FaceVertexCandidate& other) const;¶
- bool operator<(const FaceVertexCandidate& other) const;¶
Compare FaceVertexCandidate for sorting.
-
inline VectorMax4d compute_coefficients(
Eigen::ConstRef<Eigen::MatrixXd> vertices,
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const;¶ Compute the coefficients of the stencil s.t.
\(d(x) = \|\sum c_i \mathbf{x}_i\|^2\).
-
virtual VectorMax4d compute_coefficients(
Eigen::ConstRef<VectorMax12d> positions) const
= 0; Compute the coefficients of the stencil s.t.
d(x) = ‖∑ cᵢ xᵢ‖².
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Coefficients of the stencil.
-
inline double compute_distance(
Eigen::ConstRef<Eigen::MatrixXd> vertices,
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const;¶ Compute the distance of the stencil.
-
virtual double compute_distance(
Eigen::ConstRef<VectorMax12d> positions) const
= 0; Compute the distance of the stencil.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Distance of the stencil.
-
inline VectorMax12d compute_distance_gradient(
Eigen::ConstRef<Eigen::MatrixXd> vertices,
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const;¶ Compute the distance gradient of the stencil w.r.t.
the stencil’s vertex positions.
-
virtual VectorMax12d compute_distance_gradient(
Eigen::ConstRef<VectorMax12d> positions) const
= 0; Compute the distance gradient of the stencil w.r.t.
the stencil’s vertex positions.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Distance gradient of the stencil w.r.t. the stencil’s vertex positions.
-
inline MatrixMax12d compute_distance_hessian(
Eigen::ConstRef<Eigen::MatrixXd> vertices,
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const;¶ Compute the distance Hessian of the stencil w.r.t.
the stencil’s vertex positions.
-
virtual MatrixMax12d compute_distance_hessian(
Eigen::ConstRef<VectorMax12d> positions) const
= 0; Compute the distance Hessian of the stencil w.r.t.
the stencil’s vertex positions.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Distance Hessian of the stencil w.r.t. the stencil’s vertex positions.
Protected Functions¶
-
virtual VectorMax3d compute_unnormalized_normal(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the unnormalized normal of the stencil.
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Unnormalized normal of the stencil.
-
virtual MatrixMax<double, 3, 12>
compute_unnormalized_normal_jacobian(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the Jacobian of the unnormalized normal of the stencil.
- Parameters:¶
- Eigen::ConstRef<VectorMax12d> positions¶
Stencil’s vertex positions.
- Returns:¶
Jacobian of the unnormalized normal of the stencil.