Vector_eigen.hpp 6.32 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 14
#include <type_traits>
#include "util/mul_array_extents.hpp"
#include <fstream>
#include "Vector_eigen_util.hpp"
Pietro Incardona's avatar
Pietro Incardona committed
15
#include "Grid/staggered_dist_grid.hpp"
16 17
#include "Space/Ghost.hpp"
#include "FiniteDifference/util/common.hpp"
18
#include <boost/mpl/vector_c.hpp>
19
#include <unordered_map>
20

21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68
#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>
{
	// row
	long int r;

	// value
	T val;

public:

	// Get the row
	long int & row()
	{
		return r;
	}

	// Get the value
	T & value()
	{
		return val;
	}

	/*! \brief Default constructor
	 *
	 */
	rval()	{}

	/*! \brief Constructor from row, colum and value
	 *
	 * \param i row
	 * \param val value
	 *
	 */
	rval(long int i, T val)
	{
		row() = i;
		value() = val;
	}
};
69

70 71
template<typename T>
class Vector<T,Eigen::Matrix<T, Eigen::Dynamic, 1>>
72
{
73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89
	mutable Eigen::Matrix<T, Eigen::Dynamic, 1> v;

	// row value vector
	mutable openfpm::vector<rval<T,EIGEN_RVAL>> row_val;
	mutable openfpm::vector<rval<T,EIGEN_RVAL>> row_val_recv;

	// global to local map
	mutable std::unordered_map<size_t,size_t> map;

	// invalid
	T invalid;

	// Processors from where we gather
	mutable openfpm::vector<size_t> prc;

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

91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135

	/*! \brief Here we collect the full matrix on master
	 *
	 */
	void collect() const
	{
		Vcluster & vcl = *global_v_cluster;

		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;
	}

136 137
public:

138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174
	/*! \brief Copy the vector
	 *
	 * \param v vector to copy
	 *
	 */
	Vector(const Vector<T> & v)
	{
		this->operator=(v);
	}

	/*! \brief Copy the vector
	 *
	 * \param v vector to copy
	 *
	 */
	Vector(const Vector<T> && v)
	{
		this->operator=(v);
	}

	/*! \brief Create a vector with n elements
	 *
	 * \param n number of elements in the vector
	 *
	 */
	Vector(size_t n)
	{
		resize(n);
	}

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

175 176 177 178 179 180 181 182 183 184 185 186 187
	/*! \brief Resize the Vector
	 *
	 * \param row numbers of row
	 *
	 */
	void resize(size_t row)
	{
		v.resize(row);
	}

	/*! \brief Return a reference to the vector element
	 *
	 * \param i element
Pietro Incardona's avatar
Pietro Incardona committed
188
	 * \param val value
189 190
	 *
	 */
191
	void insert(size_t i, T val)
192
	{
193 194 195 196 197 198 199 200 201
		row_val.add();

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

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

Pietro Incardona's avatar
Pietro Incardona committed
202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237
	/*! \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();
	}

238 239 240 241 242 243 244 245 246
	/*! \brief Return a reference to the vector element
	 *
	 * \warning The element must exist
	 *
	 * \param i element
	 *
	 * \return reference to the element vector
	 *
	 */
Pietro Incardona's avatar
Pietro Incardona committed
247
	const T & operator()(size_t i) const
248 249 250 251 252
	{
		// Search if exist

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

Pietro Incardona's avatar
Pietro Incardona committed
253 254
		if ( it != map.end() )
			return row_val.get(it->second).value();
255

Pietro Incardona's avatar
Pietro Incardona committed
256 257
		return insert(i);
	}
258

Pietro Incardona's avatar
Pietro Incardona committed
259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274
	/*! \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() )
275 276
			return row_val.get(it->second).value();

Pietro Incardona's avatar
Pietro Incardona committed
277
		return insert(i);
278 279 280 281 282 283 284 285 286
	}

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

290 291 292 293 294 295 296 297 298 299
		return v;
	}

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

303 304 305
		return v;
	}

306 307 308 309 310 311 312 313 314
	/*! \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()
	{
315 316 317 318
		row_val_recv.clear();
		Vcluster & vcl = *global_v_cluster;

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

320 321 322 323 324 325 326 327
		// 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();
		}
328 329
	}

330 331 332 333 334 335 336 337 338
	/*! \brief Load from file
	 *
	 *
	 *
	 */
	void fromFile(std::string file)
	{
		std::ifstream inputf;
		inputf.open(file);
339

340 341 342 343
		for (size_t i = 0 ; i < v.size() ; i++)
			inputf >> v(i);

		inputf.close();
344 345

	}
346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390

	/*! \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;
	}
391 392 393 394
};


#endif /* OPENFPM_NUMERICS_SRC_VECTOR_VECTOR_EIGEN_HPP_ */