10#ifndef EIGEN_COMPRESSED_STORAGE_H
11#define EIGEN_COMPRESSED_STORAGE_H
14#include "./InternalHeaderCheck.h"
24template <
typename Scalar_,
typename StorageIndex_>
25class CompressedStorage {
27 typedef Scalar_ Scalar;
28 typedef StorageIndex_ StorageIndex;
31 typedef typename NumTraits<Scalar>::Real RealScalar;
34 CompressedStorage() : m_values(0), m_indices(0), m_size(0), m_allocatedSize(0) {}
36 explicit CompressedStorage(Index size) : m_values(0), m_indices(0), m_size(0), m_allocatedSize(0) { resize(size); }
38 CompressedStorage(
const CompressedStorage& other) : m_values(0), m_indices(0), m_size(0), m_allocatedSize(0) {
42 CompressedStorage& operator=(
const CompressedStorage& other) {
44 if (other.size() > 0) {
45 internal::smart_copy(other.m_values, other.m_values + m_size, m_values);
46 internal::smart_copy(other.m_indices, other.m_indices + m_size, m_indices);
51 void swap(CompressedStorage& other) {
52 std::swap(m_values, other.m_values);
53 std::swap(m_indices, other.m_indices);
54 std::swap(m_size, other.m_size);
55 std::swap(m_allocatedSize, other.m_allocatedSize);
58 ~CompressedStorage() {
59 conditional_aligned_delete_auto<Scalar, true>(m_values, m_allocatedSize);
60 conditional_aligned_delete_auto<StorageIndex, true>(m_indices, m_allocatedSize);
63 void reserve(Index size) {
64 Index newAllocatedSize = m_size + size;
65 if (newAllocatedSize > m_allocatedSize) reallocate(newAllocatedSize);
69 if (m_allocatedSize > m_size) reallocate(m_size);
72 void resize(Index size,
double reserveSizeFactor = 0) {
73 if (m_allocatedSize < size) {
75 using SmallerIndexType =
76 typename std::conditional<static_cast<size_t>((std::numeric_limits<Index>::max)()) <
77 static_cast<size_t>((std::numeric_limits<StorageIndex>::max)()),
78 Index, StorageIndex>::type;
80 (std::min<Index>)(NumTraits<SmallerIndexType>::highest(), size +
Index(reserveSizeFactor *
double(size)));
81 if (realloc_size < size) internal::throw_std_bad_alloc();
82 reallocate(realloc_size);
87 void append(
const Scalar& v, Index i) {
89 resize(m_size + 1, 1);
91 m_indices[id] = internal::convert_index<StorageIndex>(i);
94 inline Index size()
const {
return m_size; }
95 inline Index allocatedSize()
const {
return m_allocatedSize; }
96 inline void clear() { m_size = 0; }
98 const Scalar* valuePtr()
const {
return m_values; }
99 Scalar* valuePtr() {
return m_values; }
100 const StorageIndex* indexPtr()
const {
return m_indices; }
101 StorageIndex* indexPtr() {
return m_indices; }
103 inline Scalar& value(Index i) {
104 eigen_internal_assert(m_values != 0);
107 inline const Scalar& value(Index i)
const {
108 eigen_internal_assert(m_values != 0);
112 inline StorageIndex& index(Index i) {
113 eigen_internal_assert(m_indices != 0);
116 inline const StorageIndex& index(Index i)
const {
117 eigen_internal_assert(m_indices != 0);
122 inline Index searchLowerIndex(Index key)
const {
return searchLowerIndex(0, m_size, key); }
125 inline Index searchLowerIndex(Index start, Index end, Index key)
const {
126 return static_cast<Index
>(std::distance(m_indices, std::lower_bound(m_indices + start, m_indices + end, key)));
131 inline Scalar at(Index key,
const Scalar& defaultValue = Scalar(0))
const {
134 else if (key == m_indices[m_size - 1])
135 return m_values[m_size - 1];
138 const Index
id = searchLowerIndex(0, m_size - 1, key);
139 return ((
id < m_size) && (m_indices[
id] == key)) ? m_values[id] : defaultValue;
143 inline Scalar atInRange(Index start, Index end, Index key,
const Scalar& defaultValue = Scalar(0))
const {
146 else if (end > start && key == m_indices[end - 1])
147 return m_values[end - 1];
150 const Index
id = searchLowerIndex(start, end - 1, key);
151 return ((
id < end) && (m_indices[
id] == key)) ? m_values[id] : defaultValue;
157 inline Scalar& atWithInsertion(Index key,
const Scalar& defaultValue = Scalar(0)) {
158 Index
id = searchLowerIndex(0, m_size, key);
159 if (
id >= m_size || m_indices[
id] != key) {
160 if (m_allocatedSize < m_size + 1) {
161 Index newAllocatedSize = 2 * (m_size + 1);
162 m_values = conditional_aligned_realloc_new_auto<Scalar, true>(m_values, newAllocatedSize, m_allocatedSize);
164 conditional_aligned_realloc_new_auto<StorageIndex, true>(m_indices, newAllocatedSize, m_allocatedSize);
165 m_allocatedSize = newAllocatedSize;
168 internal::smart_memmove(m_values +
id, m_values + m_size, m_values +
id + 1);
169 internal::smart_memmove(m_indices +
id, m_indices + m_size, m_indices +
id + 1);
172 m_indices[id] = internal::convert_index<StorageIndex>(key);
173 m_values[id] = defaultValue;
178 inline void moveChunk(Index from, Index to, Index chunkSize) {
179 eigen_internal_assert(chunkSize >= 0 && to + chunkSize <= m_size);
180 internal::smart_memmove(m_values + from, m_values + from + chunkSize, m_values + to);
181 internal::smart_memmove(m_indices + from, m_indices + from + chunkSize, m_indices + to);
185 inline void reallocate(Index size) {
186#ifdef EIGEN_SPARSE_COMPRESSED_STORAGE_REALLOCATE_PLUGIN
187 EIGEN_SPARSE_COMPRESSED_STORAGE_REALLOCATE_PLUGIN
189 eigen_internal_assert(size != m_allocatedSize);
190 m_values = conditional_aligned_realloc_new_auto<Scalar, true>(m_values, size, m_allocatedSize);
191 m_indices = conditional_aligned_realloc_new_auto<StorageIndex, true>(m_indices, size, m_allocatedSize);
192 m_allocatedSize = size;
197 StorageIndex* m_indices;
199 Index m_allocatedSize;
Namespace containing all symbols from the Eigen library.
Definition Core:137
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:83