TR-mbed 1.0
Loading...
Searching...
No Matches
CompressedStorage.h
Go to the documentation of this file.
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2008-2014 Gael Guennebaud <gael.guennebaud@inria.fr>
5//
6// This Source Code Form is subject to the terms of the Mozilla
7// Public License v. 2.0. If a copy of the MPL was not distributed
8// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9
10#ifndef EIGEN_COMPRESSED_STORAGE_H
11#define EIGEN_COMPRESSED_STORAGE_H
12
13namespace Eigen {
14
15namespace internal {
16
21template<typename _Scalar,typename _StorageIndex>
23{
24 public:
25
26 typedef _Scalar Scalar;
27 typedef _StorageIndex StorageIndex;
28
29 protected:
30
32
33 public:
34
38
41 {
42 resize(size);
43 }
44
47 {
48 *this = other;
49 }
50
52 {
53 resize(other.size());
54 if(other.size()>0)
55 {
58 }
59 return *this;
60 }
61
63 {
64 std::swap(m_values, other.m_values);
65 std::swap(m_indices, other.m_indices);
66 std::swap(m_size, other.m_size);
67 std::swap(m_allocatedSize, other.m_allocatedSize);
68 }
69
71 {
72 delete[] m_values;
73 delete[] m_indices;
74 }
75
82
83 void squeeze()
84 {
87 }
88
100
101 void append(const Scalar& v, Index i)
102 {
103 Index id = m_size;
104 resize(m_size+1, 1);
105 m_values[id] = v;
107 }
108
109 inline Index size() const { return m_size; }
110 inline Index allocatedSize() const { return m_allocatedSize; }
111 inline void clear() { m_size = 0; }
112
113 const Scalar* valuePtr() const { return m_values; }
114 Scalar* valuePtr() { return m_values; }
115 const StorageIndex* indexPtr() const { return m_indices; }
117
119 inline const Scalar& value(Index i) const { eigen_internal_assert(m_values!=0); return m_values[i]; }
120
122 inline const StorageIndex& index(Index i) const { eigen_internal_assert(m_indices!=0); return m_indices[i]; }
123
125 inline Index searchLowerIndex(Index key) const
126 {
127 return searchLowerIndex(0, m_size, key);
128 }
129
131 inline Index searchLowerIndex(Index start, Index end, Index key) const
132 {
133 while(end>start)
134 {
135 Index mid = (end+start)>>1;
136 if (m_indices[mid]<key)
137 start = mid+1;
138 else
139 end = mid;
140 }
141 return start;
142 }
143
146 inline Scalar at(Index key, const Scalar& defaultValue = Scalar(0)) const
147 {
148 if (m_size==0)
149 return defaultValue;
150 else if (key==m_indices[m_size-1])
151 return m_values[m_size-1];
152 // ^^ optimization: let's first check if it is the last coefficient
153 // (very common in high level algorithms)
154 const Index id = searchLowerIndex(0,m_size-1,key);
155 return ((id<m_size) && (m_indices[id]==key)) ? m_values[id] : defaultValue;
156 }
157
159 inline Scalar atInRange(Index start, Index end, Index key, const Scalar &defaultValue = Scalar(0)) const
160 {
161 if (start>=end)
162 return defaultValue;
163 else if (end>start && key==m_indices[end-1])
164 return m_values[end-1];
165 // ^^ optimization: let's first check if it is the last coefficient
166 // (very common in high level algorithms)
167 const Index id = searchLowerIndex(start,end-1,key);
168 return ((id<end) && (m_indices[id]==key)) ? m_values[id] : defaultValue;
169 }
170
175 {
176 Index id = searchLowerIndex(0,m_size,key);
177 if (id>=m_size || m_indices[id]!=key)
178 {
180 {
181 m_allocatedSize = 2*(m_size+1);
184
185 // copy first chunk
188
189 // copy the rest
190 if(m_size>id)
191 {
194 }
195 std::swap(m_values,newValues.ptr());
196 std::swap(m_indices,newIndices.ptr());
197 }
198 else if(m_size>id)
199 {
202 }
203 m_size++;
206 }
207 return m_values[id];
208 }
209
225
226 void prune(const Scalar& reference, const RealScalar& epsilon = NumTraits<RealScalar>::dummy_precision())
227 {
228 Index k = 0;
229 Index n = size();
230 for (Index i=0; i<n; ++i)
231 {
232 if (!internal::isMuchSmallerThan(value(i), reference, epsilon))
233 {
234 value(k) = value(i);
235 index(k) = index(i);
236 ++k;
237 }
238 }
239 resize(k,0);
240 }
241
242 protected:
243
244 inline void reallocate(Index size)
245 {
246 #ifdef EIGEN_SPARSE_COMPRESSED_STORAGE_REALLOCATE_PLUGIN
248 #endif
252 Index copySize = (std::min)(size, m_size);
253 if (copySize>0) {
256 }
257 std::swap(m_values,newValues.ptr());
258 std::swap(m_indices,newIndices.ptr());
260 }
261
262 protected:
267
268};
269
270} // end namespace internal
271
272} // end namespace Eigen
273
274#endif // EIGEN_COMPRESSED_STORAGE_H
Array< int, Dynamic, 1 > v
Definition Array_initializer_list_vector_cxx11.cpp:1
int n
Definition BiCGSTAB_simple.cpp:1
int i
Definition BiCGSTAB_step_by_step.cpp:9
#define eigen_internal_assert(x)
Definition Macros.h:1043
Definition CompressedStorage.h:23
const StorageIndex & index(Index i) const
Definition CompressedStorage.h:122
Scalar at(Index key, const Scalar &defaultValue=Scalar(0)) const
Definition CompressedStorage.h:146
~CompressedStorage()
Definition CompressedStorage.h:70
void reserve(Index size)
Definition CompressedStorage.h:76
Index searchLowerIndex(Index start, Index end, Index key) const
Definition CompressedStorage.h:131
CompressedStorage(Index size)
Definition CompressedStorage.h:39
Index size() const
Definition CompressedStorage.h:109
StorageIndex * m_indices
Definition CompressedStorage.h:264
void moveChunk(Index from, Index to, Index chunkSize)
Definition CompressedStorage.h:210
Scalar * valuePtr()
Definition CompressedStorage.h:114
CompressedStorage & operator=(const CompressedStorage &other)
Definition CompressedStorage.h:51
CompressedStorage(const CompressedStorage &other)
Definition CompressedStorage.h:45
Index allocatedSize() const
Definition CompressedStorage.h:110
_Scalar Scalar
Definition CompressedStorage.h:26
Scalar atInRange(Index start, Index end, Index key, const Scalar &defaultValue=Scalar(0)) const
Definition CompressedStorage.h:159
Scalar & value(Index i)
Definition CompressedStorage.h:118
const Scalar & value(Index i) const
Definition CompressedStorage.h:119
const Scalar * valuePtr() const
Definition CompressedStorage.h:113
StorageIndex * indexPtr()
Definition CompressedStorage.h:116
void clear()
Definition CompressedStorage.h:111
void append(const Scalar &v, Index i)
Definition CompressedStorage.h:101
Index m_size
Definition CompressedStorage.h:265
StorageIndex & index(Index i)
Definition CompressedStorage.h:121
Index m_allocatedSize
Definition CompressedStorage.h:266
Scalar & atWithInsertion(Index key, const Scalar &defaultValue=Scalar(0))
Definition CompressedStorage.h:174
void squeeze()
Definition CompressedStorage.h:83
void resize(Index size, double reserveSizeFactor=0)
Definition CompressedStorage.h:89
void swap(CompressedStorage &other)
Definition CompressedStorage.h:62
void prune(const Scalar &reference, const RealScalar &epsilon=NumTraits< RealScalar >::dummy_precision())
Definition CompressedStorage.h:226
Scalar * m_values
Definition CompressedStorage.h:263
NumTraits< Scalar >::Real RealScalar
Definition CompressedStorage.h:31
_StorageIndex StorageIndex
Definition CompressedStorage.h:27
const StorageIndex * indexPtr() const
Definition CompressedStorage.h:115
Index searchLowerIndex(Index key) const
Definition CompressedStorage.h:125
void reallocate(Index size)
Definition CompressedStorage.h:244
CompressedStorage()
Definition CompressedStorage.h:35
EIGEN_DEVICE_FUNC bool isMuchSmallerThan(const Scalar &x, const OtherScalar &y, const typename NumTraits< Scalar >::Real &precision=NumTraits< Scalar >::dummy_precision())
Definition MathFunctions.h:1940
EIGEN_DEVICE_FUNC void throw_std_bad_alloc()
Definition Memory.h:67
EIGEN_DEVICE_FUNC void smart_copy(const T *start, const T *end, T *target)
Definition Memory.h:515
void smart_memmove(const T *start, const T *end, T *target)
Definition Memory.h:539
Namespace containing all symbols from the Eigen library.
Definition bench_norm.cpp:85
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:74
Definition BandTriangularSolver.h:13
#define EIGEN_SPARSE_COMPRESSED_STORAGE_REALLOCATE_PLUGIN
Definition sparse_basic.cpp:14
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition NumTraits.h:233
Definition ForwardDeclarations.h:17