OpenFPM_pdata  4.1.0
Project that contain the implementation of distributed structures
 
Loading...
Searching...
No Matches
Vector_eigen.hpp
1/*
2 * Vector_eigen.hpp
3 *
4 * Created on: Nov 27, 2015
5 * Author: i-bird
6 */
7
8#ifndef OPENFPM_NUMERICS_SRC_VECTOR_VECTOR_EIGEN_HPP_
9#define OPENFPM_NUMERICS_SRC_VECTOR_VECTOR_EIGEN_HPP_
10
11#include <type_traits>
12#include "util/mul_array_extents.hpp"
13#include <fstream>
14#include "Grid/staggered_dist_grid.hpp"
15#include "Space/Ghost.hpp"
16#include "FiniteDifference/util/common.hpp"
17#include <boost/mpl/vector_c.hpp>
18#include <unordered_map>
19#include "Vector_util.hpp"
20
21#define EIGEN_RVAL 1
22
29template<typename T>
30class rval<T,EIGEN_RVAL>
31{
33 long int r;
34
36 T val;
37
38public:
39
45 long int & row()
46 {
47 return r;
48 }
49
55 T & value()
56 {
57 return val;
58 }
59
64 :r(0)
65 {}
66
73 rval(long int i, T val)
74 {
75 row() = i;
76 value() = val;
77 }
78
84 static inline bool noPointers()
85 {
86 return true;
87 }
88};
89
90template<typename T>
91class Vector<T, EIGEN_BASE>
92{
94 mutable Eigen::Matrix<T, Eigen::Dynamic, 1> v;
95
98
101
103 mutable std::unordered_map<size_t,size_t> map;
104
107
110
111 //size of each chunk
112 mutable openfpm::vector<size_t> sz;
113
114
118 void collect() const
119 {
120 Vcluster<> & vcl = create_vcluster();
121
122 row_val_recv.clear();
123
124 // here we collect all the triplet in one array on the root node
125 vcl.SGather(row_val,row_val_recv,prc,sz,0);
126
127 if (vcl.getProcessUnitID() != 0)
128 row_val.resize(0);
129 else
130 row_val.swap(row_val_recv);
131
132 build_map();
133 }
134
139 void setEigen() const
140 {
141 // set the vector
142
143 for (size_t i = 0 ; i < row_val.size() ; i++)
144 v[row_val.get(i).row()] = row_val.get(i).value();
145 }
146
151 void build_map() const
152 {
153 map.clear();
154
155 for (size_t i = 0 ; i < row_val.size() ; i++)
156 map[row_val.get(i).row()] = i;
157 }
158
159public:
160
166 Vector(const Vector<T> & v)
167 :invalid(0)
168 {
169 this->operator=(v);
170 }
171
177 Vector(const Vector<T> && v)
178 :invalid(0)
179 {
180 this->operator=(v);
181 }
182
188 Vector(size_t n)
189 {
190 resize(n,0);
191 }
192
197 {
198 }
199
206 void resize(size_t row, size_t l_row)
207 {
208 v.resize(row);
209 }
210
217 void insert(size_t i, T val)
218 {
219 row_val.add();
220
221 // Map
222 map[i] = row_val.size()-1;
223
224 row_val.last().row() = i;
225 row_val.last().value() = val;
226 }
227
235 inline T & insert(size_t i)
236 {
237 row_val.add();
238
239 // Map
240 map[i] = row_val.size()-1;
241
242 row_val.last().row() = i;
243 return row_val.last().value();
244 }
245
253 inline const T & insert(size_t i) const
254 {
255 row_val.add();
256
257 // Map
258 map[i] = row_val.size()-1;
259
260 row_val.last().row() = i;
261 return row_val.last().value();
262 }
263
273 const T & operator()(size_t i) const
274 {
275 // Search if exist
276
277 std::unordered_map<size_t,size_t>::iterator it = map.find(i);
278
279 if ( it != map.end() )
280 return row_val.get(it->second).value();
281
282 return insert(i);
283 }
284
294 T & operator()(size_t i)
295 {
296 // Search if exist
297
298 std::unordered_map<size_t,size_t>::iterator it = map.find(i);
299
300 if ( it != map.end() )
301 return row_val.get(it->second).value();
302
303 return insert(i);
304 }
305
311 const Eigen::Matrix<T, Eigen::Dynamic, 1> & getVec() const
312 {
313 collect();
314 setEigen();
315
316 return v;
317 }
318
324 Eigen::Matrix<T, Eigen::Dynamic, 1> & getVec()
325 {
326 collect();
327 setEigen();
328
329 return v;
330 }
331
339 void scatter()
340 {
341 row_val_recv.clear();
342 Vcluster<> & vcl = create_vcluster();
343
344 vcl.SScatter(row_val,row_val_recv,prc,sz,0);
345
346 // if we do not receive anything a previous collect has not been performed
347 // and so nothing is scattered
348 if (row_val_recv.size() != 0)
349 {
350 row_val.clear();
351 row_val.add(row_val_recv);
352 build_map();
353 }
354 }
355
361 void fromFile(std::string file)
362 {
363 std::ifstream inputf;
364 inputf.open(file);
365
366 for (size_t i = 0 ; i < v.size() ; i++)
367 inputf >> v(i);
368
369 inputf.close();
370
371 }
372
379 {
380 prc = v.prc;
381 sz = v.sz;
382 map = v.map;
383 row_val = v.row_val;
384
385 return *this;
386 }
387
394 {
395 prc = v.prc;
396 sz = v.sz;
397 map = v.map;
398 row_val = v.row_val;
399
400 return *this;
401 }
402
410 Vector<T> & operator=(Eigen::Matrix<T, Eigen::Dynamic, 1> & v)
411 {
412 for (size_t i = 0 ; i < row_val.size() ; i++)
413 row_val.get(i).value() = v(row_val.get(i).row());
414
415 return *this;
416 }
417};
418
419
420#endif /* OPENFPM_NUMERICS_SRC_VECTOR_VECTOR_EIGEN_HPP_ */
size_t getProcessUnitID()
Get the process unit id.
Implementation of VCluster class.
Definition VCluster.hpp:59
bool SScatter(T &send, S &recv, openfpm::vector< size_t > &prc, openfpm::vector< size_t > &sz, size_t root)
Semantic Scatter, scatter the data from one processor to the other node.
Definition VCluster.hpp:621
bool SGather(T &send, S &recv, size_t root)
Semantic Gather, gather the data from all processors into one node.
Definition VCluster.hpp:450
Eigen::Matrix< T, Eigen::Dynamic, 1 > & getVec()
Get the Eigen Vector object.
void insert(size_t i, T val)
Return a reference to the vector element.
Vector< T > & operator=(Eigen::Matrix< T, Eigen::Dynamic, 1 > &v)
Copy the vector (it is used for special purpose)
std::unordered_map< size_t, size_t > map
global to local map
void scatter()
Scatter the vector information to the other processors.
const T & insert(size_t i) const
Return a reference to the vector element.
Vector(size_t n)
Create a vector with n elements.
openfpm::vector< rval< T, EIGEN_RVAL > > row_val
row val vector
void build_map() const
Build the map.
T & insert(size_t i)
Return a reference to the vector element.
openfpm::vector< size_t > prc
Processors from where we gather.
Vector< T > & operator=(const Vector< T > &v)
Copy the vector.
Vector< T > & operator=(const Vector< T > &&v)
Copy the vector.
void resize(size_t row, size_t l_row)
Resize the Vector.
Vector(const Vector< T > &v)
Copy the vector.
void collect() const
Here we collect the full vector on master.
Eigen::Matrix< T, Eigen::Dynamic, 1 > v
Eigen vector.
T & operator()(size_t i)
Return a reference to the vector element.
Vector()
Create a vector with 0 elements.
Vector(const Vector< T > &&v)
Copy the vector.
void setEigen() const
Set the Eigen internal vector.
void fromFile(std::string file)
Load from file.
openfpm::vector< rval< T, EIGEN_RVAL > > row_val_recv
row val vector received
const T & operator()(size_t i) const
Return a reference to the vector element.
const Eigen::Matrix< T, Eigen::Dynamic, 1 > & getVec() const
Get the Eigen Vector object.
Sparse Matrix implementation stub object when OpenFPM is compiled with no linear algebra support.
Definition Vector.hpp:40
void resize(size_t row, size_t row_n)
stub resize
Definition Vector.hpp:97
void insert(size_t i, T val)
stub insert
Definition Vector.hpp:105
Vector< T > & operator=(const Vector< T > &v)
stub operator=
Definition Vector.hpp:161
Implementation of 1-D std::vector like structure.
size_t size()
Stub size.
static bool noPointers()
Indicate that the structure has no pointer.
rval(long int i, T val)
Constructor from row, colum and value.
T & value()
Return the value.
rval()
Default constructor.
long int & row()
Return the row index.
It store one row value of a vector.
Definition Vector.hpp:22