libigl v2.5.0
Loading...
Searching...
No Matches
scaf.h
Go to the documentation of this file.
1// This file is part of libigl, a simple c++ geometry processing library.
2//
3// Copyright (C) 2018 Zhongshi Jiang <jiangzs@nyu.edu>
4//
5// This Source Code Form is subject to the terms of the Mozilla Public License
6// v. 2.0. If a copy of the MPL was not distributed with this file, You can
7// obtain one at http://mozilla.org/MPL/2.0/.
8#ifndef IGL_SCAF_H
9#define IGL_SCAF_H
10
11#include "../slim.h"
12#include "../igl_inline.h"
13#include "../MappingEnergyType.h"
14
15namespace igl
16{
17 namespace triangle
18 {
24 struct SCAFData
25 {
26 double scaffold_factor = 10;
29
30 // Output
31 int dim = 2;
35 double energy;
36
37 long mv_num = 0, mf_num = 0;
38 long sv_num = 0, sf_num = 0;
39 long v_num{}, f_num = 0;
41 Eigen::MatrixXd m_V;
43 Eigen::MatrixXi m_T;
44 // INTERNAL
46 Eigen::MatrixXd w_uv;
48 Eigen::MatrixXi s_T;
49 Eigen::MatrixXi w_T;
50
52 Eigen::VectorXd m_M;
54 Eigen::VectorXd s_M;
56 Eigen::VectorXd w_M;
58 double mesh_measure = 0;
59 double proximal_p = 0;
60
61 Eigen::VectorXi frame_ids;
62 Eigen::VectorXi fixed_ids;
63
64 std::map<int, Eigen::RowVectorXd> soft_cons;
65 double soft_const_p = 1e4;
66
67 Eigen::VectorXi internal_bnd;
68 Eigen::MatrixXd rect_frame_V;
69 // multi-chart support
70 std::vector<int> component_sizes;
71 std::vector<int> bnd_sizes;
72
73 // reweightedARAP interior variables.
74 bool has_pre_calc = false;
75 Eigen::SparseMatrix<double> Dx_s, Dy_s, Dz_s;
76 Eigen::SparseMatrix<double> Dx_m, Dy_m, Dz_m;
77 Eigen::MatrixXd Ri_m, Ji_m, Ri_s, Ji_s;
78 Eigen::MatrixXd W_m, W_s;
79 };
80
81
97 const Eigen::MatrixXd &V,
98 const Eigen::MatrixXi &F,
99 const Eigen::MatrixXd &V_init,
100 triangle::SCAFData &data,
101 MappingEnergyType slim_energy,
102 Eigen::VectorXi& b,
103 Eigen::MatrixXd& bc,
104 double soft_p);
105
110 IGL_INLINE Eigen::MatrixXd scaf_solve(triangle::SCAFData &data, int iter_num);
111
117 IGL_INLINE void scaf_system(triangle::SCAFData &s, Eigen::SparseMatrix<double> &L, Eigen::VectorXd &rhs);
118
119 namespace scaf
120 {
126 IGL_INLINE double compute_energy(SCAFData &s, const Eigen::MatrixXd &w_uv, bool whole);
127 }
128 }
129}
130
131#ifndef IGL_STATIC_LIBRARY
132# include "scaf.cpp"
133#endif
134
135#endif //IGL_SCAF_H
#define IGL_INLINE
Definition igl_inline.h:15
double compute_energy(SCAFData &s, const Eigen::MatrixXd &w_uv, bool whole)
Compute SCAF energy.
void scaf_precompute(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, const Eigen::MatrixXd &V_init, triangle::SCAFData &data, MappingEnergyType slim_energy, Eigen::VectorXi &b, Eigen::MatrixXd &bc, double soft_p)
Compute necessary information to start using SCAF.
void scaf_system(triangle::SCAFData &s, Eigen::SparseMatrix< double > &L, Eigen::VectorXd &rhs)
Set up the SCAF system L * uv = rhs, without solving it.
Eigen::MatrixXd scaf_solve(triangle::SCAFData &data, int iter_num)
Run iter_num iterations of SCAF, with precomputed data.
Definition AABB.h:17
MappingEnergyType
Energy Types used for Parameterization/Mapping.
Definition MappingEnergyType.h:16
@ SYMMETRIC_DIRICHLET
Definition MappingEnergyType.h:19
Use a similar interface to igl::slim Implement ready-to-use 2D version of the algorithm described in ...
Definition scaf.h:25
Eigen::MatrixXd Ji_s
Definition scaf.h:77
Eigen::SparseMatrix< double > Dx_m
Definition scaf.h:76
long sv_num
Definition scaf.h:38
Eigen::VectorXi internal_bnd
Definition scaf.h:67
int dim
Definition scaf.h:31
Eigen::MatrixXd Ri_s
Definition scaf.h:77
double total_energy
scaffold + isometric
Definition scaf.h:33
bool has_pre_calc
Definition scaf.h:74
Eigen::SparseMatrix< double > Dz_s
Definition scaf.h:75
igl::MappingEnergyType scaf_energy
Definition scaf.h:27
Eigen::SparseMatrix< double > Dy_s
Definition scaf.h:75
Eigen::VectorXi fixed_ids
Definition scaf.h:62
Eigen::MatrixXd m_V
input initial mesh V
Definition scaf.h:41
Eigen::MatrixXd Ri_m
Definition scaf.h:77
Eigen::MatrixXi s_T
scaffold domain tets: scaffold tets
Definition scaf.h:48
long mv_num
Definition scaf.h:37
Eigen::MatrixXd w_uv
whole domain uv: mesh + free vertices
Definition scaf.h:46
Eigen::SparseMatrix< double > Dz_m
Definition scaf.h:76
Eigen::MatrixXd rect_frame_V
Definition scaf.h:68
Eigen::MatrixXd W_s
Definition scaf.h:78
Eigen::SparseMatrix< double > Dx_s
Definition scaf.h:75
double soft_const_p
Definition scaf.h:65
std::map< int, Eigen::RowVectorXd > soft_cons
Definition scaf.h:64
Eigen::MatrixXd Ji_m
Definition scaf.h:77
Eigen::MatrixXd W_m
Definition scaf.h:78
Eigen::MatrixXi w_T
Definition scaf.h:49
Eigen::VectorXi frame_ids
Definition scaf.h:61
Eigen::VectorXd m_M
mesh area or volume
Definition scaf.h:52
double mesh_measure
area or volume
Definition scaf.h:58
long sf_num
Definition scaf.h:38
Eigen::VectorXd w_M
area/volume weights for whole
Definition scaf.h:56
long mf_num
Definition scaf.h:37
Eigen::VectorXd s_M
scaffold area or volume
Definition scaf.h:54
Eigen::MatrixXi m_T
input initial mesh F/T
Definition scaf.h:43
Eigen::SparseMatrix< double > Dy_m
Definition scaf.h:76
double energy
objective value
Definition scaf.h:35
igl::MappingEnergyType slim_energy
Definition scaf.h:28
double scaffold_factor
Definition scaf.h:26
long v_num
Definition scaf.h:39
double proximal_p
Definition scaf.h:59
std::vector< int > component_sizes
Definition scaf.h:70
std::vector< int > bnd_sizes
Definition scaf.h:71
long f_num
Definition scaf.h:39