Vector_eigen.hpp 6.69 KB
Newer Older
1 2 3 4 5 6 7 8 9 10
/*
 * Vector_eigen.hpp
 *
 *  Created on: Nov 27, 2015
 *      Author: i-bird
 */

#ifndef OPENFPM_NUMERICS_SRC_VECTOR_VECTOR_EIGEN_HPP_
#define OPENFPM_NUMERICS_SRC_VECTOR_VECTOR_EIGEN_HPP_

11 12 13
#include <type_traits>
#include "util/mul_array_extents.hpp"
#include <fstream>
14
#include "Grid/staggered_dist_grid.hpp"
15 16
#include "Space/Ghost.hpp"
#include "FiniteDifference/util/common.hpp"
17
#include <boost/mpl/vector_c.hpp>
18
#include <unordered_map>
19
#include "Vector_util.hpp"
20

21 22 23 24 25 26 27 28 29 30 31
#define EIGEN_RVAL 1

/*! \brief It store one row value of a vector
 *
 * Given a row, store a value
 *
 *
 */
template<typename T>
class rval<T,EIGEN_RVAL>
{
32
	//! row
33 34
	long int r;

35
	//! value
36 37 38 39
	T val;

public:

40 41 42 43 44
	/*! \brief Return the row index
	 *
	 * \return a reference to the row index
	 *
	 */
45 46 47 48 49
	long int & row()
	{
		return r;
	}

50 51 52 53 54
	/*! \brief Return the value
	 *
	 * \return a reference to the row value
	 *
	 */
55 56 57 58 59 60 61 62
	T & value()
	{
		return val;
	}

	/*! \brief Default constructor
	 *
	 */
63 64 65
	rval()
	:r(0)
	{}
66 67 68 69 70 71 72 73 74 75 76 77

	/*! \brief Constructor from row, colum and value
	 *
	 * \param i row
	 * \param val value
	 *
	 */
	rval(long int i, T val)
	{
		row() = i;
		value() = val;
	}
incardon's avatar
incardon committed
78 79 80 81 82 83 84 85 86 87

	/*! \brief Indicate that the structure has no pointer
	 *
	 * \return true
	 *
	 */
	static inline bool noPointers()
	{
		return true;
	}
88
};
89

90
template<typename T>
91
class Vector<T, EIGEN_BASE>
92
{
93
	//! Eigen vector
94 95
	mutable Eigen::Matrix<T, Eigen::Dynamic, 1> v;

96
	//! row val vector
97
	mutable openfpm::vector<rval<T,EIGEN_RVAL>> row_val;
98 99

	//! row val vector received
100 101
	mutable openfpm::vector<rval<T,EIGEN_RVAL>> row_val_recv;

102
	//! global to local map
103 104
	mutable std::unordered_map<size_t,size_t> map;

105
	//! invalid
106 107
	T invalid;

108
	//! Processors from where we gather
109 110 111 112
	mutable openfpm::vector<size_t> prc;

	//size of each chunk
	mutable openfpm::vector<size_t> sz;
113

114

Pietro Incardona's avatar
Pietro Incardona committed
115
	/*! \brief Here we collect the full vector on master
116 117 118 119
	 *
	 */
	void collect() const
	{
incardon's avatar
incardon committed
120
		Vcluster<> & vcl = create_vcluster();
121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158

		row_val_recv.clear();

		// here we collect all the triplet in one array on the root node
		vcl.SGather(row_val,row_val_recv,prc,sz,0);

		if (vcl.getProcessUnitID() != 0)
			row_val.resize(0);
		else
			row_val.swap(row_val_recv);

		build_map();
	}

	/*! \brief Set the Eigen internal vector
	 *
	 *
	 */
	void setEigen() const
	{
		// set the vector

		for (size_t i = 0 ; i < row_val.size() ; i++)
			v[row_val.get(i).row()] = row_val.get(i).value();
	}

	/*! \brief Build the map
	 *
	 *
	 */
	void build_map() const
	{
		map.clear();

		for (size_t i = 0 ; i < row_val.size() ; i++)
			map[row_val.get(i).row()] = i;
	}

159 160
public:

161 162 163 164 165 166
	/*! \brief Copy the vector
	 *
	 * \param v vector to copy
	 *
	 */
	Vector(const Vector<T> & v)
167
	:invalid(0)
168 169 170 171 172 173 174 175 176 177
	{
		this->operator=(v);
	}

	/*! \brief Copy the vector
	 *
	 * \param v vector to copy
	 *
	 */
	Vector(const Vector<T> && v)
178
	:invalid(0)
179 180 181 182 183 184 185 186 187 188 189
	{
		this->operator=(v);
	}

	/*! \brief Create a vector with n elements
	 *
	 * \param n number of elements in the vector
	 *
	 */
	Vector(size_t n)
	{
190
		resize(n,0);
191 192 193 194 195 196 197 198 199
	}

	/*! \brief Create a vector with 0 elements
	 *
	 */
	Vector()
	{
	}

200 201 202
	/*! \brief Resize the Vector
	 *
	 * \param row numbers of row
203
	 * \param l_row unused
204 205
	 *
	 */
206
	void resize(size_t row, size_t l_row)
207 208 209 210 211 212 213
	{
		v.resize(row);
	}

	/*! \brief Return a reference to the vector element
	 *
	 * \param i element
214
	 * \param val value
215 216
	 *
	 */
217
	void insert(size_t i, T val)
218
	{
219 220 221 222 223 224 225 226 227
		row_val.add();

		// Map
		map[i] = row_val.size()-1;

		row_val.last().row() = i;
		row_val.last().value() = val;
	}

228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263
	/*! \brief Return a reference to the vector element
	 *
	 * \param i element
	 *
	 * \return reference to the element vector
	 *
	 */
	inline T & insert(size_t i)
	{
		row_val.add();

		// Map
		map[i] = row_val.size()-1;

		row_val.last().row() = i;
		return row_val.last().value();
	}

	/*! \brief Return a reference to the vector element
	 *
	 * \param i element
	 *
	 * \return reference to the element vector
	 *
	 */
	inline const T & insert(size_t i) const
	{
		row_val.add();

		// Map
		map[i] = row_val.size()-1;

		row_val.last().row() = i;
		return row_val.last().value();
	}

264 265 266 267 268 269 270 271 272
	/*! \brief Return a reference to the vector element
	 *
	 * \warning The element must exist
	 *
	 * \param i element
	 *
	 * \return reference to the element vector
	 *
	 */
273
	const T & operator()(size_t i) const
274 275 276 277 278
	{
		// Search if exist

		std::unordered_map<size_t,size_t>::iterator it = map.find(i);

279 280
		if ( it != map.end() )
			return row_val.get(it->second).value();
281

282 283
		return insert(i);
	}
284

285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300
	/*! \brief Return a reference to the vector element
	 *
	 * \warning The element must exist
	 *
	 * \param i element
	 *
	 * \return reference to the element vector
	 *
	 */
	T & operator()(size_t i)
	{
		// Search if exist

		std::unordered_map<size_t,size_t>::iterator it = map.find(i);

		if ( it != map.end() )
301 302
			return row_val.get(it->second).value();

303
		return insert(i);
304 305 306 307 308 309 310 311 312
	}

	/*! \brief Get the Eigen Vector object
	 *
	 * \return the Eigen Vector
	 *
	 */
	const Eigen::Matrix<T, Eigen::Dynamic, 1> & getVec() const
	{
313 314 315
		collect();
		setEigen();

316 317 318 319 320 321 322 323 324 325
		return v;
	}

	/*! \brief Get the Eigen Vector object
	 *
	 * \return the Eigen Vector
	 *
	 */
	Eigen::Matrix<T, Eigen::Dynamic, 1> & getVec()
	{
326 327 328
		collect();
		setEigen();

329 330 331
		return v;
	}

332 333 334 335 336 337 338 339 340
	/*! \brief Scatter the vector information to the other processors
	 *
	 * Eigen does not have a real parallel vector, so in order to work we have to scatter
	 * the vector from one processor to the other
	 *
	 *
	 */
	void scatter()
	{
341
		row_val_recv.clear();
incardon's avatar
incardon committed
342
		Vcluster<> & vcl = create_vcluster();
343 344

		vcl.SScatter(row_val,row_val_recv,prc,sz,0);
345

346 347 348 349 350 351 352 353
		// if we do not receive anything a previous collect has not been performed
		// and so nothing is scattered
		if (row_val_recv.size() != 0)
		{
			row_val.clear();
			row_val.add(row_val_recv);
			build_map();
		}
354 355
	}

356 357 358 359 360 361 362 363 364
	/*! \brief Load from file
	 *
	 *
	 *
	 */
	void fromFile(std::string file)
	{
		std::ifstream inputf;
		inputf.open(file);
365

366 367 368 369
		for (size_t i = 0 ; i < v.size() ; i++)
			inputf >> v(i);

		inputf.close();
370 371

	}
372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416

	/*! \brief Copy the vector
	 *
	 * \param v vector to copy
	 *
	 */
	Vector<T> & operator=(const Vector<T> & v)
	{
		prc = v.prc;
		sz = v.sz;
		map = v.map;
		row_val = v.row_val;

		return *this;
	}

	/*! \brief Copy the vector
	 *
	 * \param v vector to copy
	 *
	 */
	Vector<T> & operator=(const Vector<T> && v)
	{
		prc = v.prc;
		sz = v.sz;
		map = v.map;
		row_val = v.row_val;

		return *this;
	}

	/*! \brief Copy the vector (it is used for special purpose)
	 *
	 * \warning v MUST contain at least all the elements of the vector
	 *
	 * \param v base eigen vector to copy
	 *
	 */
	Vector<T> & operator=(Eigen::Matrix<T, Eigen::Dynamic, 1> & v)
	{
		for (size_t i = 0 ; i < row_val.size() ; i++)
			row_val.get(i).value() = v(row_val.get(i).row());

		return *this;
	}
417 418 419 420
};


#endif /* OPENFPM_NUMERICS_SRC_VECTOR_VECTOR_EIGEN_HPP_ */