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

#ifndef OPENFPM_NUMERICS_SRC_MATRIX_SPARSEMATRIX_EIGEN_HPP_
#define OPENFPM_NUMERICS_SRC_MATRIX_SPARSEMATRIX_EIGEN_HPP_

#include "Vector/map_vector.hpp"
#include <boost/mpl/int.hpp>
13
#include "VCluster.hpp"
14 15 16

#define EIGEN_TRIPLET 1

Pietro Incardona's avatar
Pietro Incardona committed
17 18 19 20 21 22
/*! \brief It store one non-zero element in the sparse matrix
 *
 * Given a row, and a column, store a value
 *
 *
 */
23
template<typename T>
24
class triplet<T,EIGEN_TRIPLET>
25
{
26 27 28
	Eigen::Triplet<T,long int> triplet_ei;

public:
29 30 31

	long int & row()
	{
32
		size_t ptr = (size_t)&triplet_ei.row();
33 34 35 36 37 38

		return (long int &)*(long int *)ptr;
	}

	long int & col()
	{
39
		size_t ptr = (size_t)&triplet_ei.col();
40 41 42 43 44 45

		return (long int &)*(long int *)ptr;
	}

	T & value()
	{
46
		size_t ptr = (size_t)&triplet_ei.value();
47 48 49

		return (T &)*(T *)ptr;
	}
50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66

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

	// Default constructor
	triplet()	{};
67 68
};

Pietro Incardona's avatar
Pietro Incardona committed
69 70 71 72 73 74 75
/* ! \brief Sparse Matrix implementation, that map over Eigen
 *
 * \tparam T Type of the sparse Matrix store on each row,colums
 * \tparam id_t type of id
 * \tparam impl implementation
 *
 */
76 77 78 79 80
template<typename T, typename id_t>
class SparseMatrix<T,id_t,Eigen::SparseMatrix<T,0,id_t>>
{
public:

81
	//! Triplet implementation id
82 83
	typedef boost::mpl::int_<EIGEN_TRIPLET> triplet_impl;

84
	//! Triplet type
85 86
	typedef triplet<T,EIGEN_TRIPLET> triplet_type;

87 88 89 90
private:

	Eigen::SparseMatrix<T,0,id_t> mat;
	openfpm::vector<triplet_type> trpl;
91
	openfpm::vector<triplet_type> trpl_recv;
92 93 94 95 96

	/*! \brief Assemble the matrix
	 *
	 *
	 */
97
	void assemble()
98
	{
99
		Vcluster & vcl = create_vcluster();
100 101 102 103 104

		////// On Master and only here
		// we assemble the Matrix from the collected data
		if (vcl.getProcessingUnits() != 1)
		{
105 106 107 108
			collect();
			// only master assemble the Matrix
			if (vcl.getProcessUnitID() == 0)
				mat.setFromTriplets(trpl_recv.begin(),trpl_recv.end());
109 110
		}
		else
111
			mat.setFromTriplets(trpl.begin(),trpl.end());
112 113 114 115 116 117 118
	}

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

121
		trpl_recv.clear();
122

123 124
		// here we collect all the triplet in one array on the root node
		vcl.SGather(trpl,trpl_recv,0);
125 126 127 128 129 130
	}

public:



131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148
	/*! \brief Create an empty Matrix
	 *
	 * \param N1 number of row
	 * \param N2 number of colums
	 *
	 */
	SparseMatrix(size_t N1, size_t N2)
	:mat(N1,N2)
	{
	}

	/*! \brief Create a Matrix from a set of triplet
	 *
	 * \param N1 number of row
	 * \param N2 number of colums
	 * \param trpl triplet set
	 *
	 */
149
	SparseMatrix(size_t N1, size_t N2, const openfpm::vector<Eigen::Triplet<T,id_t>> & trpl)
150 151
	:mat(N1,N2)
	{
152
		this->trpl = trpl;
153 154 155 156 157 158 159
	}

	/*! \brief Create an empty Matrix
	 *
	 */
	SparseMatrix()	{}

160 161 162
	/*! \brief Get the Matrix triplets bugger
	 *
	 * It return a buffer that can be filled with triplets
163
	 *
164
	 * \return triplet buffer
165 166
	 *
	 */
167
	openfpm::vector<triplet_type> & getMatrixTriplets()
168
	{
169
		return this->trpl;
170 171 172 173 174 175 176 177 178
	}

	/*! \brief Get the Eigen Matrix object
	 *
	 * \return the Eigen Matrix
	 *
	 */
	const Eigen::SparseMatrix<T,0,id_t> & getMat() const
	{
179 180 181
		// Here we collect the information on master
		assemble();

182 183 184 185 186 187 188 189 190 191
		return mat;
	}

	/*! \brief Get the Eigen Matrix object
	 *
	 * \return the Eigen Matrix
	 *
	 */
	Eigen::SparseMatrix<T,0,id_t> & getMat()
	{
192 193
		assemble();

194 195 196 197 198 199 200 201 202 203 204 205 206
		return mat;
	}

	/*! \brief Resize the Sparse Matrix
	 *
	 * \param row number for row
	 * \param col number of colums
	 *
	 */
	void resize(size_t row, size_t col)
	{
		mat.resize(row,col);
	}
207 208 209 210 211 212 213 214 215 216

	/*! \brief Get the row i and the colum j of the Matrix
	 *
	 * \return the value of the matrix at row i colum j
	 *
	 */
	T operator()(id_t i, id_t j)
	{
		return mat.coeffRef(i,j);
	}
Pietro Incardona's avatar
Pietro Incardona committed
217 218 219 220 221 222 223 224 225 226 227 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 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313

	/*! \brief Save this object into file
	 *
	 * \param file filename
	 *
	 * \return true if succed
	 *
	 */
	bool save(const std::string & file) const
	{
		std::vector<size_t> pap_prp;

		Packer<decltype(trpl),HeapMemory>::packRequest(trpl,pap_prp);

		// Calculate how much preallocated memory we need to pack all the objects
		size_t req = ExtPreAlloc<HeapMemory>::calculateMem(pap_prp);

		// allocate the memory
		HeapMemory pmem;
		pmem.allocate(req);
		ExtPreAlloc<HeapMemory> mem(pap_prp,pmem);

		//Packing

		Pack_stat sts;
		Packer<decltype(trpl),HeapMemory>::pack(mem,trpl,sts);

		// Save into a binary file
	    std::ofstream dump (file, std::ios::out | std::ios::binary);
	    if (dump.is_open() == false)
	    	return false;
	    dump.write ((const char *)pmem.getPointer(), pmem.size());

	    return true;
	}

	/*! \brief Load this object from file
	 *
	 * \param file filename
	 *
	 * \return true if succed
	 *
	 */
	bool load(const std::string & file)
	{
	    std::ifstream fs (file, std::ios::in | std::ios::binary | std::ios::ate );
	    if (fs.is_open() == false)
	    	return false;

	    // take the size of the file
	    size_t sz = fs.tellg();

	    fs.close();

	    // reopen the file without ios::ate to read
	    std::ifstream input (file, std::ios::in | std::ios::binary );
	    if (input.is_open() == false)
	    	return false;

	    // Create the HeapMemory and the ExtPreAlloc memory
	    std::vector<size_t> pap_prp;
	    pap_prp.push_back(sz);
	    HeapMemory pmem;
		ExtPreAlloc<HeapMemory> mem(pap_prp,pmem);

		// read
	    input.read((char *)pmem.getPointer(), sz);

	    //close the file
	    input.close();

		//Unpacking
		Unpack_stat ps;

	 	Unpacker<decltype(trpl),HeapMemory>::unpack(mem,trpl,ps);

	 	return true;
	}

	/*! \brief Get the value from triplet
	 *
	 * \warning It is extremly slow because it do a full search across the triplets elements
	 *
	 * \param r row
	 * \param c colum
	 *
	 */
	T getValue(size_t r, size_t c)
	{
		for (size_t i = 0 ; i < trpl.size() ; i++)
		{
			if (r == (size_t)trpl.get(i).row() && c == (size_t)trpl.get(i).col())
				return trpl.get(i).value();
		}

		return 0;
	}
314 315 316 317 318
};



#endif /* OPENFPM_NUMERICS_SRC_MATRIX_SPARSEMATRIX_EIGEN_HPP_ */