eos 1.4.0
Loading...
Searching...
No Matches
Classes | Enumerations | Functions
eos::fitting Namespace Reference

Pose and shape fitting of a 3D Morphable Model. More...

Classes

struct  ContourLandmarks
 Defines which 2D landmarks comprise the right and left face contour. More...
 
struct  FittingResult
 A struct holding the result from a fitting. More...
 
struct  Frustum
 A class representing a camera viewing frustum. At the moment only fully tested with orthographic camera. More...
 
struct  KDTreeVectorOfVectorsAdaptor
 
struct  ModelContour
 Definition of the vertex indices that define the right and left model contour. More...
 
struct  NormCost
 
struct  PerspectiveProjectionLandmarkCost
 
class  RenderingParameters
 Represents a set of estimated model parameters (rotation, translation) and camera parameters (viewing frustum). More...
 
struct  ScaledOrthoProjectionParameters
 
struct  VertexColorCost
 

Enumerations

enum class  CameraType { Orthographic , Perspective }
 

Functions

std::vector< float > fit_blendshapes_to_landmarks_linear (const std::vector< morphablemodel::Blendshape > &blendshapes, const Eigen::VectorXf &face_instance, Eigen::Matrix< float, 3, 4 > affine_camera_matrix, const std::vector< Eigen::Vector2f > &landmarks, const std::vector< int > &vertex_ids, float lambda=500.0f)
 
std::vector< float > fit_blendshapes_to_landmarks_nnls (const std::vector< morphablemodel::Blendshape > &blendshapes, const Eigen::VectorXf &face_instance, Eigen::Matrix< float, 3, 4 > affine_camera_matrix, const std::vector< Eigen::Vector2f > &landmarks, const std::vector< int > &vertex_ids)
 
template<typename T >
Eigen::Vector3< T > get_shape_at_point (const eos::morphablemodel::PcaModel &shape_model, const eos::morphablemodel::Blendshapes &blendshapes, int vertex_id, Eigen::Map< const Eigen::VectorX< T > > shape_coeffs, Eigen::Map< const Eigen::VectorX< T > > blendshape_coeffs)
 
template<typename T >
Eigen::Vector3< T > get_vertex_color_at_point (const eos::morphablemodel::PcaModel &color_model, int vertex_id, Eigen::Map< const Eigen::VectorX< T > > color_coeffs)
 
std::vector< int > occluding_boundary_vertices (const core::Mesh &mesh, const morphablemodel::EdgeTopology &edge_topology, Eigen::Matrix3f R, render::ProjectionType projection_type, bool perform_self_occlusion_check=true)
 Computes the vertices that lie on occluding boundaries, given a particular pose.
 
std::pair< std::vector< Eigen::Vector2f >, std::vector< int > > find_occluding_edge_correspondences (const core::Mesh &mesh, const morphablemodel::EdgeTopology &edge_topology, const fitting::RenderingParameters &rendering_parameters, const std::vector< Eigen::Vector2f > &image_edges, float distance_threshold=64.0f, bool perform_self_occlusion_check=true)
 For a given list of 2D edge points, find corresponding 3D vertex IDs.
 
std::pair< std::vector< std::string >, std::vector< int > > select_contour (float yaw_angle, const ContourLandmarks &contour_landmarks, const ModelContour &model_contour, float frontal_range_threshold=7.5f)
 
std::tuple< std::vector< Eigen::Vector2f >, std::vector< Eigen::Vector4f >, std::vector< int > > get_nearest_contour_correspondences (const core::LandmarkCollection< Eigen::Vector2f > &landmarks, const std::vector< std::string > &landmark_contour_identifiers, const std::vector< int > &model_contour_indices, const core::Mesh &mesh, const Eigen::Matrix4f &view_model, const Eigen::Matrix4f &ortho_projection, const Eigen::Vector4f &viewport)
 
std::tuple< std::vector< Eigen::Vector2f >, std::vector< Eigen::Vector4f >, std::vector< int > > get_contour_correspondences (const core::LandmarkCollection< Eigen::Vector2f > &landmarks, const ContourLandmarks &contour_landmarks, const ModelContour &model_contour, float yaw_angle, const core::Mesh &mesh, const Eigen::Matrix4f &view_model, const Eigen::Matrix4f &ortho_projection, const Eigen::Vector4f &viewport, float frontal_range_threshold=7.5f)
 
template<class T >
auto concat (const std::vector< T > &vec_a, const std::vector< T > &vec_b)
 Concatenates two std::vector's of the same type and returns the concatenated vector. The elements of the second vector are appended after the first one.
 
Eigen::VectorXf fit_shape (Eigen::Matrix< float, 3, 4 > affine_camera_matrix, const morphablemodel::MorphableModel &morphable_model, const std::vector< morphablemodel::Blendshape > &blendshapes, const std::vector< Eigen::Vector2f > &image_points, const std::vector< int > &vertex_indices, float lambda, cpp17::optional< int > num_coefficients_to_fit, std::vector< float > &pca_shape_coefficients, std::vector< float > &blendshape_coefficients)
 
Eigen::VectorXf fit_shape (Eigen::Matrix< float, 3, 4 > affine_camera_matrix, const morphablemodel::MorphableModel &morphable_model, const std::vector< morphablemodel::Blendshape > &blendshapes, const std::vector< Eigen::Vector2f > &image_points, const std::vector< int > &vertex_indices, float lambda=3.0f, cpp17::optional< int > num_coefficients_to_fit=cpp17::optional< int >())
 
template<class T >
auto get_corresponding_pointset (const T &landmarks, const core::LandmarkMapper &landmark_mapper, const morphablemodel::MorphableModel &morphable_model)
 Takes a LandmarkCollection of 2D landmarks and, using the landmark_mapper, finds the corresponding 3D vertex indices and returns them, along with the coordinates of the 3D points.
 
std::vector< float > fit_expressions (const morphablemodel::ExpressionModel &expression_model, const Eigen::VectorXf &face_instance, const Eigen::Matrix< float, 3, 4 > &affine_camera_matrix, const std::vector< Eigen::Vector2f > &landmarks, const std::vector< int > &vertex_ids, cpp17::optional< float > lambda_expressions=cpp17::optional< float >(), cpp17::optional< int > num_expression_coefficients_to_fit=cpp17::optional< int >())
 Fits the given expression model to landmarks.
 
std::pair< core::Mesh, fitting::RenderingParametersfit_shape_and_pose (const morphablemodel::MorphableModel &morphable_model, const core::LandmarkCollection< Eigen::Vector2f > &landmarks, const core::LandmarkMapper &landmark_mapper, int image_width, int image_height, const morphablemodel::EdgeTopology &edge_topology, const fitting::ContourLandmarks &contour_landmarks, const fitting::ModelContour &model_contour, int num_iterations, cpp17::optional< int > num_shape_coefficients_to_fit, float lambda_identity, cpp17::optional< int > num_expression_coefficients_to_fit, cpp17::optional< float > lambda_expressions, std::vector< float > &pca_shape_coefficients, std::vector< float > &expression_coefficients, std::vector< Eigen::Vector2f > &fitted_image_points)
 Fit the pose (camera), shape model, and expression blendshapes to landmarks, in an iterative way.
 
std::pair< core::Mesh, fitting::RenderingParametersfit_shape_and_pose (const morphablemodel::MorphableModel &morphable_model, const core::LandmarkCollection< Eigen::Vector2f > &landmarks, const core::LandmarkMapper &landmark_mapper, int image_width, int image_height, const morphablemodel::EdgeTopology &edge_topology, const fitting::ContourLandmarks &contour_landmarks, const fitting::ModelContour &model_contour, int num_iterations=5, cpp17::optional< int > num_shape_coefficients_to_fit=cpp17::nullopt, float lambda_identity=50.0f, cpp17::optional< int > num_expression_coefficients_to_fit=cpp17::nullopt, cpp17::optional< float > lambda_expressions=cpp17::nullopt)
 Fit the pose (camera), shape model, and expression blendshapes to landmarks, in an iterative way.
 
std::pair< core::Mesh, fitting::RenderingParametersfit_shape_and_pose (const morphablemodel::MorphableModel &morphable_model, const std::vector< Eigen::Vector2f > &image_points, const std::vector< int > &vertex_indices, int image_width, int image_height, int num_iterations, cpp17::optional< int > num_shape_coefficients_to_fit, float lambda_identity, cpp17::optional< int > num_expression_coefficients_to_fit, cpp17::optional< float > lambda_expressions, std::vector< float > &pca_shape_coefficients, std::vector< float > &expression_coefficients, std::vector< Eigen::Vector2f > &fitted_image_points)
 Fit the pose (camera), shape model, and expression blendshapes to landmarks, in an iterative way.
 
std::pair< core::Mesh, fitting::RenderingParametersfit_shape_and_pose (const morphablemodel::MorphableModel &morphable_model, const std::vector< Eigen::Vector2f > &image_points, const std::vector< int > &vertex_indices, int image_width, int image_height, int num_iterations, cpp17::optional< int > num_shape_coefficients_to_fit, float lambda_identity, cpp17::optional< int > num_expression_coefficients_to_fit, cpp17::optional< float > lambda_expressions)
 Fit the pose (camera), shape model, and expression blendshapes to landmarks, in an iterative way, given fixed landmarks-to-vertex correspondences.
 
std::vector< float > fit_shape_to_landmarks_linear (const morphablemodel::PcaModel &shape_model, Eigen::Matrix< float, 3, 4 > affine_camera_matrix, const std::vector< Eigen::Vector2f > &landmarks, const std::vector< int > &vertex_ids, Eigen::VectorXf base_face=Eigen::VectorXf(), float lambda=3.0f, cpp17::optional< int > num_coefficients_to_fit=cpp17::optional< int >(), cpp17::optional< float > detector_standard_deviation=cpp17::optional< float >(), cpp17::optional< float > model_standard_deviation=cpp17::optional< float >())
 
std::vector< float > fit_shape_to_landmarks_linear_multi (const morphablemodel::PcaModel &shape_model, const std::vector< Eigen::Matrix< float, 3, 4 > > &affine_camera_matrices, const std::vector< std::vector< Eigen::Vector2f > > &landmarks, const std::vector< std::vector< int > > &vertex_ids, std::vector< Eigen::VectorXf > base_faces=std::vector< Eigen::VectorXf >(), float lambda=3.0f, cpp17::optional< int > num_coefficients_to_fit=cpp17::optional< int >(), cpp17::optional< float > detector_standard_deviation=cpp17::optional< float >(), cpp17::optional< float > model_standard_deviation=cpp17::optional< float >())
 
std::pair< std::vector< core::Mesh >, std::vector< fitting::RenderingParameters > > fit_shape_and_pose (const morphablemodel::MorphableModel &morphable_model, const std::vector< morphablemodel::Blendshape > &blendshapes, const std::vector< core::LandmarkCollection< Eigen::Vector2f > > &landmarks, const core::LandmarkMapper &landmark_mapper, std::vector< int > image_width, std::vector< int > image_height, const morphablemodel::EdgeTopology &edge_topology, const fitting::ContourLandmarks &contour_landmarks, const fitting::ModelContour &model_contour, int num_iterations, cpp17::optional< int > num_shape_coefficients_to_fit, float lambda, cpp17::optional< fitting::RenderingParameters > initial_rendering_params, std::vector< float > &pca_shape_coefficients, std::vector< std::vector< float > > &blendshape_coefficients, std::vector< std::vector< Eigen::Vector2f > > &fitted_image_points)
 Fit the pose (camera), shape model, and expression blendshapes to landmarks, in an iterative way. This function takes a set of images and landmarks and estimates per-frame pose and expressions, as well as identity shape jointly from all images.
 
std::pair< std::vector< core::Mesh >, std::vector< fitting::RenderingParameters > > fit_shape_and_pose (const morphablemodel::MorphableModel &morphable_model, const std::vector< morphablemodel::Blendshape > &blendshapes, const std::vector< core::LandmarkCollection< Eigen::Vector2f > > &landmarks, const core::LandmarkMapper &landmark_mapper, std::vector< int > image_width, std::vector< int > image_height, const morphablemodel::EdgeTopology &edge_topology, const fitting::ContourLandmarks &contour_landmarks, const fitting::ModelContour &model_contour, int num_iterations=5, cpp17::optional< int > num_shape_coefficients_to_fit=cpp17::nullopt, float lambda=30.0f)
 Fit the pose (camera), shape model, and expression blendshapes to landmarks, in an iterative way. This function takes a set of images and landmarks and estimates per-frame pose and expressions, as well as identity shape jointly from all images.
 
RenderingParameters estimate_orthographic_camera (std::vector< Eigen::Vector2f > image_points, std::vector< Eigen::Vector4f > model_points, int width, int height)
 This algorithm estimates the rotation angles and translation of the model, as well as the viewing frustum of the camera, given a set of corresponding 2D-3D points.
 
ScaledOrthoProjectionParameters estimate_orthographic_projection_linear (std::vector< Eigen::Vector2f > image_points, std::vector< Eigen::Vector4f > model_points, bool is_viewport_upsidedown, cpp17::optional< int > viewport_height=cpp17::nullopt)
 
void save_rendering_parameters (RenderingParameters rendering_parameters, std::string filename)
 
Eigen::Vector4f get_opencv_viewport (int width, int height)
 Returns a glm/OpenGL compatible viewport vector that flips y and has the origin on the top-left, like in OpenCV.
 
Eigen::Matrix< float, 3, 4 > get_3x4_affine_camera_matrix (RenderingParameters params, int width, int height)
 Creates a 3x4 affine camera matrix from given fitting parameters. The matrix transforms points directly from model-space to screen-space.
 
template<typename T >
Eigen::Matrix< T, 3, 1 > tait_bryan_angles (Eigen::Matrix< T, 3, 3 > rotation_matrix, Eigen::Index axis_0, Eigen::Index axis_1, Eigen::Index axis_2)
 Computes Tait-Bryan angles (in radians) from the given rotation matrix and axes order.
 

Detailed Description

Pose and shape fitting of a 3D Morphable Model.

Function Documentation

◆ concat()

template<class T >
auto eos::fitting::concat ( const std::vector< T > &  vec_a,
const std::vector< T > &  vec_b 
)
inline

Concatenates two std::vector's of the same type and returns the concatenated vector. The elements of the second vector are appended after the first one.

Concatenating two vectors is commonly needed when using the contour fitting, to add contour points to an existing vector of 2D/3D point correspondences.

Parameters
[in]vec_aFirst vector.
[in]vec_bSecond vector.
Returns
The concatenated vector.

◆ estimate_orthographic_camera()

RenderingParameters eos::fitting::estimate_orthographic_camera ( std::vector< Eigen::Vector2f >  image_points,
std::vector< Eigen::Vector4f >  model_points,
int  width,
int  height 
)
inline

This algorithm estimates the rotation angles and translation of the model, as well as the viewing frustum of the camera, given a set of corresponding 2D-3D points.

It assumes an orthographic camera and estimates 6 parameters, [r_x, r_y, r_z, t_x, t_y, frustum_scale], where the first five describe how to transform the model, and the last one describes the cameras viewing frustum (see CameraParameters). This 2D-3D correspondence problem is solved using Eigen's LevenbergMarquardt algorithm.

The method is slightly inspired by "Computer Vision: Models Learning and Inference", Simon J.D. Prince, 2012, but different in a lot of respects.

Eigen's LM implementation requires at least 6 data points, so we require >= 6 corresponding points.

Notes/improvements: The algorithm works reliable as it is, however, it could be improved with the following:

  • A better initial guess (see e.g. Prince)
  • We could add a parameter to pass an initial guess
  • Using the analytic derivatives instead of Eigen::NumericalDiff - they're easy to calculate.
Parameters
[in]image_pointsA list of 2D image points.
[in]model_pointsCorresponding points of a 3D model.
[in]widthWidth of the image (or viewport).
[in]heightHeight of the image (or viewport).
Returns
The estimated model and camera parameters.

◆ estimate_orthographic_projection_linear()

ScaledOrthoProjectionParameters eos::fitting::estimate_orthographic_projection_linear ( std::vector< Eigen::Vector2f >  image_points,
std::vector< Eigen::Vector4f >  model_points,
bool  is_viewport_upsidedown,
cpp17::optional< int >  viewport_height = cpp17::nullopt 
)
inline

Estimates the parameters of a scaled orthographic projection.

Given a set of 2D-3D correspondences, this algorithm estimates rotation, translation (in x and y) and a scaling parameters of the scaled orthographic projection model using a closed-form solution. It does so by first computing an affine camera matrix using algorithm [1], and then finds the closest orthonormal matrix to the estimated affine transform using SVD. This algorithm follows the original implementation 2 of William Smith, University of York.

Requires >= 4 corresponding points.

[1]: Gold Standard Algorithm for estimating an affine camera matrix from world to image correspondences, Algorithm 7.2 in Multiple View Geometry, Hartley & Zisserman, 2nd Edition, 2003.

Parameters
[in]image_pointsA list of 2D image points.
[in]model_pointsCorresponding points of a 3D model.
[in]is_viewport_upsidedownFlag to set whether the viewport of the image points is upside-down (e.g. as in OpenCV).
[in]viewport_heightHeight of the viewport of the image points (needs to be given if is_viewport_upsidedown == true).
Returns
Rotation, translation and scaling of the estimated scaled orthographic projection.

◆ find_occluding_edge_correspondences()

std::pair< std::vector< Eigen::Vector2f >, std::vector< int > > eos::fitting::find_occluding_edge_correspondences ( const core::Mesh mesh,
const morphablemodel::EdgeTopology edge_topology,
const fitting::RenderingParameters rendering_parameters,
const std::vector< Eigen::Vector2f > &  image_edges,
float  distance_threshold = 64.0f,
bool  perform_self_occlusion_check = true 
)
inline

For a given list of 2D edge points, find corresponding 3D vertex IDs.

This algorithm first computes the 3D mesh's occluding boundary vertices under the given pose. Then, for each 2D image edge point given, it searches for the closest 3D edge vertex (projected to 2D). Correspondences lying further away than distance_threshold (times a scale-factor) are discarded. It returns a list of the remaining image edge points and their corresponding 3D vertex ID.

The given rendering_parameters camery_type must be CameraType::Orthographic.

The units of distance_threshold are somewhat complicated. The function uses squared distances, and the distance_threshold is further multiplied with a face-size and image resolution dependent scale factor. It's reasonable to use correspondences that are 10 to 15 pixels away on a 1280x720 image with s=0.93. This would be a distance_threshold of around 200. 64 might be a conservative default.

Parameters
[in]meshThe 3D mesh.
[in]edge_topologyThe mesh's edge topology (used for fast computation).
[in]rendering_parametersRendering (pose) parameters of the mesh.
[in]image_edgesA list of points that are edges.
[in]distance_thresholdAll correspondences below this threshold.
[in]perform_self_occlusion_checkCheck the computed boundary vertices for self-occlusion and remove these.
Returns
A pair consisting of the used image edge points and their associated 3D vertex index.

◆ fit_blendshapes_to_landmarks_linear()

std::vector< float > eos::fitting::fit_blendshapes_to_landmarks_linear ( const std::vector< morphablemodel::Blendshape > &  blendshapes,
const Eigen::VectorXf &  face_instance,
Eigen::Matrix< float, 3, 4 >  affine_camera_matrix,
const std::vector< Eigen::Vector2f > &  landmarks,
const std::vector< int > &  vertex_ids,
float  lambda = 500.0f 
)
inline

Fits blendshape coefficients to given 2D landmarks, given a current face shape instance. It's a linear, closed-form solution fitting algorithm, with regularisation (constraining the L2-norm of the coefficients). However, there is no constraint on the coefficients, so negative coefficients are allowed, which, with linear blendshapes (offsets), will most likely not be desired. Thus, prefer the function fit_blendshapes_to_landmarks_nnls(std::vector<eos::morphablemodel::Blendshape>, cv::Mat, cv::Mat, std::vector<cv::Vec2f>, std::vector<int>).

This algorithm is very similar to the shape fitting in fit_shape_to_landmarks_linear. Instead of the PCA basis, the blendshapes are used, and instead of the mean, a current face instance is used to do the fitting from.

Parameters
[in]blendshapesA vector with blendshapes to estimate the coefficients for.
[in]face_instanceA shape instance from which the blendshape coefficients should be estimated (i.e. the current mesh without expressions, e.g. estimated from a previous PCA-model fitting). A 3m x 1 matrix.
[in]affine_camera_matrixA 3x4 affine camera matrix from model to screen-space.
[in]landmarks2D landmarks from an image to fit the blendshapes to.
[in]vertex_idsThe vertex ids in the model that correspond to the 2D points.
[in]lambdaA regularisation parameter, constraining the L2-norm of the coefficients.
Returns
The estimated blendshape-coefficients.

◆ fit_blendshapes_to_landmarks_nnls()

std::vector< float > eos::fitting::fit_blendshapes_to_landmarks_nnls ( const std::vector< morphablemodel::Blendshape > &  blendshapes,
const Eigen::VectorXf &  face_instance,
Eigen::Matrix< float, 3, 4 >  affine_camera_matrix,
const std::vector< Eigen::Vector2f > &  landmarks,
const std::vector< int > &  vertex_ids 
)
inline

Fits blendshape coefficients to given 2D landmarks, given a current face shape instance. Uses non-negative least-squares (NNLS) to solve for the coefficients. The NNLS algorithm used doesn't support any regularisation.

This algorithm is very similar to the shape fitting in fit_shape_to_landmarks_linear. Instead of the PCA basis, the blendshapes are used, and instead of the mean, a current face instance is used to do the fitting from.

Parameters
[in]blendshapesA vector with blendshapes to estimate the coefficients for.
[in]face_instanceA shape instance from which the blendshape coefficients should be estimated (i.e. the current mesh without expressions, e.g. estimated from a previous PCA-model fitting). A 3m x 1 matrix.
[in]affine_camera_matrixA 3x4 affine camera matrix from model to screen-space.
[in]landmarks2D landmarks from an image to fit the blendshapes to.
[in]vertex_idsThe vertex ids in the model that correspond to the 2D points.
Returns
The estimated blendshape-coefficients.

◆ fit_expressions()

std::vector< float > eos::fitting::fit_expressions ( const morphablemodel::ExpressionModel expression_model,
const Eigen::VectorXf &  face_instance,
const Eigen::Matrix< float, 3, 4 > &  affine_camera_matrix,
const std::vector< Eigen::Vector2f > &  landmarks,
const std::vector< int > &  vertex_ids,
cpp17::optional< float >  lambda_expressions = cpp17::optional<float>(),
cpp17::optional< int >  num_expression_coefficients_to_fit = cpp17::optional<int>() 
)
inline

Fits the given expression model to landmarks.

The function uses fit_blendshapes_to_landmarks_nnls(...) if the given expression model consists of Blendshapes, and fit_shape_to_landmarks_linear(...) if it is a PCA model.

Note that if the expression model is a PCA model, and we're doing fit_shape_to_landmarks_linear(...), we should probably pass the regularisation and num_coeffs_to_fit along as well.

Parameters
[in]expression_modelThe expression model (blendshapes or a PCA model).
[in]face_instanceA shape instance from which the expression coefficients should be estimated (i.e. the current mesh without expressions, e.g. estimated from a previous PCA-model fitting).
[in]affine_camera_matrixSecond vector.
[in]landmarks2D landmarks from an image to fit the blendshapes to.
[in]vertex_idsThe vertex ids in the model that correspond to the 2D points.
[in]lambda_expressionsThe regularisation parameter (weight of the prior towards the mean). Only used for expression-PCA fitting.
[in]num_expression_coefficients_to_fitHow many expression coefficients to fit (all others will stay 0). Should be bigger than zero, or std::nullopt to fit all coefficients. Only used for expression-PCA fitting.
Returns
A vector of fitted expression coefficients.

◆ fit_shape() [1/2]

Eigen::VectorXf eos::fitting::fit_shape ( Eigen::Matrix< float, 3, 4 >  affine_camera_matrix,
const morphablemodel::MorphableModel morphable_model,
const std::vector< morphablemodel::Blendshape > &  blendshapes,
const std::vector< Eigen::Vector2f > &  image_points,
const std::vector< int > &  vertex_indices,
float  lambda,
cpp17::optional< int >  num_coefficients_to_fit,
std::vector< float > &  pca_shape_coefficients,
std::vector< float > &  blendshape_coefficients 
)
inline

Convenience function that fits the shape model and expression blendshapes to landmarks. Makes the fitted PCA shape and blendshape coefficients accessible via the out parameters pca_shape_coefficients and blendshape_coefficients. It iterates PCA-shape and blendshape fitting until convergence (usually it converges within 5 to 10 iterations).

See fit_shape_model(cv::Mat, eos::morphablemodel::MorphableModel, std::vector<eos::morphablemodel::Blendshape>, std::vector<cv::Vec2f>, std::vector<int>, float lambda) for a simpler overload that just returns the shape instance.

Parameters
[in]affine_camera_matrixThe estimated pose as a 3x4 affine camera matrix that is used to fit the shape.
[in]morphable_modelThe 3D Morphable Model used for the shape fitting.
[in]blendshapesA vector of blendshapes that are being fit to the landmarks in addition to the PCA model.
[in]image_points2D landmarks from an image to fit the model to.
[in]vertex_indicesThe vertex indices in the model that correspond to the 2D points.
[in]lambdaRegularisation parameter of the PCA shape fitting.
[in]num_coefficients_to_fitHow many shape-coefficients to fit (all others will stay 0). Should be bigger than zero, or std::nullopt to fit all coefficients.
[out]pca_shape_coefficientsOutput parameter that will contain the resulting pca shape coefficients.
[out]blendshape_coefficientsOutput parameter that will contain the resulting blendshape coefficients.
Returns
The fitted model shape instance.

◆ fit_shape() [2/2]

Eigen::VectorXf eos::fitting::fit_shape ( Eigen::Matrix< float, 3, 4 >  affine_camera_matrix,
const morphablemodel::MorphableModel morphable_model,
const std::vector< morphablemodel::Blendshape > &  blendshapes,
const std::vector< Eigen::Vector2f > &  image_points,
const std::vector< int > &  vertex_indices,
float  lambda = 3.0f,
cpp17::optional< int >  num_coefficients_to_fit = cpp17::optional<int>() 
)
inline

Convenience function that fits the shape model and expression blendshapes to landmarks. It iterates PCA-shape and blendshape fitting until convergence (usually it converges within 5 to 10 iterations).

Parameters
[in]affine_camera_matrixThe estimated pose as a 3x4 affine camera matrix that is used to fit the shape.
[in]morphable_modelThe 3D Morphable Model used for the shape fitting.
[in]blendshapesA vector of blendshapes that are being fit to the landmarks in addition to the PCA model.
[in]image_points2D landmarks from an image to fit the model to.
[in]vertex_indicesThe vertex indices in the model that correspond to the 2D points.
[in]lambdaRegularisation parameter of the PCA shape fitting.
[in]num_coefficients_to_fitHow many shape-coefficients to fit (all others will stay 0). Should be bigger than zero, or std::nullopt to fit all coefficients.
Returns
The fitted model shape instance.

◆ fit_shape_and_pose() [1/6]

std::pair< core::Mesh, fitting::RenderingParameters > eos::fitting::fit_shape_and_pose ( const morphablemodel::MorphableModel morphable_model,
const core::LandmarkCollection< Eigen::Vector2f > &  landmarks,
const core::LandmarkMapper landmark_mapper,
int  image_width,
int  image_height,
const morphablemodel::EdgeTopology edge_topology,
const fitting::ContourLandmarks contour_landmarks,
const fitting::ModelContour model_contour,
int  num_iterations,
cpp17::optional< int >  num_shape_coefficients_to_fit,
float  lambda_identity,
cpp17::optional< int >  num_expression_coefficients_to_fit,
cpp17::optional< float >  lambda_expressions,
std::vector< float > &  pca_shape_coefficients,
std::vector< float > &  expression_coefficients,
std::vector< Eigen::Vector2f > &  fitted_image_points 
)
inline

Fit the pose (camera), shape model, and expression blendshapes to landmarks, in an iterative way.

Convenience function that fits pose (camera), the shape model, and expression blendshapes to landmarks, in an iterative (alternating) way. It fits both sides of the face contour as well.

If pca_shape_coefficients and/or blendshape_coefficients are given, they are used as starting values in the fitting. When the function returns, they contain the coefficients from the last iteration.

Use render::Mesh fit_shape_and_pose(const morphablemodel::MorphableModel&, const std::vector<morphablemodel::Blendshape>&, const core::LandmarkCollection<cv::Vec2f>&, const core::LandmarkMapper&, int, int, const morphablemodel::EdgeTopology&, const fitting::ContourLandmarks&, const fitting::ModelContour&, int, cpp17::optional<int>, float). for a simpler overload with reasonable defaults and no optional output.

num_iterations: Results are good for even a single iteration. For single-image fitting and for full convergence of all parameters, it can take up to 300 iterations. In tracking, particularly if initialising with the previous frame, it works well with as low as 1 to 5 iterations. edge_topology is used for the occluding-edge face contour fitting. contour_landmarks and model_contour are used to fit the front-facing contour.

Note: If the given morphable_model contains a PCA expression model, alternating the shape identity and expression fitting is theoretically not needed - the basis matrices could be stacked, and then both coefficients could be solved for in one go. The two bases are most likely not orthogonal though. In any case, alternating hopefully doesn't do any harm.

Todo: Add a convergence criterion.

Parameters
[in]morphable_modelThe 3D Morphable Model used for the shape fitting.
[in]landmarks2D landmarks from an image to fit the model to.
[in]landmark_mapperMapping info from the 2D landmark points to 3D vertex indices.
[in]image_widthWidth of the input image (needed for the camera model).
[in]image_heightHeight of the input image (needed for the camera model).
[in]edge_topologyPrecomputed edge topology of the 3D model, needed for fast edge-lookup.
[in]contour_landmarks2D image contour ids of left or right side (for example for ibug landmarks).
[in]model_contourThe model contour indices that should be considered to find the closest corresponding 3D vertex.
[in]num_iterationsNumber of iterations that the different fitting parts will be alternated for.
[in]num_shape_coefficients_to_fitHow many shape-coefficients to fit (all others will stay 0). Should be bigger than zero, or std::nullopt to fit all coefficients.
[in]lambda_identityRegularisation parameter of the PCA shape fitting.
[in]num_expression_coefficients_to_fitHow many shape-coefficients to fit (all others will stay 0). Should be bigger than zero, or std::nullopt to fit all coefficients. Only used for expression-PCA fitting.
[in]lambda_expressionsRegularisation parameter of the expression fitting. Only used for expression-PCA fitting.
[in,out]pca_shape_coefficientsIf given, will be used as initial PCA shape coefficients to start the fitting. Will contain the final estimated coefficients.
[in,out]expression_coefficientsIf given, will be used as initial expression blendshape coefficients to start the fitting. Will contain the final estimated coefficients.
[out]fitted_image_pointsDebug parameter: Returns all the 2D points that have been used for the fitting.
Returns
The fitted model shape instance and the final pose.
Exceptions
std::runtime_errorif an error occurs, for example more shape coefficients are given than the model contains

◆ fit_shape_and_pose() [2/6]

std::pair< core::Mesh, fitting::RenderingParameters > eos::fitting::fit_shape_and_pose ( const morphablemodel::MorphableModel morphable_model,
const core::LandmarkCollection< Eigen::Vector2f > &  landmarks,
const core::LandmarkMapper landmark_mapper,
int  image_width,
int  image_height,
const morphablemodel::EdgeTopology edge_topology,
const fitting::ContourLandmarks contour_landmarks,
const fitting::ModelContour model_contour,
int  num_iterations = 5,
cpp17::optional< int >  num_shape_coefficients_to_fit = cpp17::nullopt,
float  lambda_identity = 50.0f,
cpp17::optional< int >  num_expression_coefficients_to_fit = cpp17::nullopt,
cpp17::optional< float >  lambda_expressions = cpp17::nullopt 
)
inline

Fit the pose (camera), shape model, and expression blendshapes to landmarks, in an iterative way.

Convenience function that fits pose (camera), the shape model, and expression blendshapes to landmarks, in an iterative (alternating) way. It fits both sides of the face contour as well.

If you want to access the values of shape or blendshape coefficients, or want to set starting values for them, use the following overload to this function: std::pair<render::Mesh, fitting::RenderingParameters> fit_shape_and_pose(const morphablemodel::MorphableModel&, const std::vector<morphablemodel::Blendshape>&, const core::LandmarkCollection<cv::Vec2f>&, const core::LandmarkMapper&, int, int, const morphablemodel::EdgeTopology&, const fitting::ContourLandmarks&, const fitting::ModelContour&, int, cpp17::optional<int>, float, cpp17::optional<fitting::RenderingParameters>, std::vector<float>&, std::vector<float>&, std::vector<cv::Vec2f>&)

Todo: Add a convergence criterion.

num_iterations: Results are good for even a single iteration. For single-image fitting and for full convergence of all parameters, it can take up to 300 iterations. In tracking, particularly if initialising with the previous frame, it works well with as low as 1 to 5 iterations. edge_topology is used for the occluding-edge face contour fitting. contour_landmarks and model_contour are used to fit the front-facing contour.

Parameters
[in]morphable_modelThe 3D Morphable Model used for the shape fitting.
[in]landmarks2D landmarks from an image to fit the model to.
[in]landmark_mapperMapping info from the 2D landmark points to 3D vertex indices.
[in]image_widthWidth of the input image (needed for the camera model).
[in]image_heightHeight of the input image (needed for the camera model).
[in]edge_topologyPrecomputed edge topology of the 3D model, needed for fast edge-lookup.
[in]contour_landmarks2D image contour ids of left or right side (for example for ibug landmarks).
[in]model_contourThe model contour indices that should be considered to find the closest corresponding 3D vertex.
[in]num_iterationsNumber of iterations that the different fitting parts will be alternated for.
[in]num_shape_coefficients_to_fitHow many shape-coefficients to fit (all others will stay 0). Should be bigger than zero, or std::nullopt to fit all coefficients.
[in]lambda_identityRegularisation parameter of the PCA shape fitting.
[in]num_expression_coefficients_to_fitHow many shape-coefficients to fit (all others will stay 0). Should be bigger than zero, or std::nullopt to fit all coefficients. Only used for expression-PCA fitting.
[in]lambda_expressionsRegularisation parameter of the expression fitting. Only used for expression-PCA fitting.
Returns
The fitted model shape instance and the final pose.

◆ fit_shape_and_pose() [3/6]

std::pair< core::Mesh, fitting::RenderingParameters > eos::fitting::fit_shape_and_pose ( const morphablemodel::MorphableModel morphable_model,
const std::vector< Eigen::Vector2f > &  image_points,
const std::vector< int > &  vertex_indices,
int  image_width,
int  image_height,
int  num_iterations,
cpp17::optional< int >  num_shape_coefficients_to_fit,
float  lambda_identity,
cpp17::optional< int >  num_expression_coefficients_to_fit,
cpp17::optional< float >  lambda_expressions 
)
inline

Fit the pose (camera), shape model, and expression blendshapes to landmarks, in an iterative way, given fixed landmarks-to-vertex correspondences.

This is a convenience overload without the last 3 in/out parameters of above function.

Parameters
[in]morphable_modelThe 3D Morphable Model used for the shape fitting.
[in]image_points2D image points to fit the model to.
[in]vertex_indicesThe 3D vertex indices corresponding to the given image_points.
[in]image_widthWidth of the input image (needed for the camera model).
[in]image_heightHeight of the input image (needed for the camera model).
[in]num_iterationsNumber of iterations that the different fitting parts will be alternated for.
[in]num_shape_coefficients_to_fitHow many shape-coefficients to fit (all others will stay 0). Should be bigger than zero, or std::nullopt to fit all coefficients.
[in]lambda_identityRegularisation parameter of the PCA shape fitting.
[in]num_expression_coefficients_to_fitHow many shape-coefficients to fit (all others will stay 0). Should be bigger than zero, or std::nullopt to fit all coefficients. Only used for expression-PCA fitting.
[in]lambda_expressionsRegularisation parameter of the expression fitting. Only used for expression-PCA fitting.
Returns
The fitted model shape instance and the final pose.

◆ fit_shape_and_pose() [4/6]

std::pair< core::Mesh, fitting::RenderingParameters > eos::fitting::fit_shape_and_pose ( const morphablemodel::MorphableModel morphable_model,
const std::vector< Eigen::Vector2f > &  image_points,
const std::vector< int > &  vertex_indices,
int  image_width,
int  image_height,
int  num_iterations,
cpp17::optional< int >  num_shape_coefficients_to_fit,
float  lambda_identity,
cpp17::optional< int >  num_expression_coefficients_to_fit,
cpp17::optional< float >  lambda_expressions,
std::vector< float > &  pca_shape_coefficients,
std::vector< float > &  expression_coefficients,
std::vector< Eigen::Vector2f > &  fitted_image_points 
)
inline

Fit the pose (camera), shape model, and expression blendshapes to landmarks, in an iterative way.

Convenience function that fits pose (camera), the shape model, and expression blendshapes to landmarks, in an iterative (alternating) way. It uses the given, fixed landmarks-to-vertex correspondences, and does not use any dynamic contour fitting.

If pca_shape_coefficients and/or blendshape_coefficients are given, they are used as starting values in the fitting. When the function returns, they contain the coefficients from the last iteration.

num_iterations: Results are good for even a single iteration. For single-image fitting and for full convergence of all parameters, it can take up to 300 iterations. In tracking, particularly if initialising with the previous frame, it works well with as low as 1 to 5 iterations. edge_topology is used for the occluding-edge face contour fitting. contour_landmarks and model_contour are used to fit the front-facing contour.

Note: If the given morphable_model contains a PCA expression model, alternating the shape identity and expression fitting is theoretically not needed - the basis matrices could be stacked, and then both coefficients could be solved for in one go. The two bases are most likely not orthogonal though. In any case, alternating hopefully doesn't do any harm.

Todo: Add a convergence criterion.

Parameters
[in]morphable_modelThe 3D Morphable Model used for the shape fitting.
[in]image_points2D image points to fit the model to.
[in]vertex_indicesThe 3D vertex indices corresponding to the given image_points.
[in]image_widthWidth of the input image (needed for the camera model).
[in]image_heightHeight of the input image (needed for the camera model).
[in]num_iterationsNumber of iterations that the different fitting parts will be alternated for.
[in]num_shape_coefficients_to_fitHow many shape-coefficients to fit (all others will stay 0). Should be bigger than zero, or std::nullopt to fit all coefficients.
[in]lambda_identityRegularisation parameter of the PCA shape fitting.
[in]num_expression_coefficients_to_fitHow many shape-coefficients to fit (all others will stay 0). Should be bigger than zero, or std::nullopt to fit all coefficients. Only used for expression-PCA fitting.
[in]lambda_expressionsRegularisation parameter of the expression fitting. Only used for expression-PCA fitting.
[in,out]pca_shape_coefficientsIf given, will be used as initial PCA shape coefficients to start the fitting. Will contain the final estimated coefficients.
[in,out]expression_coefficientsIf given, will be used as initial expression blendshape coefficients to start the fitting. Will contain the final estimated coefficients.
[out]fitted_image_pointsDebug parameter: Returns all the 2D points that have been used for the fitting.
Returns
The fitted model shape instance and the final pose.

◆ fit_shape_and_pose() [5/6]

std::pair< std::vector< core::Mesh >, std::vector< fitting::RenderingParameters > > eos::fitting::fit_shape_and_pose ( const morphablemodel::MorphableModel morphable_model,
const std::vector< morphablemodel::Blendshape > &  blendshapes,
const std::vector< core::LandmarkCollection< Eigen::Vector2f > > &  landmarks,
const core::LandmarkMapper landmark_mapper,
std::vector< int >  image_width,
std::vector< int >  image_height,
const morphablemodel::EdgeTopology edge_topology,
const fitting::ContourLandmarks contour_landmarks,
const fitting::ModelContour model_contour,
int  num_iterations,
cpp17::optional< int >  num_shape_coefficients_to_fit,
float  lambda,
cpp17::optional< fitting::RenderingParameters initial_rendering_params,
std::vector< float > &  pca_shape_coefficients,
std::vector< std::vector< float > > &  blendshape_coefficients,
std::vector< std::vector< Eigen::Vector2f > > &  fitted_image_points 
)
inline

Fit the pose (camera), shape model, and expression blendshapes to landmarks, in an iterative way. This function takes a set of images and landmarks and estimates per-frame pose and expressions, as well as identity shape jointly from all images.

Convenience function that fits pose (camera), the shape model, and expression blendshapes to landmarks, in an iterative (alternating) way. It fits both sides of the face contour as well.

If pca_shape_coefficients and/or blendshape_coefficients are given, they are used as starting values in the fitting. When the function returns, they contain the coefficients from the last iteration.

Use render::Mesh fit_shape_and_pose(const morphablemodel::MorphableModel&, const std::vector<morphablemodel::Blendshape>&, const core::LandmarkCollection<Eigen::Vector2f>&, const core::LandmarkMapper&, int, int, const morphablemodel::EdgeTopology&, const fitting::ContourLandmarks&, const fitting::ModelContour&, int, cpp17::optional<int>, float). for a simpler overload with reasonable defaults and no optional output.

num_iterations: Results are good for even a few iterations. For full convergence of all parameters, it can take up to 300 iterations. In tracking, particularly if initialising with the previous frame, it works well with as low as 1 to 5 iterations. edge_topology is used for the occluding-edge face contour fitting. contour_landmarks and model_contour are used to fit the front-facing contour.

Todo: Add a convergence criterion.

Parameters
[in]morphable_modelThe 3D Morphable Model used for the shape fitting.
[in]blendshapesA vector of blendshapes that are being fit to the landmarks in addition to the PCA model.
[in]landmarks2D landmarks from an image to fit the model to.
[in]landmark_mapperMapping info from the 2D landmark points to 3D vertex indices.
[in]image_widthWidth of the input image (needed for the camera model).
[in]image_heightHeight of the input image (needed for the camera model).
[in]edge_topologyPrecomputed edge topology of the 3D model, needed for fast edge-lookup.
[in]contour_landmarks2D image contour ids of left or right side (for example for ibug landmarks).
[in]model_contourThe model contour indices that should be considered to find the closest corresponding 3D vertex.
[in]num_iterationsNumber of iterations that the different fitting parts will be alternated for.
[in]num_shape_coefficients_to_fitHow many shape-coefficients to fit (all others will stay 0). Should be bigger than zero, or std::nullopt to fit all coefficients.
[in]lambdaRegularisation parameter of the PCA shape fitting.
[in]initial_rendering_paramsCurrently ignored (not used).
[in,out]pca_shape_coefficientsIf given, will be used as initial PCA shape coefficients to start the fitting. Will contain the final estimated coefficients.
[in,out]blendshape_coefficientsIf given, will be used as initial expression blendshape coefficients to start the fitting. Will contain the final estimated coefficients.
[out]fitted_image_pointsDebug parameter: Returns all the 2D points that have been used for the fitting.
Returns
The fitted model shape instance and the final pose.

◆ fit_shape_and_pose() [6/6]

std::pair< std::vector< core::Mesh >, std::vector< fitting::RenderingParameters > > eos::fitting::fit_shape_and_pose ( const morphablemodel::MorphableModel morphable_model,
const std::vector< morphablemodel::Blendshape > &  blendshapes,
const std::vector< core::LandmarkCollection< Eigen::Vector2f > > &  landmarks,
const core::LandmarkMapper landmark_mapper,
std::vector< int >  image_width,
std::vector< int >  image_height,
const morphablemodel::EdgeTopology edge_topology,
const fitting::ContourLandmarks contour_landmarks,
const fitting::ModelContour model_contour,
int  num_iterations = 5,
cpp17::optional< int >  num_shape_coefficients_to_fit = cpp17::nullopt,
float  lambda = 30.0f 
)
inline

Fit the pose (camera), shape model, and expression blendshapes to landmarks, in an iterative way. This function takes a set of images and landmarks and estimates per-frame pose and expressions, as well as identity shape jointly from all images.

Convenience function that fits pose (camera), the shape model, and expression blendshapes to landmarks, in an iterative (alternating) way. It fits both sides of the face contour as well.

If pca_shape_coefficients and/or blendshape_coefficients are given, they are used as starting values in the fitting. When the function returns, they contain the coefficients from the last iteration.

If you want to access the values of shape or blendshape coefficients, or want to set starting values for them, use the following overload to this function: std::pair<render::Mesh, fitting::RenderingParameters> fit_shape_and_pose(const morphablemodel::MorphableModel&, const std::vector<morphablemodel::Blendshape>&, const core::LandmarkCollection<Eigen::Vector2f>&, const core::LandmarkMapper&, int, int, const morphablemodel::EdgeTopology&, const fitting::ContourLandmarks&, const fitting::ModelContour&, int, cpp17::optional<int>, float, cpp17::optional<fitting::RenderingParameters>, std::vector<float>&, std::vector<float>&, std::vector<Eigen::Vector2f>&)

num_iterations: Results are good for even a few iterations. For full convergence of all parameters, it can take up to 300 iterations. In tracking, particularly if initialising with the previous frame, it works well with as low as 1 to 5 iterations. edge_topology is used for the occluding-edge face contour fitting. contour_landmarks and model_contour are used to fit the front-facing contour.

Parameters
[in]morphable_modelThe 3D Morphable Model used for the shape fitting.
[in]blendshapesA vector of blendshapes that are being fit to the landmarks in addition to the PCA model.
[in]landmarks2D landmarks from an image to fit the model to.
[in]landmark_mapperMapping info from the 2D landmark points to 3D vertex indices.
[in]image_widthWidth of the input image (needed for the camera model).
[in]image_heightHeight of the input image (needed for the camera model).
[in]edge_topologyPrecomputed edge topology of the 3D model, needed for fast edge-lookup.
[in]contour_landmarks2D image contour ids of left or right side (for example for ibug landmarks).
[in]model_contourThe model contour indices that should be considered to find the closest corresponding 3D vertex.
[in]num_iterationsNumber of iterations that the different fitting parts will be alternated for.
[in]num_shape_coefficients_to_fitHow many shape-coefficients to fit (all others will stay 0). Should be bigger than zero, or std::nullopt to fit all coefficients.
[in]lambdaRegularisation parameter of the PCA shape fitting.
Returns
The fitted model shape instance and the final pose.

◆ fit_shape_to_landmarks_linear()

std::vector< float > eos::fitting::fit_shape_to_landmarks_linear ( const morphablemodel::PcaModel shape_model,
Eigen::Matrix< float, 3, 4 >  affine_camera_matrix,
const std::vector< Eigen::Vector2f > &  landmarks,
const std::vector< int > &  vertex_ids,
Eigen::VectorXf  base_face = Eigen::VectorXf(),
float  lambda = 3.0f,
cpp17::optional< int >  num_coefficients_to_fit = cpp17::optional<int>(),
cpp17::optional< float >  detector_standard_deviation = cpp17::optional<float>(),
cpp17::optional< float >  model_standard_deviation = cpp17::optional<float>() 
)
inline

Fits the shape of a Morphable Model to given 2D landmarks (i.e. estimates the maximum likelihood solution of the shape coefficients) as proposed in [1]. It's a linear, closed-form solution fitting of the shape, with regularisation (prior towards the mean).

[1] O. Aldrian & W. Smith, Inverse Rendering of Faces with a 3D Morphable Model, PAMI 2013.

Note: Using less than the maximum number of coefficients to fit is not thoroughly tested yet and may contain an error. Note: Returns coefficients following standard normal distribution (i.e. all have similar magnitude). Why? Because we fit using the normalised basis? Note: The standard deviations given should be a vector, i.e. different for each landmark. This is not implemented yet.

Parameters
[in]shape_modelThe Morphable Model whose shape (coefficients) are estimated.
[in]affine_camera_matrixA 3x4 affine camera matrix from model to screen-space.
[in]landmarks2D landmarks from an image to fit the model to.
[in]vertex_idsThe vertex ids in the model that correspond to the 2D points.
[in]base_faceThe base or reference face from where the fitting is started. Usually this would be the models mean face, which is what will be used if the parameter is not explicitly specified.
[in]lambdaThe regularisation parameter (weight of the prior towards the mean).
[in]num_coefficients_to_fitHow many shape coefficients to fit (all others will stay 0). Should be bigger than zero, or std::nullopt to fit all coefficients.
[in]detector_standard_deviationThe standard deviation of the 2D landmarks given (e.g. of the detector used), in pixels.
[in]model_standard_deviationThe standard deviation of the 3D vertex points in the 3D model, projected to 2D (so the value is in pixels).
Returns
The estimated shape coefficients (alphas).

◆ fit_shape_to_landmarks_linear_multi()

std::vector< float > eos::fitting::fit_shape_to_landmarks_linear_multi ( const morphablemodel::PcaModel shape_model,
const std::vector< Eigen::Matrix< float, 3, 4 > > &  affine_camera_matrices,
const std::vector< std::vector< Eigen::Vector2f > > &  landmarks,
const std::vector< std::vector< int > > &  vertex_ids,
std::vector< Eigen::VectorXf >  base_faces = std::vector<Eigen::VectorXf>(),
float  lambda = 3.0f,
cpp17::optional< int >  num_coefficients_to_fit = cpp17::optional<int>(),
cpp17::optional< float >  detector_standard_deviation = cpp17::optional<float>(),
cpp17::optional< float >  model_standard_deviation = cpp17::optional<float>() 
)
inline

Fits the shape of a Morphable Model to given 2D landmarks from multiple images.

See the documentation of the single-image version, fit_shape_to_landmarks_linear(...).

Note: I think fit_shape_to_landmarks_linear(...) got updated since Philipp wrote this, so maybe we could copy some of these updates from above.

Parameters
[in]shape_modelThe Morphable Model whose shape (coefficients) are estimated.
[in]affine_camera_matricesA 3x4 affine camera matrix from model to screen-space (should probably be of type CV_32FC1 as all our calculations are done with float).
[in]landmarks2D landmarks from an image to fit the model to.
[in]vertex_idsThe vertex ids in the model that correspond to the 2D points.
[in]base_facesThe base or reference face from where the fitting is started. Usually this would be the models mean face, which is what will be used if the parameter is not explicitly specified.
[in]lambdaThe regularisation parameter (weight of the prior towards the mean). Gets normalized by the number of images given.
[in]num_coefficients_to_fitHow many shape-coefficients to fit (all others will stay 0). Should be bigger than zero, or boost::none to fit all coefficients.
[in]detector_standard_deviationThe standard deviation of the 2D landmarks given (e.g. of the detector used), in pixels.
[in]model_standard_deviationThe standard deviation of the 3D vertex points in the 3D model, projected to 2D (so the value is in pixels).
Returns
The estimated shape-coefficients (alphas).

◆ get_3x4_affine_camera_matrix()

Eigen::Matrix< float, 3, 4 > eos::fitting::get_3x4_affine_camera_matrix ( RenderingParameters  params,
int  width,
int  height 
)
inline

Creates a 3x4 affine camera matrix from given fitting parameters. The matrix transforms points directly from model-space to screen-space.

This function is mainly used since the linear shape fitting fitting::fit_shape_to_landmarks_linear expects one of these 3x4 affine camera matrices, as well as render::extract_texture.

◆ get_contour_correspondences()

std::tuple< std::vector< Eigen::Vector2f >, std::vector< Eigen::Vector4f >, std::vector< int > > eos::fitting::get_contour_correspondences ( const core::LandmarkCollection< Eigen::Vector2f > &  landmarks,
const ContourLandmarks contour_landmarks,
const ModelContour model_contour,
float  yaw_angle,
const core::Mesh mesh,
const Eigen::Matrix4f &  view_model,
const Eigen::Matrix4f &  ortho_projection,
const Eigen::Vector4f &  viewport,
float  frontal_range_threshold = 7.5f 
)
inline

Given a set of 2D image landmarks, finds the closest (in a L2 sense) 3D vertex from a list of vertices pre-defined in model_contour. landmarks can contain all landmarks, and the function will sub-select the relevant contour landmarks with the help of the given contour_landmarks. This function choses the front-facing contour and only fits this contour to the 3D model, since these correspondences are approximately static and do not move with changing pose-angle.

It's the main contour fitting function that calls all other functions.

Note: Maybe rename to find_contour_correspondences, to highlight that there is (potentially a lot) computational cost involved? Note: Does ortho_projection have to be specifically orthographic? Otherwise, if it works with perspective too, rename to just "projection".

Parameters
[in]landmarksAll image landmarks.
[in]contour_landmarks2D image contour ids of left or right side (for example for ibug landmarks).
[in]model_contourThe model contour indices that should be considered to find the closest corresponding 3D vertex.
[in]yaw_angleYaw angle of the current fitting, in degrees. The front-facing contour will be chosen depending on this yaw angle.
[in]meshThe mesh that's used to find the nearest contour points.
[in]view_modelModel-view matrix of the current fitting to project the 3D model vertices to 2D.
[in]ortho_projectionProjection matrix to project the 3D model vertices to 2D.
[in]viewportCurrent viewport to use.
[in]frontal_range_thresholdIf the yaw angle is between +- this threshold (in degrees), both contours will be selected.
Returns
A tuple with the 2D contour landmark points, the corresponding points in the 3D shape model and their vertex indices.

◆ get_corresponding_pointset()

template<class T >
auto eos::fitting::get_corresponding_pointset ( const T &  landmarks,
const core::LandmarkMapper landmark_mapper,
const morphablemodel::MorphableModel morphable_model 
)
inline

Takes a LandmarkCollection of 2D landmarks and, using the landmark_mapper, finds the corresponding 3D vertex indices and returns them, along with the coordinates of the 3D points.

The function only returns points which the landmark mapper was able to convert, and skips all points for which there is no mapping. Thus, the number of returned points might be smaller than the number of input points. All three output vectors have the same size and contain the points in the same order. landmarks can be an eos::core::LandmarkCollection<cv::Vec2f> or an rcr::LandmarkCollection<cv::Vec2f>.

Notes:

  • Split into two functions, one which maps from 2D LMs to vtx_idx and returns a reduced vec of 2D LMs. And then the other one to go from vtx_idx to a vector<Vec4f>.
  • Place in a potentially more appropriate header (shape-fitting?).
  • Could move to detail namespace or forward-declare.
  • landmarks has to be a collection of LMs, with size(), [] and Vec2f ::coordinates.
  • Probably model_points would better be a Vector3f and not in homogeneous coordinates?
Parameters
[in]landmarksA LandmarkCollection of 2D landmarks.
[in]landmark_mapperA mapper which maps the 2D landmark identifiers to 3D model vertex indices.
[in]morphable_modelModel to get the 3D point coordinates from.
Returns
A tuple of [image_points, model_points, vertex_indices].

◆ get_nearest_contour_correspondences()

std::tuple< std::vector< Eigen::Vector2f >, std::vector< Eigen::Vector4f >, std::vector< int > > eos::fitting::get_nearest_contour_correspondences ( const core::LandmarkCollection< Eigen::Vector2f > &  landmarks,
const std::vector< std::string > &  landmark_contour_identifiers,
const std::vector< int > &  model_contour_indices,
const core::Mesh mesh,
const Eigen::Matrix4f &  view_model,
const Eigen::Matrix4f &  ortho_projection,
const Eigen::Vector4f &  viewport 
)
inline

Given a set of 2D image landmarks, finds the closest (in a L2 sense) 3D vertex from a list of vertices pre-defined in model_contour. Assumes to be given contour correspondences of the front-facing contour.

Todo: Return Vector3f here instead of Vector4f?

Note: Maybe rename to find_nearest_contour_points, to highlight that there is (potentially a lot) computational cost involved? Note: Does ortho_projection have to be specifically orthographic? Otherwise, if it works with perspective too, rename to just "projection". Note: Actually, only return the vertex id, not the model point as well? Same with get_corresponding_pointset?

Parameters
[in]landmarksAll image landmarks.
[in]landmark_contour_identifiers2D image contour ids of left or right side (for example for ibug landmarks).
[in]model_contour_indicesThe model contour indices that should be considered to find the closest corresponding 3D vertex.
[in]meshThe mesh that's projected to find the nearest contour vertex.
[in]view_modelModel-view matrix of the current fitting to project the 3D model vertices to 2D.
[in]ortho_projectionProjection matrix to project the 3D model vertices to 2D.
[in]viewportCurrent viewport to use.
Returns
A tuple with the 2D contour landmark points, the corresponding points in the 3D shape model and their vertex indices.

◆ get_shape_at_point()

template<typename T >
Eigen::Vector3< T > eos::fitting::get_shape_at_point ( const eos::morphablemodel::PcaModel shape_model,
const eos::morphablemodel::Blendshapes blendshapes,
int  vertex_id,
Eigen::Map< const Eigen::VectorX< T > >  shape_coeffs,
Eigen::Map< const Eigen::VectorX< T > >  blendshape_coeffs 
)

Returns the 3D position of a single point of the 3D shape generated by the parameters given.

Parameters
[in]shape_modelA PCA 3D shape model.
[in]blendshapesA set of 3D blendshapes.
[in]vertex_idVertex id of the 3D model that should be projected.
[in]shape_coeffsA set of PCA shape coefficients used to generate the point.
[in]blendshape_coeffsA set of blendshape coefficients used to generate the point.
Returns
The 3D point.

◆ get_vertex_color_at_point()

template<typename T >
Eigen::Vector3< T > eos::fitting::get_vertex_color_at_point ( const eos::morphablemodel::PcaModel color_model,
int  vertex_id,
Eigen::Map< const Eigen::VectorX< T > >  color_coeffs 
)

Returns the colour value of a single point of the 3D model generated by the parameters given.

Parameters
[in]color_modelA PCA 3D colour (albedo) model.
[in]vertex_idVertex id of the 3D model whose colour is to be returned.
[in]color_coeffsA set of PCA colour coefficients.
Returns
The colour. As RGB? In [0, 1]?

◆ occluding_boundary_vertices()

std::vector< int > eos::fitting::occluding_boundary_vertices ( const core::Mesh mesh,
const morphablemodel::EdgeTopology edge_topology,
Eigen::Matrix3f  R,
render::ProjectionType  projection_type,
bool  perform_self_occlusion_check = true 
)
inline

Computes the vertices that lie on occluding boundaries, given a particular pose.

This algorithm computes the edges that lie on occluding boundaries of the mesh. It performs a visibility test of each vertex, and returns a list of the (unique) vertices that make up the boundary edges. An edge is defined as the line whose two adjacent faces normals flip the sign.

Parameters
[in]meshThe mesh to use.
[in]edge_topologyThe edge topology of the given mesh.
[in]RThe rotation (pose) under which the occluding boundaries should be computed.
[in]projection_typeIndicates whether the projection used is orthographic or perspective.
[in]perform_self_occlusion_checkCheck the computed boundary vertices for self-occlusion and remove these.
Returns
A vector with unique vertex id's making up the edges.

◆ save_rendering_parameters()

void eos::fitting::save_rendering_parameters ( RenderingParameters  rendering_parameters,
std::string  filename 
)
inline

Saves the rendering parameters for an image to a json file.

Parameters
[in]rendering_parametersAn instance of class RenderingParameters.
[in]filenameThe file to write.
Exceptions
std::runtime_errorif unable to open the given file for writing.

◆ select_contour()

std::pair< std::vector< std::string >, std::vector< int > > eos::fitting::select_contour ( float  yaw_angle,
const ContourLandmarks contour_landmarks,
const ModelContour model_contour,
float  frontal_range_threshold 
)
inline

Takes a set of 2D and 3D contour landmarks and a yaw angle and returns two vectors with either the right or the left 2D and 3D contour indices. This function does not establish correspondence between the 2D and 3D landmarks, it just selects the front-facing contour. The two returned vectors can thus have different size. Correspondence can be established using get_nearest_contour_correspondences().

If the yaw angle is between +-7.5°, both contours will be selected.

Parameters
[in]yaw_angleYaw angle in degrees.
[in]contour_landmarks2D image contour ids of left or right side (for example for ibug landmarks).
[in]model_contourThe model contour indices that should be used/considered to find the closest corresponding 3D vertex.
[in]frontal_range_thresholdIf the yaw angle is between +- this threshold (in degrees), both contours will be selected.
Returns
A pair with two vectors containing the selected 2D image contour landmark ids and the 3D model contour indices.

◆ tait_bryan_angles()

template<typename T >
Eigen::Matrix< T, 3, 1 > eos::fitting::tait_bryan_angles ( Eigen::Matrix< T, 3, 3 >  rotation_matrix,
Eigen::Index  axis_0,
Eigen::Index  axis_1,
Eigen::Index  axis_2 
)
inline

Computes Tait-Bryan angles (in radians) from the given rotation matrix and axes order.

Calls Eigen::Matrix3::eulerAngles(), but then swap the solution for the one where the middle (pitch) axis is constrained to -PI/2 to PI/2.