Tangent

Tangent Basis

MatrixMax<double, 3, 2> ipc::point_point_tangent_basis(
    
const Eigen::Ref<const VectorMax3d>p0,
    
const Eigen::Ref<const VectorMax3d>p1);

Compute a basis for the space tangent to the point-point pair.

Parameters:
const Eigen::Ref<const VectorMax3d> &p0

First point

const Eigen::Ref<const VectorMax3d> &p1

Second point

Returns:

A 3x2 matrix whose columns are the basis vectors.

MatrixMax<double, 3, 2> ipc::point_edge_tangent_basis(
    
const Eigen::Ref<const VectorMax3d>p,
    
const Eigen::Ref<const VectorMax3d>e0,
    
const Eigen::Ref<const VectorMax3d>e1);

Compute a basis for the space tangent to the point-edge pair.

Parameters:
const Eigen::Ref<const VectorMax3d> &p

Point

const Eigen::Ref<const VectorMax3d> &e0

First edge point

const Eigen::Ref<const VectorMax3d> &e1

Second edge point

Returns:

A 3x2 matrix whose columns are the basis vectors.

Eigen::Matrix<double, 3, 2> ipc::edge_edge_tangent_basis(
    
const Eigen::Ref<const Eigen::Vector3d>ea0,
    
const Eigen::Ref<const Eigen::Vector3d>ea1,
    
const Eigen::Ref<const Eigen::Vector3d>eb0,
    
const Eigen::Ref<const Eigen::Vector3d>eb1);

Compute a basis for the space tangent to the edge-edge pair.

Parameters:
const Eigen::Ref<const Eigen::Vector3d> &ea0

First point of the first edge

const Eigen::Ref<const Eigen::Vector3d> &ea1

Second point of the first edge

const Eigen::Ref<const Eigen::Vector3d> &eb0

First point of the second edge

const Eigen::Ref<const Eigen::Vector3d> &eb1

Second point of the second edge

Returns:

A 3x2 matrix whose columns are the basis vectors.

Eigen::Matrix<double, 3, 2> ipc::point_triangle_tangent_basis(
    
const Eigen::Ref<const Eigen::Vector3d>p,
    
const Eigen::Ref<const Eigen::Vector3d>t0,
    
const Eigen::Ref<const Eigen::Vector3d>t1,
    
const Eigen::Ref<const Eigen::Vector3d>t2);

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} \]

Parameters:
const Eigen::Ref<const Eigen::Vector3d> &p

Point

const Eigen::Ref<const Eigen::Vector3d> &t0

Triangle’s first vertex

const Eigen::Ref<const Eigen::Vector3d> &t1

Triangle’s second vertex

const Eigen::Ref<const Eigen::Vector3d> &t2

Triangle’s third vertex

Returns:

A 3x2 matrix whose columns are the basis vectors.

Tangent Basis Jacobians

MatrixMax<double, 18, 2> ipc::point_point_tangent_basis_jacobian(
    
const Eigen::Ref<const VectorMax3d>p0,
    
const Eigen::Ref<const VectorMax3d>p1);

Compute the Jacobian of the tangent basis for the point-point pair.

Parameters:
const Eigen::Ref<const VectorMax3d> &p0

First point

const Eigen::Ref<const VectorMax3d> &p1

Second point

Returns:

A 6*3x2 matrix whose columns are the basis vectors.

MatrixMax<double, 27, 2> ipc::point_edge_tangent_basis_jacobian(
    
const Eigen::Ref<const VectorMax3d>p,
    
const Eigen::Ref<const VectorMax3d>e0,
    
const Eigen::Ref<const VectorMax3d>e1);

Compute the Jacobian of the tangent basis for the point-edge pair.

Parameters:
const Eigen::Ref<const VectorMax3d> &p

Point

const Eigen::Ref<const VectorMax3d> &e0

First edge point

const Eigen::Ref<const VectorMax3d> &e1

Second edge point

Returns:

A 9*3x2 matrix whose columns are the basis vectors.

Eigen::Matrix<double, 36, 2> ipc::edge_edge_tangent_basis_jacobian(
    
const Eigen::Ref<const Eigen::Vector3d>ea0,
    
const Eigen::Ref<const Eigen::Vector3d>ea1,
    
const Eigen::Ref<const Eigen::Vector3d>eb0,
    
const Eigen::Ref<const Eigen::Vector3d>eb1);

Compute the Jacobian of the tangent basis for the edge-edge pair.

Parameters:
const Eigen::Ref<const Eigen::Vector3d> &ea0

First point of the first edge

const Eigen::Ref<const Eigen::Vector3d> &ea1

Second point of the first edge

const Eigen::Ref<const Eigen::Vector3d> &eb0

First point of the second edge

const Eigen::Ref<const Eigen::Vector3d> &eb1

Second point of the second edge

Returns:

A 12*3x2 matrix whose columns are the basis vectors.

Eigen::Matrix<double, 36, 2>
ipc::point_triangle_tangent_basis_jacobian(
    
const Eigen::Ref<const Eigen::Vector3d>p,
    
const Eigen::Ref<const Eigen::Vector3d>t0,
    
const Eigen::Ref<const Eigen::Vector3d>t1,
    
const Eigen::Ref<const Eigen::Vector3d>t2);

Compute the Jacobian of the tangent basis for the point-triangle pair.

Parameters:
const Eigen::Ref<const Eigen::Vector3d> &p

Point

const Eigen::Ref<const Eigen::Vector3d> &t0

Triangle’s first vertex

const Eigen::Ref<const Eigen::Vector3d> &t1

Triangle’s second vertex

const Eigen::Ref<const Eigen::Vector3d> &t2

Triangle’s third vertex

Returns:

A 12*3x2 matrix whose columns are the basis vectors.

Relative Velocity

VectorMax3d ipc::point_point_relative_velocity(
    
const Eigen::Ref<const VectorMax3d>dp0,
    
const Eigen::Ref<const VectorMax3d>dp1);

Compute the relative velocity of two points.

Parameters:
const Eigen::Ref<const VectorMax3d> &dp0

Velocity of the first point

const Eigen::Ref<const VectorMax3d> &dp1

Velocity of the second point

Returns:

The relative velocity of the two points

VectorMax3d ipc::point_edge_relative_velocity(
    
const Eigen::Ref<const VectorMax3d>dp,
    
const Eigen::Ref<const VectorMax3d>de0,
    
const Eigen::Ref<const VectorMax3d>de1const double alpha);

Compute the relative velocity of a point and an edge.

Parameters:
const Eigen::Ref<const VectorMax3d> &dp

Velocity of the point

const Eigen::Ref<const VectorMax3d> &de0

Velocity of the first endpoint of the edge

const Eigen::Ref<const VectorMax3d> &de1

Velocity of the second endpoint of the edge

const double alpha

Parametric coordinate of the closest point on the edge

Returns:

The relative velocity of the point and the edge

Eigen::Vector3d ipc::edge_edge_relative_velocity(
    
const Eigen::Ref<const Eigen::Vector3d>dea0,
    
const Eigen::Ref<const Eigen::Vector3d>dea1,
    
const Eigen::Ref<const Eigen::Vector3d>deb0,
    
const Eigen::Ref<const Eigen::Vector3d>deb1,
    
const Eigen::Ref<const Eigen::Vector2d>coords);

Compute the relative velocity of the edges.

Parameters:
const Eigen::Ref<const Eigen::Vector3d> &dea0

Velocity of the first endpoint of the first edge

const Eigen::Ref<const Eigen::Vector3d> &dea1

Velocity of the second endpoint of the first edge

const Eigen::Ref<const Eigen::Vector3d> &deb0

Velocity of the first endpoint of the second edge

const Eigen::Ref<const Eigen::Vector3d> &deb1

Velocity of the second endpoint of the second edge

const Eigen::Ref<const Eigen::Vector2d> &coords

Two parametric coordinates of the closest points on the edges

Returns:

The relative velocity of the edges

Eigen::Vector3d ipc::point_triangle_relative_velocity(
    
const Eigen::Ref<const Eigen::Vector3d>dp,
    
const Eigen::Ref<const Eigen::Vector3d>dt0,
    
const Eigen::Ref<const Eigen::Vector3d>dt1,
    
const Eigen::Ref<const Eigen::Vector3d>dt2,
    
const Eigen::Ref<const Eigen::Vector2d>coords);

Compute the relative velocity of the point to the triangle.

Parameters:
const Eigen::Ref<const Eigen::Vector3d> &dp

Velocity of the point

const Eigen::Ref<const Eigen::Vector3d> &dt0

Velocity of the first vertex of the triangle

const Eigen::Ref<const Eigen::Vector3d> &dt1

Velocity of the second vertex of the triangle

const Eigen::Ref<const Eigen::Vector3d> &dt2

Velocity of the third vertex of the triangle

const Eigen::Ref<const Eigen::Vector2d> &coords

Baricentric coordinates of the closest point on the triangle

Returns:

The relative velocity of the point to the triangle

Relative Velocity as Multiplier Matricies

MatrixMax<double, 3, 6> ipc::point_point_relative_velocity_matrix(
    
const int dim);

Compute the point-point relative velocity premultiplier matrix.

Parameters:
const int dim

Dimension (2 or 3)

Returns:

The relative velocity premultiplier matrix

MatrixMax<double, 3, 9> ipc::point_edge_relative_velocity_matrix(
    
const int dimconst double alpha);

Compute the point-edge relative velocity premultiplier matrix.

Parameters:
const int dim

Dimension (2 or 3)

const double alpha

Parametric coordinate of the closest point on the edge

Returns:

The relative velocity premultiplier matrix

MatrixMax<double, 3, 12> ipc::edge_edge_relative_velocity_matrix(
    
const int dimconst Eigen::Ref<const Eigen::Vector2d>coords);

Compute the edge-edge relative velocity matrix.

Parameters:
const int dim

Dimension (2 or 3)

const Eigen::Ref<const Eigen::Vector2d> &coords

Two parametric coordinates of the closest points on the edges

Returns:

The relative velocity matrix

MatrixMax<double, 3, 12>
ipc::point_triangle_relative_velocity_matrix(
    
const int dimconst Eigen::Ref<const Eigen::Vector2d>coords);

Compute the point-triangle relative velocity matrix.

Parameters:
const int dim

Dimension (2 or 3)

const Eigen::Ref<const Eigen::Vector2d> &coords

Baricentric coordinates of the closest point on the triangle

Returns:

The relative velocity matrix

Relative Velocity Matrix Jacobians

MatrixMax<double, 3, 6>
ipc::point_point_relative_velocity_matrix_jacobian(const int dim);

Compute the Jacobian of the relative velocity premultiplier matrix.

Parameters:
const int dim

Dimension (2 or 3)

Returns:

The Jacobian of the relative velocity premultiplier matrix

MatrixMax<double, 3, 9>
ipc::point_edge_relative_velocity_matrix_jacobian(
    
const int dimconst double alpha);

Compute the Jacobian of the relative velocity premultiplier matrix.

Parameters:
const int dim

Dimension (2 or 3)

const double alpha

Parametric coordinate of the closest point on the edge

Returns:

The Jacobian of the relative velocity premultiplier matrix

MatrixMax<double, 6, 12>
ipc::edge_edge_relative_velocity_matrix_jacobian(
    
const int dimconst Eigen::Ref<const Eigen::Vector2d>coords);

Compute the Jacobian of the edge-edge relative velocity matrix.

Parameters:
const int dim

Dimension (2 or 3)

const Eigen::Ref<const Eigen::Vector2d> &coords

Two parametric coordinates of the closest points on the edges

Returns:

The Jacobian of the relative velocity matrix

MatrixMax<double, 6, 12>
ipc::point_triangle_relative_velocity_matrix_jacobian(
    
const int dimconst Eigen::Ref<const Eigen::Vector2d>coords);

Compute the Jacobian of the point-triangle relative velocity matrix.

Parameters:
const int dim

Dimension (2 or 3)

const Eigen::Ref<const Eigen::Vector2d> &coords

Baricentric coordinates of the closest point on the triangle

Returns:

The Jacobian of the relative velocity matrix

Closet Points

double ipc::point_edge_closest_point(
    
const Eigen::Ref<const VectorMax3d>p,
    
const Eigen::Ref<const VectorMax3d>e0,
    
const Eigen::Ref<const VectorMax3d>e1);

Compute the baricentric coordinate of the closest point on the edge.

Parameters:
const Eigen::Ref<const VectorMax3d> &p

Point

const Eigen::Ref<const VectorMax3d> &e0

First edge point

const Eigen::Ref<const VectorMax3d> &e1

Second edge point

Returns:

barycentric coordinates of the closest point

Eigen::Vector2d ipc::edge_edge_closest_point(
    
const Eigen::Ref<const Eigen::Vector3d>ea0,
    
const Eigen::Ref<const Eigen::Vector3d>ea1,
    
const Eigen::Ref<const Eigen::Vector3d>eb0,
    
const Eigen::Ref<const Eigen::Vector3d>eb1);

Compute the barycentric coordinates of the closest points between two edges.

Parameters:
const Eigen::Ref<const Eigen::Vector3d> &ea0

First point of the first edge

const Eigen::Ref<const Eigen::Vector3d> &ea1

Second point of the first edge

const Eigen::Ref<const Eigen::Vector3d> &eb0

First point of the second edge

const Eigen::Ref<const Eigen::Vector3d> &eb1

Second point of the second edge

Returns:

Barycentric coordinates of the closest points

Eigen::Vector2d ipc::point_triangle_closest_point(
    
const Eigen::Ref<const Eigen::Vector3d>p,
    
const Eigen::Ref<const Eigen::Vector3d>t0,
    
const Eigen::Ref<const Eigen::Vector3d>t1,
    
const Eigen::Ref<const Eigen::Vector3d>t2);

Compute the barycentric coordinates of the closest point on the triangle.

Parameters:
const Eigen::Ref<const Eigen::Vector3d> &p

Point

const Eigen::Ref<const Eigen::Vector3d> &t0

Triangle’s first vertex

const Eigen::Ref<const Eigen::Vector3d> &t1

Triangle’s second vertex

const Eigen::Ref<const Eigen::Vector3d> &t2

Triangle’s third vertex

Returns:

Barycentric coordinates of the closest point

Closet Points Jacobians

VectorMax9d ipc::point_edge_closest_point_jacobian(
    
const Eigen::Ref<const VectorMax3d>p,
    
const Eigen::Ref<const VectorMax3d>e0,
    
const Eigen::Ref<const VectorMax3d>e1);

Compute the Jacobian of the closest point on the edge.

Parameters:
const Eigen::Ref<const VectorMax3d> &p

Point

const Eigen::Ref<const VectorMax3d> &e0

First edge point

const Eigen::Ref<const VectorMax3d> &e1

Second edge point

Returns:

Jacobian of the closest point

Eigen::Matrix<double, 2, 12> ipc::edge_edge_closest_point_jacobian(
    
const Eigen::Ref<const Eigen::Vector3d>ea0,
    
const Eigen::Ref<const Eigen::Vector3d>ea1,
    
const Eigen::Ref<const Eigen::Vector3d>eb0,
    
const Eigen::Ref<const Eigen::Vector3d>eb1);

Compute the Jacobian of the closest points between two edges.

Parameters:
const Eigen::Ref<const Eigen::Vector3d> &ea0

First point of the first edge

const Eigen::Ref<const Eigen::Vector3d> &ea1

Second point of the first edge

const Eigen::Ref<const Eigen::Vector3d> &eb0

First point of the second edge

const Eigen::Ref<const Eigen::Vector3d> &eb1

Second point of the second edge

Returns:

Jacobian of the closest points

Eigen::Matrix<double, 2, 12>
ipc::point_triangle_closest_point_jacobian(
    
const Eigen::Ref<const Eigen::Vector3d>p,
    
const Eigen::Ref<const Eigen::Vector3d>t0,
    
const Eigen::Ref<const Eigen::Vector3d>t1,
    
const Eigen::Ref<const Eigen::Vector3d>t2);

Compute the Jacobian of the closest point on the triangle.

Parameters:
const Eigen::Ref<const Eigen::Vector3d> &p

Point

const Eigen::Ref<const Eigen::Vector3d> &t0

Triangle’s first vertex

const Eigen::Ref<const Eigen::Vector3d> &t1

Triangle’s second vertex

const Eigen::Ref<const Eigen::Vector3d> &t2

Triangle’s third vertex

Returns:

Jacobian of the closest point


Last update: Feb 18, 2025