QPALM
1.2.3
Proximal Augmented Lagrangian method for Quadratic Programs
|
#include <Eigen/Sparse>
#include <Eigen/src/Core/Map.h>
#include <Eigen/src/Core/Matrix.h>
#include <ladel_constants.h>
#include <ladel_types.h>
#include <qpalm/global_opts.h>
#include <memory>
#include <qpalm/qpalm_cxx-export.hpp>
Go to the source code of this file.
Data Structures | |
struct | qpalm::alloc::ladel_sparse_matrix_deleter |
Namespaces | |
namespace | qpalm |
namespace | qpalm::alloc |
RAII-based wrappers for the allocation and deallocation functions of the C API. | |
Typedefs | |
using | qpalm::index_t = Eigen::Index |
Index types for vectors and matrices. | |
using | qpalm::sp_index_t = ladel_int |
Index type for sparse matrices representation. | |
using | qpalm::sparse_mat_t = Eigen::SparseMatrix< c_float, Eigen::ColMajor, sp_index_t > |
Owning sparse matrix type. | |
using | qpalm::sparse_mat_view_t = Eigen::Map< const sparse_mat_t > |
Read-only view on a sparse matrix. | |
using | qpalm::sparse_mat_ref_t = Eigen::Ref< const sparse_mat_t > |
Read-only reference to a sparse matrix. | |
using | qpalm::triplet_t = Eigen::Triplet< c_float, sp_index_t > |
Type for (row, column, value) triplets for initializing sparse matrices. | |
using | qpalm::vec_t = Eigen::Matrix< c_float, Eigen::Dynamic, 1 > |
Owning dense vector type. | |
using | qpalm::borrowed_vec_t = Eigen::Map< vec_t > |
Borrowed dense vector type (vector view). | |
using | qpalm::const_borrowed_vec_t = Eigen::Map< const vec_t > |
Read-only borrowed dense vector type (vector view). | |
using | qpalm::ref_vec_t = Eigen::Ref< vec_t > |
Reference to a dense vector (vector view). | |
using | qpalm::const_ref_vec_t = Eigen::Ref< const vec_t > |
Read-only reference to a dense vector (vector view). | |
using | qpalm::ladel_sparse_matrix_ptr = std::unique_ptr< ladel_sparse_matrix, alloc::ladel_sparse_matrix_deleter > |
Smart pointer that automatically cleans up an owning ladel_sparse_matrix object. | |
Functions | |
ladel_sparse_matrix | qpalm::eigen_to_ladel (sparse_mat_t &mat, ladel_int symmetry=UNSYMMETRIC) |
Convert an Eigen sparse matrix to a LADEL sparse matrix, without creating a copy. | |
ladel_sparse_matrix_ptr | qpalm::ladel_sparse_create (index_t rows, index_t cols, index_t nnz, ladel_int symmetry, bool values=true, bool nonzeros=false) |
Create an LADEL sparse matrix of the given dimensions. | |
ladel_sparse_matrix_ptr | qpalm::eigen_to_ladel_copy (const sparse_mat_ref_t &mat, ladel_int symmetry=UNSYMMETRIC) |
Similar to eigen_to_ladel, but creates a copy of all data, in such a way that the returned matrix is completely decoupled from mat , and such that it can be reallocated and deallocated by the ladel_sparse_free and similar functions. | |