OpenFPM_pdata  4.1.0
Project that contain the implementation of distributed structures
Vector_petsc.hpp
1 /*
2  * Vector_petsc.hpp
3  *
4  * Created on: Apr 29, 2016
5  * Author: i-bird
6  */
7 
8 #ifndef OPENFPM_NUMERICS_SRC_VECTOR_VECTOR_PETSC_HPP_
9 #define OPENFPM_NUMERICS_SRC_VECTOR_VECTOR_PETSC_HPP_
10 
11 #include "Vector/map_vector.hpp"
12 #include "Vector/vector_def.hpp"
13 #include <boost/mpl/int.hpp>
14 #include <petscvec.h>
15 #include "util/petsc_util.hpp"
16 #include <unordered_map>
17 #include "VTKWriter/VTKWriter.hpp"
18 #include "CSVWriter/CSVWriter.hpp"
19 
20 #define PETSC_RVAL 2
21 
28 template<typename T>
29 class rval<T,PETSC_RVAL>
30 {
31 public:
32 
34  typedef boost::fusion::vector<PetscInt,T> type;
35 
38 
40  static const unsigned int row = 0;
41 
43  static const unsigned int value = 1;
44 
46  static const unsigned int max_prop = 2;
47 
53  long int & rw()
54  {
55  return boost::fusion::at_c<row>(data);
56  }
57 
63  T & val()
64  {
65  return boost::fusion::at_c<value>(data);
66  }
67 
71  rval() {}
72 
79  rval(long int i, T val)
80  {
81  rw() = i;
82  val() = val;
83  }
84 
90  static inline bool noPointers()
91  {
92  return true;
93  }
94 };
95 
96 constexpr unsigned int row_id = 0;
97 constexpr unsigned int val_id = 1;
98 
99 
105 template<typename T>
106 class Vector<T,PETSC_BASE>
107 {
109  size_t n_row;
110 
112  size_t n_row_local;
113 
115  mutable bool v_created = false;
116 
118  mutable Vec v;
119 
122 
124  mutable std::unordered_map<size_t,size_t> map;
125 
128 
133  void setPetsc() const
134  {
135  if (v_created == false)
136  {PETSC_SAFE_CALL(VecSetType(v,VECMPI));}
137 
138  // set the vector
139 
140  if (row_val.size() != 0)
141  {PETSC_SAFE_CALL(VecSetValues(v,row_val.size(),&row_val.template get<row_id>(0),&row_val.template get<val_id>(0),INSERT_VALUES))}
142 
143  PETSC_SAFE_CALL(VecAssemblyBegin(v));
144  PETSC_SAFE_CALL(VecAssemblyEnd(v));
145 
146  v_created = true;
147  }
148 
149 public:
150 
157  :Vector()
158  {
159  this->operator=(v);
160  }
161 
168  :Vector()
169  {
170  this->operator=(v);
171  }
172 
178  {
179  if (is_openfpm_init() == true)
180  {PETSC_SAFE_CALL(VecDestroy(&v));}
181  }
182 
189  Vector(size_t n, size_t n_row_local)
190  :n_row_local(n_row_local),v(NULL),invalid(0)
191  {
192  // Create the vector
193  PETSC_SAFE_CALL(VecCreate(PETSC_COMM_WORLD,&v));
194 
195  resize(n,n_row_local);
196  }
197 
202  :n_row(0),n_row_local(0),invalid(0)
203  {
204  // Create the vector
205  PETSC_SAFE_CALL(VecCreate(PETSC_COMM_WORLD,&v));
206  }
207 
214  void resize(size_t row, size_t l_row)
215  {
216  n_row = row;
217  n_row_local = l_row;
218 
219  PETSC_SAFE_CALL(VecSetSizes(v,n_row_local,n_row));
220  }
221 
228  void insert(size_t i, T val)
229  {
230  row_val.add();
231 
232  // Map
233  map[i] = row_val.size()-1;
234 
235  row_val.last().template get<row_id>() = i;
236  row_val.last().template get<val_id>() = val;
237  }
238 
246  inline PetscScalar & insert(size_t i)
247  {
248  row_val.add();
249 
250  // Map
251  map[i] = row_val.size()-1;
252 
253  row_val.last().template get<row_id>() = i;
254  return row_val.last().template get<val_id>();
255  }
256 
264  inline const PetscScalar & insert(size_t i) const
265  {
266  row_val.add();
267 
268  // Map
269  map[i] = row_val.size()-1;
270 
271  row_val.last().template get<row_id>() = i;
272  return row_val.last().template get<val_id>();
273  }
274 
284  const PetscScalar & operator()(size_t i) const
285  {
286  // Search if exist
287 
288  std::unordered_map<size_t,size_t>::iterator it = map.find(i);
289 
290  if ( it != map.end() )
291  return row_val.template get<val_id>(it->second);
292 
293  return insert(i);
294  }
295 
305  PetscScalar & operator()(size_t i)
306  {
307  // Search if exist
308 
309  std::unordered_map<size_t,size_t>::iterator it = map.find(i);
310 
311  if ( it != map.end() )
312  return row_val.template get<val_id>(it->second);
313 
314  return insert(i);
315  }
316 
322  const Vec & getVec() const
323  {
324  setPetsc();
325 
326  return v;
327  }
328 
334  Vec & getVec()
335  {
336  setPetsc();
337 
338  return v;
339  }
340 
344  void update()
345  {
346  PetscInt n_row;
347  PetscInt n_row_local;
348 
349  // Get the size of the vector from PETSC
350  VecGetSize(v,&n_row);
351  VecGetLocalSize(v,&n_row_local);
352 
353  this->n_row = n_row;
354  this->n_row_local = n_row_local;
355 
356  row_val.resize(n_row_local);
357 
358  //
359 
360  PetscInt low;
361  PetscInt high;
362 
363  VecGetOwnershipRange(v,&low,&high);
364 
365  // Fill the index and construct the map
366 
367  size_t k = 0;
368  for (size_t i = low ; i < (size_t)high ; i++)
369  {
370  row_val.template get<row_id>(k) = i;
371  map[i] = k;
372  k++;
373  }
374 
375  PETSC_SAFE_CALL(VecGetValues(v,row_val.size(),&row_val.template get<row_id>(0),&row_val.template get<val_id>(0)))
376  }
377 
384  {
385  map = v.map;
386  row_val = v.row_val;
387 
388  return *this;
389  }
390 
397  {
398  map.swap(v.map);
399  row_val.swap(v.row_val);
400 
401  return *this;
402  }
403 
408  void setZero()
409  {
410  if (v_created == false)
411  {PETSC_SAFE_CALL(VecSetType(v,VECMPI));}
412 
413  v_created = true;
414  }
415 
416  /* Write vector on vtk
417  *
418  * \param out file to write into
419  *
420  */
421  bool write(std::string out, size_t opt = VTK_WRITER)
422  {
423  Vcluster<> & v_cl = create_vcluster();
424 
427 
428  row_col.resize(n_row_local);
429  values.resize(n_row_local);
430 
431  int i = 0;
432  for (auto it = map.begin() ; it != map.end() ; it++, i++)
433  {
434  row_col.template get<0>(i)[1] = it->first;
435  row_col.template get<0>(i)[0] = 0.0;
436 
437  values.template get<0>(i) = row_val.template get<1>(it->second);
438  }
439 
440  if (opt == VTK_WRITER)
441  {
442  auto ft = file_type::ASCII;
443 
444  // VTKWriter for a set of points
447  VECTOR_POINTS> vtk_writer;
448 
449  vtk_writer.add(row_col,values,row_col.size());
450 
451  std::string output = std::to_string(out + "_" + std::to_string(v_cl.getProcessUnitID()) + std::to_string(".vtk"));
452 
454  prp_names.add("value");
455 
456  // Write the VTK file
457  return vtk_writer.write(output,prp_names,"vector","",ft);
458  }
459  else
460  {
461  // CSVWriter test
464 
465  std::string output = std::to_string(out + "_" + std::to_string(v_cl.getProcessUnitID()) + std::to_string(".csv"));
466 
467  // Write the CSV
468  return csv_writer.write(output,row_col,values);
469  }
470  }
471 };
472 
473 
474 #endif /* OPENFPM_NUMERICS_SRC_VECTOR_VECTOR_EIGEN_HPP_ */
475 
PetscScalar & operator()(size_t i)
Return a reference to the vector element.
void insert(size_t i, T val)
stub insert
Definition: Vector.hpp:105
const PetscScalar & insert(size_t i) const
Return a reference to the vector element.
void update()
Update the Vector with the PETSC object.
const PetscScalar & operator()(size_t i) const
Return a reference to the vector element.
size_t n_row
Number of row the petsc vector has.
rval(long int i, T val)
Constructor from row, column and value.
Vector(const Vector< T, PETSC_BASE > &v)
Copy the vector.
std::unordered_map< size_t, size_t > map
Global to local map.
size_t getProcessUnitID()
Get the process unit id.
static bool noPointers()
Indicate that the structure has no pointer.
void insert(size_t i, T val)
Return a reference to the vector element.
void resize(size_t row, size_t l_row)
Resize the Vector.
T & val()
Get the value.
Sparse Matrix implementation stub object when OpenFPM is compiled with no linear algebra support.
Definition: Vector.hpp:39
PETSC vector for linear algebra.
Vector(size_t n, size_t n_row_local)
Create a vector with n elements.
Vec & getVec()
Get the PETSC Vector object.
openfpm::vector< rval< PetscScalar, PETSC_RVAL >, HeapMemory, memory_traits_inte > row_val
Mutable row value vector.
size_t size()
Stub size.
Definition: map_vector.hpp:211
PetscScalar & insert(size_t i)
Return a reference to the vector element.
~Vector()
Destroy the vector.
CSV Writer.
Definition: CSVWriter.hpp:163
This class allocate, and destroy CPU memory.
Definition: HeapMemory.hpp:39
Transform the boost::fusion::vector into memory specification (memory_traits)
Definition: memory_conf.hpp:83
Implementation of VCluster class.
Definition: VCluster.hpp:58
void setPetsc() const
Set the Eigen internal vector.
Vector< T > & operator=(const Vector< T > &v)
stub operator=
Definition: Vector.hpp:161
bool write(std::string file, v_pos &v, v_prp &prp, size_t offset=0)
It write a CSV file.
Definition: CSVWriter.hpp:248
long int & rw()
Get the row.
type data
structure that store the data of the point
void resize(size_t row, size_t row_n)
stub resize
Definition: Vector.hpp:97
const Vec & getVec() const
Get the PETSC Vector object.
Vector(Vector< T, PETSC_BASE > &&v)
Copy the vector.
It store one row value of a vector.
Definition: Vector.hpp:21
size_t n_row_local
Number of local rows.
boost::fusion::vector< PetscInt, T > type
boost fusion that store the point
rval()
Default constructor.
Vector()
Create a vector with 0 elements.
Vec v
Mutable vector.
Vector< T, PETSC_BASE > & operator=(const Vector< T, PETSC_BASE > &v)
Copy the vector.
Implementation of 1-D std::vector like structure.
Definition: map_vector.hpp:202
Vector< T, PETSC_BASE > & operator=(Vector< T, PETSC_BASE > &&v)
Copy the vector.
void setZero()
Set to zero all the entries.