libigl v2.5.0
Loading...
Searching...
No Matches
min_quad_with_fixed.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) 2014 Alec Jacobson <alecjacobson@gmail.com>
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_MIN_QUAD_WITH_FIXED_H
9#define IGL_MIN_QUAD_WITH_FIXED_H
10#include "igl_inline.h"
11
12#define EIGEN_YES_I_KNOW_SPARSE_MODULE_IS_NOT_STABLE_YET
13#include <Eigen/Core>
14#include <Eigen/Dense>
15#include <Eigen/Sparse>
16// Bug in unsupported/Eigen/SparseExtra needs iostream first
17#include <iostream>
18#include <unsupported/Eigen/SparseExtra>
19
20namespace igl
21{
22 template <typename T>
23 struct min_quad_with_fixed_data;
49 template <typename T, typename Derivedknown>
51 const Eigen::SparseMatrix<T>& A,
52 const Eigen::MatrixBase<Derivedknown> & known,
53 const Eigen::SparseMatrix<T>& Aeq,
54 const bool pd,
55 min_quad_with_fixed_data<T> & data
56 );
69 template <
70 typename T,
71 typename DerivedB,
72 typename DerivedY,
73 typename DerivedBeq,
74 typename DerivedZ,
75 typename Derivedsol>
77 const min_quad_with_fixed_data<T> & data,
78 const Eigen::MatrixBase<DerivedB> & B,
79 const Eigen::MatrixBase<DerivedY> & Y,
80 const Eigen::MatrixBase<DerivedBeq> & Beq,
81 Eigen::PlainObjectBase<DerivedZ> & Z,
82 Eigen::PlainObjectBase<Derivedsol> & sol);
84 template <
85 typename T,
86 typename DerivedB,
87 typename DerivedY,
88 typename DerivedBeq,
89 typename DerivedZ>
91 const min_quad_with_fixed_data<T> & data,
92 const Eigen::MatrixBase<DerivedB> & B,
93 const Eigen::MatrixBase<DerivedY> & Y,
94 const Eigen::MatrixBase<DerivedBeq> & Beq,
95 Eigen::PlainObjectBase<DerivedZ> & Z);
99 template <
100 typename T,
101 typename Derivedknown,
102 typename DerivedB,
103 typename DerivedY,
104 typename DerivedBeq,
105 typename DerivedZ>
107 const Eigen::SparseMatrix<T>& A,
108 const Eigen::MatrixBase<DerivedB> & B,
109 const Eigen::MatrixBase<Derivedknown> & known,
110 const Eigen::MatrixBase<DerivedY> & Y,
111 const Eigen::SparseMatrix<T>& Aeq,
112 const Eigen::MatrixBase<DerivedBeq> & Beq,
113 const bool pd,
114 Eigen::PlainObjectBase<DerivedZ> & Z);
135 template <typename Scalar, int n, int m, bool Hpd=true>
136 IGL_INLINE Eigen::Matrix<Scalar,n,1> min_quad_with_fixed(
137 const Eigen::Matrix<Scalar,n,n> & H,
138 const Eigen::Matrix<Scalar,n,1> & f,
139 const Eigen::Array<bool,n,1> & k,
140 const Eigen::Matrix<Scalar,n,1> & bc,
141 const Eigen::Matrix<Scalar,m,n> & A,
142 const Eigen::Matrix<Scalar,m,1> & b);
144 template <typename Scalar, int n, bool Hpd=true>
145 IGL_INLINE Eigen::Matrix<Scalar,n,1> min_quad_with_fixed(
146 const Eigen::Matrix<Scalar,n,n> & H,
147 const Eigen::Matrix<Scalar,n,1> & f,
148 const Eigen::Array<bool,n,1> & k,
149 const Eigen::Matrix<Scalar,n,1> & bc);
155 template <typename Scalar, int n, int kcount, bool Hpd/*no default*/>
156 IGL_INLINE Eigen::Matrix<Scalar,n,1> min_quad_with_fixed(
157 const Eigen::Matrix<Scalar,n,n> & H,
158 const Eigen::Matrix<Scalar,n,1> & f,
159 const Eigen::Array<bool,n,1> & k,
160 const Eigen::Matrix<Scalar,n,1> & bc);
161}
162
164template <typename T>
166{
168 int n;
170 bool Auu_pd;
174 Eigen::VectorXi known;
176 Eigen::VectorXi unknown;
178 Eigen::VectorXi lagrange;
180 Eigen::VectorXi unknown_lagrange;
182 Eigen::SparseMatrix<T> preY;
193 Eigen::SimplicialLLT <Eigen::SparseMatrix<T > > llt;
194 Eigen::SimplicialLDLT<Eigen::SparseMatrix<T > > ldlt;
195 Eigen::SparseLU<Eigen::SparseMatrix<T, Eigen::ColMajor>, Eigen::COLAMDOrdering<int> > lu;
198 bool Aeq_li;
200 int neq;
201 Eigen::SparseQR<Eigen::SparseMatrix<T>, Eigen::COLAMDOrdering<int> > AeqTQR;
202 Eigen::SparseMatrix<T> Aeqk;
203 Eigen::SparseMatrix<T> Aequ;
204 Eigen::SparseMatrix<T> Auu;
205 Eigen::SparseMatrix<T> AeqTQ1;
206 Eigen::SparseMatrix<T> AeqTQ1T;
207 Eigen::SparseMatrix<T> AeqTQ2;
208 Eigen::SparseMatrix<T> AeqTQ2T;
209 Eigen::SparseMatrix<T> AeqTR1;
210 Eigen::SparseMatrix<T> AeqTR1T;
211 Eigen::SparseMatrix<T> AeqTE;
212 Eigen::SparseMatrix<T> AeqTET;
214 Eigen::SparseMatrix<T> NA;
215 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> NB;
216};
217
218#ifndef IGL_STATIC_LIBRARY
220#endif
221
222#endif
#define IGL_INLINE
Definition igl_inline.h:15
Definition AABB.h:17
bool min_quad_with_fixed(const Eigen::SparseMatrix< T > &A, const Eigen::MatrixBase< DerivedB > &B, const Eigen::MatrixBase< Derivedknown > &known, const Eigen::MatrixBase< DerivedY > &Y, const Eigen::SparseMatrix< T > &Aeq, const Eigen::MatrixBase< DerivedBeq > &Beq, const bool pd, Eigen::PlainObjectBase< DerivedZ > &Z)
This is an overloaded member function, provided for convenience. It differs from the above function o...
Definition min_quad_with_fixed.impl.h:567
bool min_quad_with_fixed_solve(const min_quad_with_fixed_data< T > &data, const Eigen::MatrixBase< DerivedB > &B, const Eigen::MatrixBase< DerivedY > &Y, const Eigen::MatrixBase< DerivedBeq > &Beq, Eigen::PlainObjectBase< DerivedZ > &Z, Eigen::PlainObjectBase< Derivedsol > &sol)
Solves a system previously factored using min_quad_with_fixed_precompute.
Definition min_quad_with_fixed.impl.h:394
bool min_quad_with_fixed_precompute(const Eigen::SparseMatrix< T > &A, const Eigen::MatrixBase< Derivedknown > &known, const Eigen::SparseMatrix< T > &Aeq, const bool pd, min_quad_with_fixed_data< T > &data)
Minimize a convex quadratic energy subject to fixed value and linear equality constraints.
Definition min_quad_with_fixed.impl.h:30
Parameters and precomputed values for min_quad_with_fixed.
Definition min_quad_with_fixed.h:166
Eigen::VectorXi known
Indices of known variables.
Definition min_quad_with_fixed.h:174
Eigen::SparseQR< Eigen::SparseMatrix< T >, Eigen::COLAMDOrdering< int > > AeqTQR
Definition min_quad_with_fixed.h:201
SolverType
Type of solver used.
Definition min_quad_with_fixed.h:185
@ LU
Definition min_quad_with_fixed.h:188
@ LDLT
Definition min_quad_with_fixed.h:187
@ QR_LLT
Definition min_quad_with_fixed.h:189
@ NUM_SOLVER_TYPES
Definition min_quad_with_fixed.h:190
@ LLT
Definition min_quad_with_fixed.h:186
Eigen::SimplicialLDLT< Eigen::SparseMatrix< T > > ldlt
Definition min_quad_with_fixed.h:194
Eigen::SparseMatrix< T > preY
Matrix multiplied against Y when constructing right hand side.
Definition min_quad_with_fixed.h:182
Eigen::VectorXi unknown
Indices of unknown variables.
Definition min_quad_with_fixed.h:176
Eigen::SimplicialLLT< Eigen::SparseMatrix< T > > llt
Solver data (factorization)
Definition min_quad_with_fixed.h:193
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > NB
Definition min_quad_with_fixed.h:215
Eigen::SparseMatrix< T > AeqTQ2
Definition min_quad_with_fixed.h:207
Eigen::SparseMatrix< T > AeqTQ1
Definition min_quad_with_fixed.h:205
Eigen::SparseLU< Eigen::SparseMatrix< T, Eigen::ColMajor >, Eigen::COLAMDOrdering< int > > lu
Definition min_quad_with_fixed.h:195
int n
Size of original system: number of unknowns + number of knowns.
Definition min_quad_with_fixed.h:168
bool Aeq_li
QR factorization Are rows of Aeq linearly independent?
Definition min_quad_with_fixed.h:198
Eigen::VectorXi lagrange
Indices of lagrange variables.
Definition min_quad_with_fixed.h:178
bool Auu_sym
Whether A(unknown,unknown) is symmetric.
Definition min_quad_with_fixed.h:172
int neq
Columns of Aeq corresponding to unknowns.
Definition min_quad_with_fixed.h:200
Eigen::SparseMatrix< T > AeqTR1
Definition min_quad_with_fixed.h:209
Eigen::VectorXi unknown_lagrange
Indices of unknown variable followed by Indices of lagrange variables.
Definition min_quad_with_fixed.h:180
Eigen::SparseMatrix< T > AeqTET
Definition min_quad_with_fixed.h:212
bool Auu_pd
Whether A(unknown,unknown) is positive definite.
Definition min_quad_with_fixed.h:170
Eigen::SparseMatrix< T > AeqTR1T
Definition min_quad_with_fixed.h:210
Eigen::SparseMatrix< T > AeqTQ2T
Definition min_quad_with_fixed.h:208
Eigen::SparseMatrix< T > Aeqk
Definition min_quad_with_fixed.h:202
Eigen::SparseMatrix< T > Auu
Definition min_quad_with_fixed.h:204
enum igl::min_quad_with_fixed_data::SolverType solver_type
Eigen::SparseMatrix< T > Aequ
Definition min_quad_with_fixed.h:203
Eigen::SparseMatrix< T > AeqTE
Definition min_quad_with_fixed.h:211
Eigen::SparseMatrix< T > AeqTQ1T
Definition min_quad_with_fixed.h:206