libigl v2.5.0
Loading...
Searching...
No Matches
arap_dof.h File Reference

"Fast Automatic Skinning Transformations" [Jacobson et al. 2012] More...

#include "igl_inline.h"
#include <Eigen/Dense>
#include <Eigen/Sparse>
#include "ARAPEnergyType.h"
#include <vector>
#include "arap_dof.cpp"

Go to the source code of this file.

Classes

struct  igl::ArapDOFData< LbsMatrixType, SSCALAR >
 Structure that contains fields for all precomputed data or data that needs to be remembered at update. More...
 

Namespaces

namespace  igl
 

Functions

template<typename LbsMatrixType , typename SSCALAR >
bool igl::arap_dof_precomputation (const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, const LbsMatrixType &M, const Eigen::Matrix< int, Eigen::Dynamic, 1 > &G, ArapDOFData< LbsMatrixType, SSCALAR > &data)
 Precomputes the system to optimize for "Fast Automatic Skinning Transformations" [Jacobson et al. 2012] skinning degrees of freedom optimization using as-rigid-as-possible energy.
 
template<typename LbsMatrixType , typename SSCALAR >
bool igl::arap_dof_recomputation (const Eigen::Matrix< int, Eigen::Dynamic, 1 > &fixed_dim, const Eigen::SparseMatrix< double > &A_eq, ArapDOFData< LbsMatrixType, SSCALAR > &data)
 Should always be called after arap_dof_precomputation, but may be called in between successive calls to arap_dof_update, recomputes precomputation given that there are only changes in free and fixed.
 
template<typename LbsMatrixType , typename SSCALAR >
bool igl::arap_dof_update (const ArapDOFData< LbsMatrixType, SSCALAR > &data, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &B_eq, const Eigen::MatrixXd &L0, const int max_iters, const double tol, Eigen::MatrixXd &L)
 Optimizes the transformations attached to each weight function based on precomputed system.
 

Detailed Description

"Fast Automatic Skinning Transformations" [Jacobson et al. 2012]

Arap DOF precomputation consists of two parts the computation. The first is that which depends solely on the mesh (V,F), the linear blend skinning weights (M) and the groups G. Then there's the part that depends on the previous precomputation and the list of free and fixed vertices.

Caller example:

Once:
  arap_dof_precomputation(...)

Each frame:
  while(not satisfied)
    arap_dof_update(...)
  end

The code and variables differ from the description in Section 3 of "Fast Automatic Skinning Transformations" by [Jacobson et al. 2012]

Here is a useful conversion table:

[article]                             [code]
S = \tilde{K} T                       S = CSM * Lsep
S --> R                               S --> R --shuffled--> Rxyz
Gamma_solve RT = Pi_1 \tilde{K} RT    L_part1xyz = CSolveBlock1 * Rxyz 
Pi_1 \tilde{K}                        CSolveBlock1
Peq = [T_full; P_pos]                 
T_full                                B_eq_fix <--- L0
P_pos                                 B_eq
Pi_2 * P_eq =                         Lpart2and3 = Lpart2 + Lpart3
  Pi_2_left T_full +                  Lpart3 = M_fullsolve(right) * B_eq_fix
  Pi_2_right P_pos                    Lpart2 = M_fullsolve(left) * B_eq
T = [Pi_1 Pi_2] [\tilde{K}TRT P_eq]   L = Lpart1 + Lpart2and3