Prusa Slicer 2.6.0
Loading...
Searching...
No Matches
igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF > Class Template Reference
+ Collaboration diagram for igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >:

Public Member Functions

IGL_INLINE void SolvePoisson (Eigen::VectorXd Stiffness, double vector_field_scale=0.1f, double grid_res=1.f, bool direct_round=true, int localIter=0, bool _integer_rounding=true, bool _singularity_rounding=true, std::vector< int > roundVertices=std::vector< int >(), std::vector< std::vector< int > > hardFeatures=std::vector< std::vector< int > >())
 
IGL_INLINE PoissonSolver (const Eigen::PlainObjectBase< DerivedV > &_V, const Eigen::PlainObjectBase< DerivedF > &_F, const Eigen::PlainObjectBase< DerivedV > &_Vcut, const Eigen::PlainObjectBase< DerivedF > &_Fcut, const Eigen::PlainObjectBase< DerivedF > &_TT, const Eigen::PlainObjectBase< DerivedF > &_TTi, const Eigen::PlainObjectBase< DerivedV > &_PD1, const Eigen::PlainObjectBase< DerivedV > &_PD2, const Eigen::Matrix< int, Eigen::Dynamic, 1 > &_Handle_Singular, const MeshSystemInfo &_Handle_SystemInfo)
 
IGL_INLINE std::complex< double > GetRotationComplex (int interval)
 START COMMON MATH FUNCTIONS return the complex encoding the rotation for a given missmatch interval.
 
IGL_INLINE void AddFixedVertex (int v)
 END COMMON MATH FUNCTIONS.
 
IGL_INLINE void FindFixedVertField ()
 find vertex to fix in case we're using a vector field NB: multiple components not handled
 
IGL_INLINE void FindFixedVert ()
 find hard constraint depending if using or not a vector field
 
IGL_INLINE int GetFirstVertexIndex (int v)
 
IGL_INLINE void FixBlockedVertex ()
 fix the vertices which are flagged as fixed
 
IGL_INLINE void AddSingularityRound ()
 END FIXING VERTICES.
 
IGL_INLINE void AddToRoundVertices (std::vector< int > ids)
 
IGL_INLINE void BuildLaplacianMatrix (double vfscale=1)
 START GENERIC SYSTEM FUNCTIONS.
 
IGL_INLINE void FindSizes ()
 find different sized of the system
 
IGL_INLINE void AllocateSystem ()
 
IGL_INLINE void InitMatrix ()
 intitialize the whole matrix
 
IGL_INLINE void MapCoords ()
 map back coordinates after that the system has been solved
 
IGL_INLINE void BuildSeamConstraintsExplicitTranslation ()
 END GENERIC SYSTEM FUNCTIONS.
 
IGL_INLINE void BuildUserDefinedConstraints ()
 set the constraints for the inter-range cuts
 
IGL_INLINE void MixedIntegerSolve (double cone_grid_res=1, bool direct_round=true, int localIter=0)
 call of the mixed integer solver
 
IGL_INLINE void clearUserConstraint ()
 
IGL_INLINE void addSharpEdgeConstraint (int fid, int vid)
 

Public Attributes

const Eigen::PlainObjectBase< DerivedV > & V
 
const Eigen::PlainObjectBase< DerivedF > & F
 
const Eigen::PlainObjectBase< DerivedV > & Vcut
 
const Eigen::PlainObjectBase< DerivedF > & Fcut
 
const Eigen::PlainObjectBase< DerivedF > & TT
 
const Eigen::PlainObjectBase< DerivedF > & TTi
 
const Eigen::PlainObjectBase< DerivedV > & PD1
 
const Eigen::PlainObjectBase< DerivedV > & PD2
 
const Eigen::Matrix< int, Eigen::Dynamic, 1 > & Handle_Singular
 
const MeshSystemInfoHandle_SystemInfo
 
Eigen::VectorXd Handle_Stiffness
 
std::vector< std::vector< int > > VF
 
std::vector< std::vector< int > > VFi
 
Eigen::MatrixXd UV
 
Eigen::MatrixXd WUV
 
Eigen::MatrixXd UV_out
 
Eigen::SparseMatrix< double > Lhs
 
Eigen::SparseMatrix< double > Constraints
 
Eigen::VectorXd rhs
 
Eigen::VectorXd constraints_rhs
 
std::vector< double > X
 vector of unknowns
 
unsigned int n_fixed_vars
 number of fixed vertex
 
unsigned int n_vert_vars
 the number of REAL variables for vertices
 
unsigned int num_total_vars
 total number of variables of the system, do not consider constraints, but consider integer vars
 
unsigned int n_integer_vars
 the total number of integer variables
 
unsigned int num_cut_constraint
 CONSTRAINT PART number of cuts constraints.
 
unsigned int num_userdefined_constraint
 
unsigned int num_constraint_equations
 total number of constraints equations
 
std::vector< int > Hard_constraints
 vector of blocked vertices
 
std::vector< int > ids_to_round
 vector of indexes to round
 
std::vector< std::vector< int > > userdefined_constraints
 vector of indexes to round
 
bool integer_rounding
 boolean that is true if rounding to integer is needed
 

Detailed Description

template<typename DerivedV, typename DerivedF>
class igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >

Constructor & Destructor Documentation

◆ PoissonSolver()

template<typename DerivedV , typename DerivedF >
IGL_INLINE igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::PoissonSolver ( const Eigen::PlainObjectBase< DerivedV > &  _V,
const Eigen::PlainObjectBase< DerivedF > &  _F,
const Eigen::PlainObjectBase< DerivedV > &  _Vcut,
const Eigen::PlainObjectBase< DerivedF > &  _Fcut,
const Eigen::PlainObjectBase< DerivedF > &  _TT,
const Eigen::PlainObjectBase< DerivedF > &  _TTi,
const Eigen::PlainObjectBase< DerivedV > &  _PD1,
const Eigen::PlainObjectBase< DerivedV > &  _PD2,
const Eigen::Matrix< int, Eigen::Dynamic, 1 > &  _Handle_Singular,
const MeshSystemInfo _Handle_SystemInfo 
)
679 :
680V(_V),
681F(_F),
682Vcut(_Vcut),
683Fcut(_Fcut),
684TT(_TT),
685TTi(_TTi),
686PD1(_PD1),
687PD2(_PD2),
688Handle_Singular(_Handle_Singular),
689Handle_SystemInfo(_Handle_SystemInfo)
690{
691 UV = Eigen::MatrixXd(V.rows(),2);
692 WUV = Eigen::MatrixXd(F.rows(),6);
693 UV_out = Eigen::MatrixXd(Vcut.rows(),2);
695}
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Index rows() const
Definition PlainObjectBase.h:151
const Eigen::PlainObjectBase< DerivedV > & Vcut
Definition miq.cpp:163
std::vector< std::vector< int > > VF
Definition miq.cpp:175
Eigen::MatrixXd WUV
Definition miq.cpp:181
Eigen::MatrixXd UV
Definition miq.cpp:177
const Eigen::PlainObjectBase< DerivedF > & Fcut
Definition miq.cpp:164
const Eigen::PlainObjectBase< DerivedV > & V
Definition miq.cpp:161
const Eigen::PlainObjectBase< DerivedV > & PD1
Definition miq.cpp:167
const MeshSystemInfo & Handle_SystemInfo
Definition miq.cpp:171
const Eigen::PlainObjectBase< DerivedF > & TT
Definition miq.cpp:165
const Eigen::Matrix< int, Eigen::Dynamic, 1 > & Handle_Singular
Definition miq.cpp:169
const Eigen::PlainObjectBase< DerivedF > & F
Definition miq.cpp:162
Eigen::MatrixXd UV_out
Definition miq.cpp:183
const Eigen::PlainObjectBase< DerivedF > & TTi
Definition miq.cpp:166
std::vector< std::vector< int > > VFi
Definition miq.cpp:176
const Eigen::PlainObjectBase< DerivedV > & PD2
Definition miq.cpp:168
IGL_INLINE void vertex_triangle_adjacency(const typename DerivedF::Scalar n, const Eigen::MatrixBase< DerivedF > &F, std::vector< std::vector< VFType > > &VF, std::vector< std::vector< VFiType > > &VFi)
Definition vertex_triangle_adjacency.cpp:12

References igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::F, Eigen::PlainObjectBase< Derived >::rows(), igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::UV, igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::UV_out, igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::V, igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::Vcut, igl::vertex_triangle_adjacency(), igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::VF, igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::VFi, and igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::WUV.

+ Here is the call graph for this function:

Member Function Documentation

◆ AddFixedVertex()

template<typename DerivedV , typename DerivedF >
IGL_INLINE void igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::AddFixedVertex ( int  v)

END COMMON MATH FUNCTIONS.

START FIXING VERTICES set a given vertex as fixed

720{
721 n_fixed_vars++;
722 Hard_constraints.push_back(v);
723}
unsigned int n_fixed_vars
number of fixed vertex
Definition miq.cpp:195
std::vector< int > Hard_constraints
vector of blocked vertices
Definition miq.cpp:219

◆ addSharpEdgeConstraint()

template<typename DerivedV , typename DerivedF >
IGL_INLINE void igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::addSharpEdgeConstraint ( int  fid,
int  vid 
)
1134{
1135 // prepare constraint
1136 std::vector<int> c(Handle_SystemInfo.num_vert_variables*2 + 1);
1137
1138 for (size_t i = 0; i < c.size(); ++i)
1139 {
1140 c[i] = 0;
1141 }
1142
1143 int v1 = Fcut(fid,vid);
1144 int v2 = Fcut(fid,(vid+1)%3);
1145
1147 e = e.normalized();
1148
1149 double d1 = fabs(e.dot(PD1.row(fid).normalized()));
1150 double d2 = fabs(e.dot(PD2.row(fid).normalized()));
1151
1152 int offset = 0;
1153
1154 if (d1>d2)
1155 offset = 1;
1156
1157 ids_to_round.push_back((v1 * 2) + offset);
1158 ids_to_round.push_back((v2 * 2) + offset);
1159
1160 // add constraint
1161 c[(v1 * 2) + offset] = 1;
1162 c[(v2 * 2) + offset] = -1;
1163
1164 // add to the user-defined constraints
1166 userdefined_constraints.push_back(c);
1167
1168}
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:180
unsigned int num_userdefined_constraint
Definition miq.cpp:213
std::vector< int > ids_to_round
vector of indexes to round
Definition miq.cpp:222
std::vector< std::vector< int > > userdefined_constraints
vector of indexes to round
Definition miq.cpp:225
int num_vert_variables
Definition miq.cpp:68
void offset(Slic3r::ExPolygon &sh, coord_t distance, const PolygonTag &)
Definition geometries.hpp:132

◆ AddSingularityRound()

template<typename DerivedV , typename DerivedF >
IGL_INLINE void igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::AddSingularityRound

END FIXING VERTICES.

HANDLING SINGULARITY

805{
806 for (unsigned int v=0;v<V.rows();v++)
807 {
808 if (Handle_Singular(v))
809 {
810 int index0=GetFirstVertexIndex(v);
811 ids_to_round.push_back( index0*2 );
812 ids_to_round.push_back((index0*2)+1);
813 }
814 }
815}
IGL_INLINE int GetFirstVertexIndex(int v)
Definition miq.cpp:760

◆ AddToRoundVertices()

template<typename DerivedV , typename DerivedF >
IGL_INLINE void igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::AddToRoundVertices ( std::vector< int >  ids)
819{
820 for (size_t i = 0; i < ids.size(); ++i)
821 {
822 if (ids[i] < 0 || ids[i] >= V.rows())
823 std::cerr << "WARNING: Ignored round vertex constraint, vertex " << ids[i] << " does not exist in the mesh." << std::endl;
824 int index0 = GetFirstVertexIndex(ids[i]);
825 ids_to_round.push_back( index0*2 );
826 ids_to_round.push_back((index0*2)+1);
827 }
828}

◆ AllocateSystem()

template<typename DerivedV , typename DerivedF >
IGL_INLINE void igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::AllocateSystem
909{
912 rhs.resize(n_vert_vars * 2);
914
915 printf("\n INITIALIZED SPARSE MATRIX OF %d x %d \n",n_vert_vars*2, n_vert_vars*2);
916 printf("\n INITIALIZED SPARSE MATRIX OF %d x %d \n",num_constraint_equations, num_total_vars);
917 printf("\n INITIALIZED VECTOR OF %d x 1 \n",n_vert_vars*2);
918 printf("\n INITIALIZED VECTOR OF %d x 1 \n",num_constraint_equations);
919}
void resize(Index rows, Index cols)
Definition SparseMatrix.h:621
Eigen::SparseMatrix< double > Constraints
Definition miq.cpp:187
Eigen::VectorXd rhs
Definition miq.cpp:188
unsigned int n_vert_vars
the number of REAL variables for vertices
Definition miq.cpp:198
Eigen::VectorXd constraints_rhs
Definition miq.cpp:189
unsigned int num_total_vars
total number of variables of the system, do not consider constraints, but consider integer vars
Definition miq.cpp:202
unsigned int num_constraint_equations
total number of constraints equations
Definition miq.cpp:216
Eigen::SparseMatrix< double > Lhs
Definition miq.cpp:186

◆ BuildLaplacianMatrix()

template<typename DerivedV , typename DerivedF >
IGL_INLINE void igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::BuildLaplacianMatrix ( double  vfscale = 1)

START GENERIC SYSTEM FUNCTIONS.

Compute LHS

Compute RHS

833{
834 Eigen::VectorXi idx = igl::LinSpaced<Eigen::VectorXi >(Vcut.rows(), 0, 2*Vcut.rows()-2);
835 Eigen::VectorXi idx2 = igl::LinSpaced<Eigen::VectorXi >(Vcut.rows(), 1, 2*Vcut.rows()-1);
836
837 // get gradient matrix
839 igl::grad(Vcut, Fcut, G);
840
841 // get triangle weights
842 Eigen::VectorXd dblA(Fcut.rows());
843 igl::doublearea(Vcut, Fcut, dblA);
844
845 // compute intermediate result
847 G2 = G.transpose() * dblA.replicate<3,1>().asDiagonal() * Handle_Stiffness.replicate<3,1>().asDiagonal();
848
851 Cotmatrix = 0.5 * G2 * G;
852 igl::slice_into(Cotmatrix, idx, idx, Lhs);
853 igl::slice_into(Cotmatrix, idx2, idx2, Lhs);
854
856 // reshape nrosy vectors
857 const Eigen::MatrixXd u = Eigen::Map<const Eigen::MatrixXd>(PD1.data(),Fcut.rows()*3,1); // this mimics a reshape at the cost of a copy.
858 const Eigen::MatrixXd v = Eigen::Map<const Eigen::MatrixXd>(PD2.data(),Fcut.rows()*3,1); // this mimics a reshape at the cost of a copy.
859
860 // multiply with weights
861 Eigen::VectorXd rhs1 = G2 * u * 0.5 * vfscale;
862 Eigen::VectorXd rhs2 = -G2 * v * 0.5 * vfscale;
863 igl::slice_into(rhs1, idx, 1, rhs);
864 igl::slice_into(rhs2, idx2, 1, rhs);
865}
A matrix or vector expression mapping an existing array of data.
Definition Map.h:96
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const Scalar * data() const
Definition PlainObjectBase.h:255
TransposeReturnType transpose()
Definition SparseMatrixBase.h:349
A versatible sparse matrix representation.
Definition SparseMatrix.h:98
Eigen::VectorXd Handle_Stiffness
Definition miq.cpp:174
IGL_INLINE void slice_into(const Eigen::SparseMatrix< T > &X, const Eigen::Matrix< int, Eigen::Dynamic, 1 > &R, const Eigen::Matrix< int, Eigen::Dynamic, 1 > &C, Eigen::SparseMatrix< T > &Y)
Definition slice_into.cpp:16
IGL_INLINE void doublearea(const Eigen::MatrixBase< DerivedV > &V, const Eigen::MatrixBase< DerivedF > &F, Eigen::PlainObjectBase< DeriveddblA > &dblA)
Definition doublearea.cpp:17
IGL_INLINE void grad(const Eigen::PlainObjectBase< DerivedV > &V, const Eigen::PlainObjectBase< DerivedF > &F, Eigen::SparseMatrix< typename DerivedV::Scalar > &G, bool uniform=false)
Definition grad.cpp:227

References igl::doublearea(), igl::grad(), igl::slice_into(), and Eigen::SparseMatrixBase< Derived >::transpose().

+ Here is the call graph for this function:

◆ BuildSeamConstraintsExplicitTranslation()

template<typename DerivedV , typename DerivedF >
IGL_INLINE void igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::BuildSeamConstraintsExplicitTranslation

END GENERIC SYSTEM FUNCTIONS.

set the constraints for the inter-range cuts

current constraint row

get the integer variable

963{
965 int constr_row = 0;
966
967 for (unsigned int i=0; i<num_cut_constraint; i++)
968 {
969 unsigned char interval = Handle_SystemInfo.EdgeSeamInfo[i].MMatch;
970 if (interval==1)
971 interval=3;
972 else
973 if(interval==3)
974 interval=1;
975
976 int p0 = Handle_SystemInfo.EdgeSeamInfo[i].v0;
977 int p0p = Handle_SystemInfo.EdgeSeamInfo[i].v0p;
978
979 std::complex<double> rot = GetRotationComplex(interval);
980
982 int integerVar = n_vert_vars + Handle_SystemInfo.EdgeSeamInfo[i].integerVar;
983
985 {
986 ids_to_round.push_back(integerVar*2);
987 ids_to_round.push_back(integerVar*2+1);
988 }
989
990 // cross boundary compatibility conditions
991 Constraints.coeffRef(constr_row, 2*p0) += rot.real();
992 Constraints.coeffRef(constr_row, 2*p0+1) += -rot.imag();
993 Constraints.coeffRef(constr_row+1, 2*p0) += rot.imag();
994 Constraints.coeffRef(constr_row+1, 2*p0+1) += rot.real();
995
996 Constraints.coeffRef(constr_row, 2*p0p) += -1;
997 Constraints.coeffRef(constr_row+1, 2*p0p+1) += -1;
998
999 Constraints.coeffRef(constr_row, 2*integerVar) += 1;
1000 Constraints.coeffRef(constr_row+1, 2*integerVar+1) += 1;
1001
1002 constraints_rhs[constr_row] = 0;
1003 constraints_rhs[constr_row+1] = 0;
1004
1005 constr_row += 2;
1006 }
1007
1008}
Scalar & coeffRef(Index row, Index col)
Definition SparseMatrix.h:206
IGL_INLINE std::complex< double > GetRotationComplex(int interval)
START COMMON MATH FUNCTIONS return the complex encoding the rotation for a given missmatch interval.
Definition miq.cpp:701
bool integer_rounding
boolean that is true if rounding to integer is needed
Definition miq.cpp:228
unsigned int num_cut_constraint
CONSTRAINT PART number of cuts constraints.
Definition miq.cpp:210
std::vector< SeamInfo > EdgeSeamInfo
this are used for drawing purposes
Definition miq.cpp:72

◆ BuildUserDefinedConstraints()

template<typename DerivedV , typename DerivedF >
IGL_INLINE void igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::BuildUserDefinedConstraints

set the constraints for the inter-range cuts

the user defined constraints are at the end

current constraint row

1013{
1015 int offset_row = num_cut_constraint*2 + n_fixed_vars*2;
1016
1018 int constr_row = offset_row;
1019
1021
1022 for (unsigned int i=0; i<num_userdefined_constraint; i++)
1023 {
1024 for (unsigned int j=0; j<userdefined_constraints[i].size()-1; ++j)
1025 {
1026 Constraints.coeffRef(constr_row, j) = userdefined_constraints[i][j];
1027 }
1028
1030 constr_row +=1;
1031 }
1032}

◆ clearUserConstraint()

template<typename DerivedV , typename DerivedF >
IGL_INLINE void igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::clearUserConstraint
1127{
1130}

◆ FindFixedVert()

template<typename DerivedV , typename DerivedF >
IGL_INLINE void igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::FindFixedVert

find hard constraint depending if using or not a vector field

754{
755 Hard_constraints.clear();
757}
IGL_INLINE void FindFixedVertField()
find vertex to fix in case we're using a vector field NB: multiple components not handled
Definition miq.cpp:728

◆ FindFixedVertField()

template<typename DerivedV , typename DerivedF >
IGL_INLINE void igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::FindFixedVertField

find vertex to fix in case we're using a vector field NB: multiple components not handled

if anything fixed fix the first

729{
730 Hard_constraints.clear();
731
732 n_fixed_vars=0;
733 //fix the first singularity
734 for (unsigned int v=0;v<V.rows();v++)
735 {
736 if (Handle_Singular(v))
737 {
739 UV.row(v) << 0,0;
740 return;
741 }
742 }
743
746 UV.row(0) << 0,0;
747 std::cerr << "No vertices to fix, I am fixing the first vertex to the origin!" << std::endl;
748}
IGL_INLINE void AddFixedVertex(int v)
END COMMON MATH FUNCTIONS.
Definition miq.cpp:719

◆ FindSizes()

template<typename DerivedV , typename DerivedF >
IGL_INLINE void igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::FindSizes

find different sized of the system

find the vertex that need to be fixed

REAL PART

INTEGER PART the total number of integer variables

CONSTRAINT PART

total variable of the system

initialize matrix size

870{
873
876
880
883
885
888
890
891 if (DEBUGPRINT) printf("\n*** SYSTEM VARIABLES *** \n");
892 if (DEBUGPRINT) printf("* NUM REAL VERTEX VARIABLES %d \n",n_vert_vars);
893
894 if (DEBUGPRINT) printf("\n*** INTEGER VARIABLES *** \n");
895 if (DEBUGPRINT) printf("* NUM INTEGER VARIABLES %d \n",(int)n_integer_vars);
896
897 if (DEBUGPRINT) printf("\n*** CONSTRAINTS *** \n ");
898 if (DEBUGPRINT) printf("* NUM FIXED CONSTRAINTS %d\n",n_fixed_vars);
899 if (DEBUGPRINT) printf("* NUM CUTS CONSTRAINTS %d\n",num_cut_constraint);
900 if (DEBUGPRINT) printf("* NUM USER DEFINED CONSTRAINTS %d\n",num_userdefined_constraint);
901
902 if (DEBUGPRINT) printf("\n*** TOTAL SIZE *** \n");
903 if (DEBUGPRINT) printf("* TOTAL VARIABLE SIZE (WITH INTEGER TRASL) %d \n",num_total_vars);
904 if (DEBUGPRINT) printf("* TOTAL CONSTRAINTS %d \n",num_constraint_equations);
905}
unsigned int n_integer_vars
the total number of integer variables
Definition miq.cpp:206
IGL_INLINE void FindFixedVert()
find hard constraint depending if using or not a vector field
Definition miq.cpp:753
#define DEBUGPRINT
Definition miq.cpp:45
int num_integer_cuts
num of integer for cuts
Definition miq.cpp:70

References DEBUGPRINT.

◆ FixBlockedVertex()

template<typename DerivedV , typename DerivedF >
IGL_INLINE void igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::FixBlockedVertex

fix the vertices which are flagged as fixed

get first index of the vertex that must blocked

multiply times 2 because of uv

find the first free row to add the constraint

add fixing constraint LHS

add fixing constraint RHS

768{
769 int offset_row = num_cut_constraint*2;
770
771 unsigned int constr_num = 0;
772 for (unsigned int i=0;i<Hard_constraints.size();i++)
773 {
774 int v = Hard_constraints[i];
775
777 //int index=v->vertex_index[0];
778 int index = GetFirstVertexIndex(v);
779
781 int indexvert = index*2;
782
784 int indexRow = (offset_row+constr_num*2);
785 int indexCol = indexRow;
786
788 Constraints.coeffRef(indexRow, indexvert) += 1;
789 Constraints.coeffRef(indexRow+1,indexvert+1) += 1;
790
792 constraints_rhs[indexCol] = UV(v,0);
793 constraints_rhs[indexCol+1] = UV(v,1);
794
795 constr_num++;
796 }
797 assert(constr_num==n_fixed_vars);
798}

◆ GetFirstVertexIndex()

template<typename DerivedV , typename DerivedF >
IGL_INLINE int igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::GetFirstVertexIndex ( int  v)
761{
762 return Fcut(VF[v][0],VFi[v][0]);
763}

◆ GetRotationComplex()

template<typename DerivedV , typename DerivedF >
IGL_INLINE std::complex< double > igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::GetRotationComplex ( int  interval)

START COMMON MATH FUNCTIONS return the complex encoding the rotation for a given missmatch interval.

702{
703 assert((interval>=0)&&(interval<4));
704
705 switch(interval)
706 {
707 case 0:return std::complex<double>(1,0);
708 case 1:return std::complex<double>(0,1);
709 case 2:return std::complex<double>(-1,0);
710 default:return std::complex<double>(0,-1);
711 }
712}

◆ InitMatrix()

template<typename DerivedV , typename DerivedF >
IGL_INLINE void igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::InitMatrix

intitialize the whole matrix

924{
925 FindSizes();
927}
IGL_INLINE void FindSizes()
find different sized of the system
Definition miq.cpp:869
IGL_INLINE void AllocateSystem()
Definition miq.cpp:908

◆ MapCoords()

template<typename DerivedV , typename DerivedF >
IGL_INLINE void igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::MapCoords

map back coordinates after that the system has been solved

map coords to faces

then get U and V coords

933{
935 for (unsigned int f=0;f<Fcut.rows();f++)
936 {
937
938 for (int k=0;k<3;k++)
939 {
940 //get the index of the variable in the system
941 int indexUV = Fcut(f,k);
943 double U=X[indexUV*2];
944 double V=X[indexUV*2+1];
945
946 WUV(f,k*2 + 0) = U;
947 WUV(f,k*2 + 1) = V;
948 }
949
950 }
951
952 for(int i = 0; i < Vcut.rows(); i++){
953 UV_out(i,0) = X[i*2];
954 UV_out(i,1) = X[i*2+1];
955 }
956}
std::vector< double > X
vector of unknowns
Definition miq.cpp:191

◆ MixedIntegerSolve()

template<typename DerivedV , typename DerivedF >
IGL_INLINE void igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::MixedIntegerSolve ( double  cone_grid_res = 1,
bool  direct_round = true,
int  localIter = 0 
)

call of the mixed integer solver

variables part

matrix A

constraints part

add penalization term for integer variables

1039{
1040 X = std::vector<double>((n_vert_vars+n_integer_vars)*2);
1041 if (DEBUGPRINT)
1042 printf("\n ALLOCATED X \n");
1043
1045 int ScalarSize = n_vert_vars*2;
1046 int SizeMatrix = (n_vert_vars+n_integer_vars)*2;
1047
1049 gmm::col_matrix< gmm::wsvector< double > > A(SizeMatrix,SizeMatrix); // lhs matrix variables
1050
1052 int CsizeX = num_constraint_equations;
1053 int CsizeY = SizeMatrix+1;
1054 gmm::row_matrix< gmm::wsvector< double > > C(CsizeX,CsizeY); // constraints
1055
1056 if (DEBUGPRINT)
1057 printf("\n ALLOCATED QMM STRUCTURES \n");
1058
1059 std::vector<double> B(SizeMatrix,0); // rhs
1060
1061 if (DEBUGPRINT)
1062 printf("\n ALLOCATED RHS STRUCTURES \n");
1063
1065 for (int k=0; k < Lhs.outerSize(); ++k){
1066 for (Eigen::SparseMatrix<double>::InnerIterator it(Lhs,k); it; ++it){
1067 int row = it.row();
1068 int col = it.col();
1069 A(row, col) += it.value();
1070 }
1071 }
1073 for (int k=0; k < Constraints.outerSize(); ++k){
1075 int row = it.row();
1076 int col = it.col();
1077 C(row, col) += it.value();
1078 }
1079 }
1080
1081 if (DEBUGPRINT)
1082 printf("\n SET %d INTEGER VALUES \n",n_integer_vars);
1083
1085 double penalization = 0.000001;
1086 int offline_index = ScalarSize;
1087 for(unsigned int i = 0; i < (n_integer_vars)*2; ++i)
1088 {
1089 int index=offline_index+i;
1090 A(index,index)=penalization;
1091 }
1092
1093 if (DEBUGPRINT)
1094 printf("\n SET RHS \n");
1095
1096 // copy RHS
1097 for(int i = 0; i < (int)ScalarSize; ++i)
1098 {
1099 B[i] = rhs[i] * cone_grid_res;
1100 }
1101
1102 // copy constraint RHS
1103 if (DEBUGPRINT)
1104 printf("\n SET %d CONSTRAINTS \n",num_constraint_equations);
1105
1106 for(unsigned int i = 0; i < num_constraint_equations; ++i)
1107 {
1108 C(i, SizeMatrix) = -constraints_rhs[i] * cone_grid_res;
1109 }
1110
1111 COMISO::ConstrainedSolver solver;
1112
1113 solver.misolver().set_local_iters(localIter);
1114
1115 solver.misolver().set_direct_rounding(direct_round);
1116
1117 std::sort(ids_to_round.begin(),ids_to_round.end());
1118 std::vector<int>::iterator new_end=std::unique(ids_to_round.begin(),ids_to_round.end());
1119 int dist=distance(ids_to_round.begin(),new_end);
1120 ids_to_round.resize(dist);
1121
1122 solver.solve( C, A, X, B, ids_to_round, 0.0, false, false);
1123}
EIGEN_DEVICE_FUNC RowXpr row(Index i)
This is the const version of row(). *‍/.
Definition BlockMethods.h:859
EIGEN_DEVICE_FUNC ColXpr col(Index i)
This is the const version of col().
Definition BlockMethods.h:838
Definition SparseCompressedBase.h:137
Index outerSize() const
Definition SparseMatrix.h:143
T dist(const boost::polygon::point_data< T > &p1, const boost::polygon::point_data< T > &p2)
Definition Geometry.cpp:280
double distance(const P &p1, const P &p2)
Definition geometry_traits.hpp:329

References col(), DEBUGPRINT, and row().

+ Here is the call graph for this function:

◆ SolvePoisson()

template<typename DerivedV , typename DerivedF >
IGL_INLINE void igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::SolvePoisson ( Eigen::VectorXd  Stiffness,
double  vector_field_scale = 0.1f,
double  grid_res = 1.f,
bool  direct_round = true,
int  localIter = 0,
bool  _integer_rounding = true,
bool  _singularity_rounding = true,
std::vector< int >  roundVertices = std::vector<int>(),
std::vector< std::vector< int > >  hardFeatures = std::vector<std::vector<int> >() 
)

Initializing Matrix

initialize the matrix ALLOCATING SPACE

build the laplacian system

605{
606 Handle_Stiffness = Stiffness;
607
608 //initialization of flags and data structures
609 integer_rounding=_integer_rounding;
610
611 ids_to_round.clear();
612
614 // copy the user constraints number
615 for (size_t i = 0; i < hardFeatures.size(); ++i)
616 {
617 addSharpEdgeConstraint(hardFeatures[i][0],hardFeatures[i][1]);
618 }
619
621
622 int t0=clock();
623
625 InitMatrix();
626 if (DEBUGPRINT)
627 printf("\n ALLOCATED THE MATRIX \n");
628
630 BuildLaplacianMatrix(vector_field_scale);
631
632 // add seam constraints
634
635 // add user defined constraints
637
640
641 if (DEBUGPRINT)
642 printf("\n BUILT THE MATRIX \n");
643
645 AddToRoundVertices(roundVertices);
646
647 if (_singularity_rounding)
649
650 int t1=clock();
651 if (DEBUGPRINT) printf("\n time:%d \n",t1-t0);
652 if (DEBUGPRINT) printf("\n SOLVING \n");
653
654 MixedIntegerSolve(grid_res,direct_round,localIter);
655
656 int t2=clock();
657 if (DEBUGPRINT) printf("\n time:%d \n",t2-t1);
658 if (DEBUGPRINT) printf("\n ASSIGNING COORDS \n");
659
660 MapCoords();
661
662 int t3=clock();
663 if (DEBUGPRINT) printf("\n time:%d \n",t3-t2);
664 if (DEBUGPRINT) printf("\n FINISHED \n");
665}
IGL_INLINE void AddToRoundVertices(std::vector< int > ids)
Definition miq.cpp:818
IGL_INLINE void BuildLaplacianMatrix(double vfscale=1)
START GENERIC SYSTEM FUNCTIONS.
Definition miq.cpp:832
IGL_INLINE void addSharpEdgeConstraint(int fid, int vid)
Definition miq.cpp:1133
IGL_INLINE void BuildUserDefinedConstraints()
set the constraints for the inter-range cuts
Definition miq.cpp:1012
IGL_INLINE void AddSingularityRound()
END FIXING VERTICES.
Definition miq.cpp:804
IGL_INLINE void BuildSeamConstraintsExplicitTranslation()
END GENERIC SYSTEM FUNCTIONS.
Definition miq.cpp:962
IGL_INLINE void MixedIntegerSolve(double cone_grid_res=1, bool direct_round=true, int localIter=0)
call of the mixed integer solver
Definition miq.cpp:1036
IGL_INLINE void clearUserConstraint()
Definition miq.cpp:1126
IGL_INLINE void FixBlockedVertex()
fix the vertices which are flagged as fixed
Definition miq.cpp:767
IGL_INLINE void MapCoords()
map back coordinates after that the system has been solved
Definition miq.cpp:932
IGL_INLINE void InitMatrix()
intitialize the whole matrix
Definition miq.cpp:923

References DEBUGPRINT.

Referenced by igl::copyleft::comiso::MIQ_class< DerivedV, DerivedF, DerivedU >::MIQ_class().

+ Here is the caller graph for this function:

Member Data Documentation

◆ Constraints

template<typename DerivedV , typename DerivedF >
Eigen::SparseMatrix<double> igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::Constraints

◆ constraints_rhs

template<typename DerivedV , typename DerivedF >
Eigen::VectorXd igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::constraints_rhs

◆ F

template<typename DerivedV , typename DerivedF >
const Eigen::PlainObjectBase<DerivedF>& igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::F

◆ Fcut

template<typename DerivedV , typename DerivedF >
const Eigen::PlainObjectBase<DerivedF>& igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::Fcut

◆ Handle_Singular

template<typename DerivedV , typename DerivedF >
const Eigen::Matrix<int, Eigen::Dynamic, 1>& igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::Handle_Singular

◆ Handle_Stiffness

template<typename DerivedV , typename DerivedF >
Eigen::VectorXd igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::Handle_Stiffness

◆ Handle_SystemInfo

template<typename DerivedV , typename DerivedF >
const MeshSystemInfo& igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::Handle_SystemInfo

◆ Hard_constraints

template<typename DerivedV , typename DerivedF >
std::vector<int> igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::Hard_constraints

vector of blocked vertices

◆ ids_to_round

template<typename DerivedV , typename DerivedF >
std::vector<int> igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::ids_to_round

vector of indexes to round

◆ integer_rounding

template<typename DerivedV , typename DerivedF >
bool igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::integer_rounding

boolean that is true if rounding to integer is needed

◆ Lhs

template<typename DerivedV , typename DerivedF >
Eigen::SparseMatrix<double> igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::Lhs

◆ n_fixed_vars

template<typename DerivedV , typename DerivedF >
unsigned int igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::n_fixed_vars

number of fixed vertex

◆ n_integer_vars

template<typename DerivedV , typename DerivedF >
unsigned int igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::n_integer_vars

the total number of integer variables

◆ n_vert_vars

template<typename DerivedV , typename DerivedF >
unsigned int igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::n_vert_vars

the number of REAL variables for vertices

◆ num_constraint_equations

template<typename DerivedV , typename DerivedF >
unsigned int igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::num_constraint_equations

total number of constraints equations

◆ num_cut_constraint

template<typename DerivedV , typename DerivedF >
unsigned int igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::num_cut_constraint

CONSTRAINT PART number of cuts constraints.

◆ num_total_vars

template<typename DerivedV , typename DerivedF >
unsigned int igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::num_total_vars

total number of variables of the system, do not consider constraints, but consider integer vars

◆ num_userdefined_constraint

template<typename DerivedV , typename DerivedF >
unsigned int igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::num_userdefined_constraint

◆ PD1

template<typename DerivedV , typename DerivedF >
const Eigen::PlainObjectBase<DerivedV>& igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::PD1

◆ PD2

template<typename DerivedV , typename DerivedF >
const Eigen::PlainObjectBase<DerivedV>& igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::PD2

◆ rhs

template<typename DerivedV , typename DerivedF >
Eigen::VectorXd igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::rhs

◆ TT

template<typename DerivedV , typename DerivedF >
const Eigen::PlainObjectBase<DerivedF>& igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::TT

◆ TTi

template<typename DerivedV , typename DerivedF >
const Eigen::PlainObjectBase<DerivedF>& igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::TTi

◆ userdefined_constraints

template<typename DerivedV , typename DerivedF >
std::vector<std::vector<int > > igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::userdefined_constraints

vector of indexes to round

◆ UV

template<typename DerivedV , typename DerivedF >
Eigen::MatrixXd igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::UV

◆ UV_out

◆ V

template<typename DerivedV , typename DerivedF >
const Eigen::PlainObjectBase<DerivedV>& igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::V

◆ Vcut

template<typename DerivedV , typename DerivedF >
const Eigen::PlainObjectBase<DerivedV>& igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::Vcut

◆ VF

template<typename DerivedV , typename DerivedF >
std::vector<std::vector<int> > igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::VF

◆ VFi

template<typename DerivedV , typename DerivedF >
std::vector<std::vector<int> > igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::VFi

◆ WUV

◆ X

template<typename DerivedV , typename DerivedF >
std::vector< double > igl::copyleft::comiso::PoissonSolver< DerivedV, DerivedF >::X

vector of unknowns


The documentation for this class was generated from the following file: