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 
13 namespace Eigen {
14 
15 namespace internal {
16 
21 template<typename _Scalar,typename _StorageIndex>
23 {
24  public:
25 
26  typedef _Scalar Scalar;
27  typedef _StorageIndex StorageIndex;
28 
29  protected:
30 
32 
33  public:
34 
36  : m_values(0), m_indices(0), m_size(0), m_allocatedSize(0)
37  {}
38 
40  : m_values(0), m_indices(0), m_size(0), m_allocatedSize(0)
41  {
42  resize(size);
43  }
44 
46  : m_values(0), m_indices(0), m_size(0), m_allocatedSize(0)
47  {
48  *this = other;
49  }
50 
52  {
53  resize(other.size());
54  if(other.size()>0)
55  {
58  }
59  return *this;
60  }
61 
62  void swap(CompressedStorage& other)
63  {
64  std::swap(m_values, other.m_values);
66  std::swap(m_size, other.m_size);
68  }
69 
71  {
72  delete[] m_values;
73  delete[] m_indices;
74  }
75 
77  {
78  Index newAllocatedSize = m_size + size;
79  if (newAllocatedSize > m_allocatedSize)
80  reallocate(newAllocatedSize);
81  }
82 
83  void squeeze()
84  {
87  }
88 
89  void resize(Index size, double reserveSizeFactor = 0)
90  {
91  if (m_allocatedSize<size)
92  {
93  Index realloc_size = (std::min<Index>)(NumTraits<StorageIndex>::highest(), size + Index(reserveSizeFactor*double(size)));
94  if(realloc_size<size)
96  reallocate(realloc_size);
97  }
98  m_size = size;
99  }
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;
106  m_indices[id] = internal::convert_index<StorageIndex>(i);
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; }
116  StorageIndex* indexPtr() { return m_indices; }
117 
118  inline Scalar& value(Index i) { eigen_internal_assert(m_values!=0); return m_values[i]; }
119  inline const Scalar& value(Index i) const { eigen_internal_assert(m_values!=0); return m_values[i]; }
120 
121  inline StorageIndex& index(Index i) { eigen_internal_assert(m_indices!=0); return m_indices[i]; }
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 
174  inline Scalar& atWithInsertion(Index key, const Scalar& defaultValue = Scalar(0))
175  {
176  Index id = searchLowerIndex(0,m_size,key);
177  if (id>=m_size || m_indices[id]!=key)
178  {
179  if (m_allocatedSize<m_size+1)
180  {
181  m_allocatedSize = 2*(m_size+1);
184 
185  // copy first chunk
186  internal::smart_copy(m_values, m_values +id, newValues.ptr());
187  internal::smart_copy(m_indices, m_indices+id, newIndices.ptr());
188 
189  // copy the rest
190  if(m_size>id)
191  {
192  internal::smart_copy(m_values +id, m_values +m_size, newValues.ptr() +id+1);
193  internal::smart_copy(m_indices+id, m_indices+m_size, newIndices.ptr()+id+1);
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++;
204  m_indices[id] = internal::convert_index<StorageIndex>(key);
205  m_values[id] = defaultValue;
206  }
207  return m_values[id];
208  }
209 
210  void prune(const Scalar& reference, const RealScalar& epsilon = NumTraits<RealScalar>::dummy_precision())
211  {
212  Index k = 0;
213  Index n = size();
214  for (Index i=0; i<n; ++i)
215  {
216  if (!internal::isMuchSmallerThan(value(i), reference, epsilon))
217  {
218  value(k) = value(i);
219  index(k) = index(i);
220  ++k;
221  }
222  }
223  resize(k,0);
224  }
225 
226  protected:
227 
228  inline void reallocate(Index size)
229  {
230  #ifdef EIGEN_SPARSE_COMPRESSED_STORAGE_REALLOCATE_PLUGIN
231  EIGEN_SPARSE_COMPRESSED_STORAGE_REALLOCATE_PLUGIN
232  #endif
234  internal::scoped_array<Scalar> newValues(size);
235  internal::scoped_array<StorageIndex> newIndices(size);
236  Index copySize = (std::min)(size, m_size);
237  if (copySize>0) {
238  internal::smart_copy(m_values, m_values+copySize, newValues.ptr());
239  internal::smart_copy(m_indices, m_indices+copySize, newIndices.ptr());
240  }
241  std::swap(m_values,newValues.ptr());
242  std::swap(m_indices,newIndices.ptr());
244  }
245 
246  protected:
247  Scalar* m_values;
248  StorageIndex* m_indices;
251 
252 };
253 
254 } // end namespace internal
255 
256 } // end namespace Eigen
257 
258 #endif // EIGEN_COMPRESSED_STORAGE_H
double epsilon
const Scalar & value(Index i) const
CompressedStorage & operator=(const CompressedStorage &other)
Definition: LDLT.h:16
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:150
const StorageIndex & index(Index i) const
EIGEN_DEVICE_FUNC void smart_copy(const T *start, const T *end, T *target)
Definition: Memory.h:485
void append(const Scalar &v, Index i)
Index searchLowerIndex(Index key) const
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
CompressedStorage(const CompressedStorage &other)
void resize(Index size, double reserveSizeFactor=0)
Scalar & atWithInsertion(Index key, const Scalar &defaultValue=Scalar(0))
void prune(const Scalar &reference, const RealScalar &epsilon=NumTraits< RealScalar >::dummy_precision())
Scalar atInRange(Index start, Index end, Index key, const Scalar &defaultValue=Scalar(0)) const
const StorageIndex * indexPtr() const
NumTraits< Scalar >::Real RealScalar
EIGEN_DEVICE_FUNC void throw_std_bad_alloc()
Definition: Memory.h:67
#define eigen_internal_assert(x)
Definition: Macros.h:583
Index searchLowerIndex(Index start, Index end, Index key) const
void smart_memmove(const T *start, const T *end, T *target)
Definition: Memory.h:508
void swap(CompressedStorage &other)
void swap(mpfr::mpreal &x, mpfr::mpreal &y)
Definition: mpreal.h:2986
Scalar at(Index key, const Scalar &defaultValue=Scalar(0)) const


hebiros
Author(s): Xavier Artache , Matthew Tesch
autogenerated on Thu Sep 3 2020 04:08:04