OpenFPM_pdata  1.1.0
Project that contain the implementation of distributed structures
 All Data Structures Namespaces Functions Variables Typedefs Enumerations Friends Pages
FDScheme.hpp
1 /*
2  * FiniteDifferences.hpp
3  *
4  * Created on: Sep 17, 2015
5  * Author: i-bird
6  */
7 
8 #ifndef OPENFPM_NUMERICS_SRC_FINITEDIFFERENCE_FDSCHEME_HPP_
9 #define OPENFPM_NUMERICS_SRC_FINITEDIFFERENCE_FDSCHEME_HPP_
10 
11 #include "../Matrix/SparseMatrix.hpp"
12 #include "Grid/grid_dist_id.hpp"
13 #include "Grid/Iterators/grid_dist_id_iterator_sub.hpp"
14 #include "eq.hpp"
15 #include "NN/CellList/CellDecomposer.hpp"
16 #include "Grid/staggered_dist_grid_util.hpp"
17 #include "Grid/grid_dist_id.hpp"
18 #include "Vector/Vector_util.hpp"
19 #include "Grid/staggered_dist_grid.hpp"
20 
123 template<typename Sys_eqs>
124 class FDScheme
125 {
126 public:
127 
129  typedef grid_dist_id<Sys_eqs::dims,typename Sys_eqs::stype,aggregate<size_t>,typename Sys_eqs::b_grid::decomposition::extended_type> g_map_type;
130 
132  typedef Sys_eqs Sys_eqs_typ;
133 
134 private:
135 
137  struct constant_b
138  {
140  typename Sys_eqs::stype scal;
141 
147  constant_b(typename Sys_eqs::stype scal)
148  {
149  this->scal = scal;
150  }
151 
161  typename Sys_eqs::stype get(grid_dist_key_dx<Sys_eqs::dims> & key)
162  {
163  return scal;
164  }
165  };
166 
168  template<typename grid, unsigned int prp>
169  struct grid_b
170  {
172  grid & gr;
173 
180  :gr(gr)
181  {}
182 
190  typename Sys_eqs::stype get(grid_dist_key_dx<Sys_eqs::dims> & key)
191  {
192  return gr.template get<prp>(key);
193  }
194  };
195 
198 
200  typedef typename Sys_eqs::SparseMatrix_type::triplet_type triplet;
201 
203  typename Sys_eqs::Vector_type b;
204 
207 
209  typename Sys_eqs::stype spacing[Sys_eqs::dims];
210 
213 
215  size_t row;
216 
218  size_t row_b;
219 
222 
224  comb<Sys_eqs::dims> s_pos[Sys_eqs::nvar];
225 
230  size_t s_pnt;
231 
235  struct key_and_eq
236  {
239 
241  size_t eq;
242  };
243 
252  {
253  key_and_eq ke;
254  auto it = g_map.getDomainIterator();
255 
256  while (it.isNext())
257  {
258  size_t row_low = g_map.template get<0>(it.get());
259 
260  if (row >= row_low * Sys_eqs::nvar && row < row_low * Sys_eqs::nvar + Sys_eqs::nvar)
261  {
262  ke.eq = row - row_low * Sys_eqs::nvar;
263  ke.key = g_map.getGKey(it.get());
264  return ke;
265  }
266 
267  ++it;
268  }
269  std::cerr << "Error: " << __FILE__ << ":" << __LINE__ << " the row does not map to any position" << "\n";
270 
271  return ke;
272  }
273 
282  inline const std::vector<size_t> padding( const size_t (& sz)[Sys_eqs::dims], Padding<Sys_eqs::dims> & pd)
283  {
284  std::vector<size_t> g_sz_pad(Sys_eqs::dims);
285 
286  for (size_t i = 0 ; i < Sys_eqs::dims ; i++)
287  g_sz_pad[i] = sz[i] + pd.getLow(i) + pd.getHigh(i);
288 
289  return g_sz_pad;
290  }
291 
295  void consistency()
296  {
297  openfpm::vector<triplet> & trpl = A.getMatrixTriplets();
298 
299  // A and B must have the same rows
300  if (row != row_b)
301  std::cerr << "Error " << __FILE__ << ":" << __LINE__ << "the term B and the Matrix A for Ax=B must contain the same number of rows\n";
302 
303  // Indicate all the non zero rows
305  nz_rows.resize(row_b);
306 
307  for (size_t i = 0 ; i < trpl.size() ; i++)
308  nz_rows.get(trpl.get(i).row() - s_pnt*Sys_eqs::nvar) = true;
309 
310  // Indicate all the non zero colums
311  // This check can be done only on single processor
312 
313  Vcluster & v_cl = create_vcluster();
314  if (v_cl.getProcessingUnits() == 1)
315  {
317  nz_cols.resize(row_b);
318 
319  for (size_t i = 0 ; i < trpl.size() ; i++)
320  nz_cols.get(trpl.get(i).col()) = true;
321 
322  // all the rows must have a non zero element
323  for (size_t i = 0 ; i < nz_rows.size() ; i++)
324  {
325  if (nz_rows.get(i) == false)
326  {
327  key_and_eq ke = from_row_to_key(i);
328  std::cerr << "Error: " << __FILE__ << ":" << __LINE__ << " Ill posed matrix row " << i << " is not filled, position " << ke.key.to_string() << " equation: " << ke.eq << "\n";
329  }
330  }
331 
332  // all the colums must have a non zero element
333  for (size_t i = 0 ; i < nz_cols.size() ; i++)
334  {
335  if (nz_cols.get(i) == false)
336  std::cerr << "Error: " << __FILE__ << ":" << __LINE__ << " Ill posed matrix colum " << i << " is not filled\n";
337  }
338  }
339  }
340 
351  template<typename Vct, typename Grid_dst ,unsigned int ... pos> void copy_staggered(Vct & v, Grid_dst & g_dst)
352  {
353  // check that g_dst is staggered
354  if (g_dst.is_staggered() == false)
355  std::cerr << __FILE__ << ":" << __LINE__ << " The destination grid must be staggered " << std::endl;
356 
357 #ifdef SE_CLASS1
358 
359  if (g_map.getLocalDomainSize() != g_dst.getLocalDomainSize())
360  std::cerr << __FILE__ << ":" << __LINE__ << " The staggered and destination grid in size does not match " << std::endl;
361 #endif
362 
363  // sub-grid iterator over the grid map
364  auto g_map_it = g_map.getDomainIterator();
365 
366  // Iterator over the destination grid
367  auto g_dst_it = g_dst.getDomainIterator();
368 
369  while (g_map_it.isNext() == true)
370  {
371  typedef typename to_boost_vmpl<pos...>::type vid;
372  typedef boost::mpl::size<vid> v_size;
373 
374  auto key_src = g_map_it.get();
375  size_t lin_id = g_map.template get<0>(key_src);
376 
377  // destination point
378  auto key_dst = g_dst_it.get();
379 
380  // Transform this id into an id for the Eigen vector
381 
382  copy_ele<Sys_eqs_typ, Grid_dst,Vct> cp(key_dst,g_dst,v,lin_id,g_map.size());
383 
384  boost::mpl::for_each_ref<boost::mpl::range_c<int,0,v_size::value>>(cp);
385 
386  ++g_map_it;
387  ++g_dst_it;
388  }
389  }
390 
401  template<typename Vct, typename Grid_dst ,unsigned int ... pos> void copy_normal(Vct & v, Grid_dst & g_dst)
402  {
403  // check that g_dst is staggered
404  if (g_dst.is_staggered() == true)
405  std::cerr << __FILE__ << ":" << __LINE__ << " The destination grid must be normal " << std::endl;
406 
409 
410  for (size_t i = 0 ; i < Grid_dst::dims ; i++)
411  {
412  start.set_d(i,pd.getLow(i));
413  stop.set_d(i,g_map.size(i) - pd.getHigh(i));
414  }
415 
416  // sub-grid iterator over the grid map
417  auto g_map_it = g_map.getSubDomainIterator(start,stop);
418 
419  // Iterator over the destination grid
420  auto g_dst_it = g_dst.getDomainIterator();
421 
422  while (g_dst_it.isNext() == true)
423  {
424  typedef typename to_boost_vmpl<pos...>::type vid;
425  typedef boost::mpl::size<vid> v_size;
426 
427  auto key_src = g_map_it.get();
428  size_t lin_id = g_map.template get<0>(key_src);
429 
430  // destination point
431  auto key_dst = g_dst_it.get();
432 
433  // Transform this id into an id for the vector
434 
435  copy_ele<Sys_eqs_typ, Grid_dst,Vct> cp(key_dst,g_dst,v,lin_id,g_map.size());
436 
437  boost::mpl::for_each_ref<boost::mpl::range_c<int,0,v_size::value>>(cp);
438 
439  ++g_map_it;
440  ++g_dst_it;
441  }
442  }
443 
457  template<typename bop, typename iterator>
458  void impose_dit_b(bop num,
459  long int id ,
460  const iterator & it_d)
461  {
462 
463  auto it = it_d;
465 
466  // iterate all the grid points
467  while (it.isNext())
468  {
469  // get the position
470  auto key = it.get();
471 
472  b(g_map.template get<0>(key)*Sys_eqs::nvar + id) = num.get(key);
473 
474  // if SE_CLASS1 is defined check the position
475 #ifdef SE_CLASS1
476 // T::position(key,gs,s_pos);
477 #endif
478 
479  ++row_b;
480  ++it;
481  }
482  }
483 
499  template<typename T, typename bop, typename iterator> void impose_dit(const T & op,
500  bop num,
501  long int id ,
502  const iterator & it_d)
503  {
504  openfpm::vector<triplet> & trpl = A.getMatrixTriplets();
505 
506  auto it = it_d;
508 
509  std::unordered_map<long int,typename Sys_eqs::stype> cols;
510 
511  // iterate all the grid points
512  while (it.isNext())
513  {
514  // get the position
515  auto key = it.get();
516 
517  // Add padding
518  for (size_t i = 0 ; i < Sys_eqs::dims ; i++)
519  key.getKeyRef().set_d(i,key.getKeyRef().get(i) + pd.getLow(i));
520 
521  // Calculate the non-zero colums
522  T::value(g_map,key,gs,spacing,cols,1.0);
523 
524  // indicate if the diagonal has been set
525  bool is_diag = false;
526 
527  // create the triplet
528  for ( auto it = cols.begin(); it != cols.end(); ++it )
529  {
530  trpl.add();
531  trpl.last().row() = g_map.template get<0>(key)*Sys_eqs::nvar + id;
532  trpl.last().col() = it->first;
533  trpl.last().value() = it->second;
534 
535  if (trpl.last().row() == trpl.last().col())
536  {is_diag = true;}
537 
538  }
539 
540  // If does not have a diagonal entry put it to zero
541  if (is_diag == false)
542  {
543  trpl.add();
544  trpl.last().row() = g_map.template get<0>(key)*Sys_eqs::nvar + id;
545  trpl.last().col() = g_map.template get<0>(key)*Sys_eqs::nvar + id;
546  trpl.last().value() = 0.0;
547  }
548 
549  b(g_map.template get<0>(key)*Sys_eqs::nvar + id) = num.get(key);
550 
551  cols.clear();
552 
553  // if SE_CLASS1 is defined check the position
554 #ifdef SE_CLASS1
555 // T::position(key,gs,s_pos);
556 #endif
557 
558  ++row;
559  ++row_b;
560  ++it;
561  }
562  }
563 
579  template<typename T, typename bop, typename iterator> void impose_git(const T & op ,
580  bop num,
581  long int id ,
582  const iterator & it_d)
583  {
584  openfpm::vector<triplet> & trpl = A.getMatrixTriplets();
585 
586  auto it = it_d;
588 
589  std::unordered_map<long int,float> cols;
590 
591  // iterate all the grid points
592  while (it.isNext())
593  {
594  // get the position
595  auto key = it.get();
596 
597  // Calculate the non-zero colums
598  T::value(g_map,key,gs,spacing,cols,1.0);
599 
600  // indicate if the diagonal has been set
601  bool is_diag = false;
602 
603  // create the triplet
604  for ( auto it = cols.begin(); it != cols.end(); ++it )
605  {
606  trpl.add();
607  trpl.last().row() = g_map.template get<0>(key)*Sys_eqs::nvar + id;
608  trpl.last().col() = it->first;
609  trpl.last().value() = it->second;
610 
611  if (trpl.last().row() == trpl.last().col())
612  is_diag = true;
613 
614 // std::cout << "(" << trpl.last().row() << "," << trpl.last().col() << "," << trpl.last().value() << ")" << "\n";
615  }
616 
617  // If does not have a diagonal entry put it to zero
618  if (is_diag == false)
619  {
620  trpl.add();
621  trpl.last().row() = g_map.template get<0>(key)*Sys_eqs::nvar + id;
622  trpl.last().col() = g_map.template get<0>(key)*Sys_eqs::nvar + id;
623  trpl.last().value() = 0.0;
624  }
625 
626  b(g_map.template get<0>(key)*Sys_eqs::nvar + id) = num.get(key);
627 
628  cols.clear();
629 
630  // if SE_CLASS1 is defined check the position
631 #ifdef SE_CLASS1
632 // T::position(key,gs,s_pos);
633 #endif
634 
635  ++row;
636  ++row_b;
637  ++it;
638  }
639  }
640 
645  {
646  Vcluster & v_cl = create_vcluster();
647 
648  // Calculate the size of the local domain
649  size_t sz = g_map.getLocalDomainSize();
650 
651  // Get the total size of the local grids on each processors
652  v_cl.allGather(sz,pnt);
653  v_cl.execute();
654  s_pnt = 0;
655 
656  // calculate the starting point for this processor
657  for (size_t i = 0 ; i < v_cl.getProcessUnitID() ; i++)
658  s_pnt += pnt.get(i);
659 
660  // resize b if needed
661  b.resize(Sys_eqs::nvar * g_map.size(),Sys_eqs::nvar * sz);
662 
663  // Calculate the starting point
664 
665  // Counter
666  size_t cnt = 0;
667 
668  // Create the re-mapping grid
669  auto it = g_map.getDomainIterator();
670 
671  while (it.isNext())
672  {
673  auto key = it.get();
674 
675  g_map.template get<0>(key) = cnt + s_pnt;
676 
677  ++cnt;
678  ++it;
679  }
680 
681  // sync the ghost
682  g_map.template ghost_get<0>();
683  }
684 
691  {
692  construct_gmap();
693 
694  // Create a CellDecomposer and calculate the spacing
695 
696  size_t sz_g[Sys_eqs::dims];
697  for (size_t i = 0 ; i < Sys_eqs::dims ; i++)
698  {
699  if (Sys_eqs::boundary[i] == NON_PERIODIC)
700  sz_g[i] = gs.getSize()[i] - 1;
701  else
702  sz_g[i] = gs.getSize()[i];
703  }
704 
705  CellDecomposer_sm<Sys_eqs::dims,typename Sys_eqs::stype> cd(domain,sz_g,0);
706 
707  for (size_t i = 0 ; i < Sys_eqs::dims ; i++)
708  spacing[i] = cd.getCellBox().getHigh(i);
709  }
710 
711 public:
712 
718  void setStagPos(comb<Sys_eqs::dims> (& sp)[Sys_eqs::nvar])
719  {
720  for (size_t i = 0 ; i < Sys_eqs::nvar ; i++)
721  s_pos[i] = sp[i];
722  }
723 
731  void computeStag()
732  {
733  typedef typename Sys_eqs::b_grid::value_type prp_type;
734 
735  openfpm::vector<comb<Sys_eqs::dims>> c_prp[prp_type::max_prop];
736 
738 
739  boost::mpl::for_each_ref< boost::mpl::range_c<int,0,prp_type::max_prop> >(ssp);
740  }
741 
748  {
749  return pd;
750  }
751 
759  const g_map_type & getMap()
760  {
761  return g_map;
762  }
763 
775  const typename Sys_eqs::b_grid & b_g)
776  :pd({0,0,0},{0,0,0}),gs(b_g.getGridInfoVoid()),g_map(b_g,stencil,pd),row(0),row_b(0)
777  {
778  Initialize(domain);
779  }
780 
792  const Ghost<Sys_eqs::dims,long int> & stencil,
794  const typename Sys_eqs::b_grid & b_g)
795  :pd(pd),gs(b_g.getGridInfoVoid()),g_map(b_g,stencil,pd),row(0),row_b(0)
796  {
797  Initialize(domain);
798  }
799 
817  template<typename T> void impose(const T & op,
818  typename Sys_eqs::stype num,
819  long int id,
820  const long int (& start)[Sys_eqs::dims],
821  const long int (& stop)[Sys_eqs::dims],
822  bool skip_first = false)
823  {
826 
827  bool increment = false;
828  if (skip_first == true)
829  {
830  start_k = grid_key_dx<Sys_eqs::dims>(start);
831  stop_k = grid_key_dx<Sys_eqs::dims>(start);
832 
833  auto it = g_map.getSubDomainIterator(start_k,stop_k);
834 
835  if (it.isNext() == true)
836  increment = true;
837  }
838 
839  // add padding to start and stop
840  start_k = grid_key_dx<Sys_eqs::dims>(start);
841  stop_k = grid_key_dx<Sys_eqs::dims>(stop);
842 
843  auto it = g_map.getSubDomainIterator(start_k,stop_k);
844 
845  if (increment == true)
846  ++it;
847 
848  constant_b b(num);
849 
850  impose_git(op,b,id,it);
851 
852  }
853 
857  void new_b()
858  {row_b = 0;}
859 
864  void new_A()
865  {row = 0;}
866 
882  template<unsigned int prp, typename b_term, typename iterator>
883  void impose_dit_b(b_term & b_t,
884  const iterator & it_d,
885  long int id = 0)
886  {
887  grid_b<b_term,prp> b(b_t);
888 
889  impose_dit_b(b,id,it_d);
890  }
891 
907  template<typename T> void impose_dit(const T & op ,
908  typename Sys_eqs::stype num,
909  long int id ,
911  {
912  constant_b b(num);
913 
914  impose_dit(op,b,id,it_d);
915  }
916 
932  template<unsigned int prp, typename T, typename b_term, typename iterator> void impose_dit(const T & op ,
933  b_term & b_t,
934  const iterator & it_d,
935  long int id = 0)
936  {
937  grid_b<b_term,prp> b(b_t);
938 
939  impose_dit(op,b,id,it_d);
940  }
941 
943  typename Sys_eqs::SparseMatrix_type A;
944 
950  typename Sys_eqs::SparseMatrix_type & getA()
951  {
952 #ifdef SE_CLASS1
953  consistency();
954 #endif
955  A.resize(g_map.size()*Sys_eqs::nvar,g_map.size()*Sys_eqs::nvar,
956  g_map.getLocalDomainSize()*Sys_eqs::nvar,
957  g_map.getLocalDomainSize()*Sys_eqs::nvar);
958 
959  return A;
960 
961  }
962 
968  typename Sys_eqs::Vector_type & getB()
969  {
970 #ifdef SE_CLASS1
971  consistency();
972 #endif
973 
974  return b;
975  }
976 
992  template<unsigned int ... pos, typename Vct, typename Grid_dst> void copy(Vct & v,const long int (& start)[Sys_eqs_typ::dims], const long int (& stop)[Sys_eqs_typ::dims], Grid_dst & g_dst)
993  {
995  {
996  if (g_dst.is_staggered() == true)
997  copy_staggered<Vct,Grid_dst,pos...>(v,g_dst);
998  else
999  {
1000  // Create a temporal staggered grid and copy the data there
1001  auto & g_map = this->getMap();
1002 
1003  // Convert the ghost in grid units
1004 
1006  for (size_t i = 0 ; i < Grid_dst::dims ; i++)
1007  {
1008  g_int.setLow(i,g_map.getDecomposition().getGhost().getLow(i) / g_map.spacing(i));
1009  g_int.setHigh(i,g_map.getDecomposition().getGhost().getHigh(i) / g_map.spacing(i));
1010  }
1011 
1013  stg.setDefaultStagPosition();
1014  copy_staggered<Vct,decltype(stg),pos...>(v,stg);
1015 
1016  // sync the ghost and interpolate to the normal grid
1017  stg.template ghost_get<pos...>();
1018  stg.template to_normal<Grid_dst,pos...>(g_dst,this->getPadding(),start,stop);
1019  }
1020  }
1021  else
1022  {
1023  copy_normal<Vct,Grid_dst,pos...>(v,g_dst);
1024  }
1025  }
1026 
1040  template<unsigned int ... pos, typename Vct, typename Grid_dst> void copy(Vct & v, Grid_dst & g_dst)
1041  {
1042  long int start[Sys_eqs_typ::dims];
1043  long int stop[Sys_eqs_typ::dims];
1044 
1045  for (size_t i = 0 ; i < Sys_eqs_typ::dims ; i++)
1046  {
1047  start[i] = 0;
1048  stop[i] = g_dst.size(i);
1049  }
1050 
1051  this->copy<pos...>(v,start,stop,g_dst);
1052  }
1053 };
1054 
1055 #define EQS_FIELDS 0
1056 #define EQS_SPACE 1
1057 
1058 
1059 #endif /* OPENFPM_NUMERICS_SRC_FINITEDIFFERENCE_FDSCHEME_HPP_ */
Padding< Sys_eqs::dims > pd
Padding.
Definition: FDScheme.hpp:197
void setDefaultStagPosition()
It set all the properties defined to be staggered on the default location.
grid_key_dx< Sys_eqs::dims > key
space position
Definition: FDScheme.hpp:238
void copy_staggered(Vct &v, Grid_dst &g_dst)
Copy a given solution vector in a staggered grid.
Definition: FDScheme.hpp:351
size_t row
row of the matrix
Definition: FDScheme.hpp:215
void computeStag()
compute the staggered position for each property
Definition: FDScheme.hpp:731
auto get(const grid_dist_key_dx< dim > &v1) const -> typename std::add_lvalue_reference< decltype(loc_grid.get(v1.getSub()).template get< p >(v1.getKey()))>::type
Get the reference of the selected element.
grid_dist_iterator< dim, device_grid, FREE > getDomainIterator() const
It return an iterator that span the full grid domain (each processor span its local domain) ...
Sys_eqs::Vector_type b
Vector b.
Definition: FDScheme.hpp:203
size_t row_b
row on b
Definition: FDScheme.hpp:218
FDScheme(const Ghost< Sys_eqs::dims, long int > &stencil, const Box< Sys_eqs::dims, typename Sys_eqs::stype > &domain, const typename Sys_eqs::b_grid &b_g)
Constructor.
Definition: FDScheme.hpp:773
g_map_type g_map
mapping grid
Definition: FDScheme.hpp:212
T getLow(int i) const
get the i-coordinate of the low bound interval of the box
Definition: Box.hpp:479
size_t getProcessUnitID()
Get the process unit id.
FDScheme(Padding< Sys_eqs::dims > &pd, const Ghost< Sys_eqs::dims, long int > &stencil, const Box< Sys_eqs::dims, typename Sys_eqs::stype > &domain, const typename Sys_eqs::b_grid &b_g)
Constructor.
Definition: FDScheme.hpp:791
const grid_sm< dim, void > & getGridInfoVoid() const
Get an object containing the grid informations without type.
void execute()
Execute all the requests.
comb< Sys_eqs::dims > s_pos[Sys_eqs::nvar]
Staggered position for each property.
Definition: FDScheme.hpp:224
Sys_eqs::stype spacing[Sys_eqs::dims]
Get the grid spacing.
Definition: FDScheme.hpp:209
grid_b(grid &gr)
gr grid that encapsulate
Definition: FDScheme.hpp:179
const grid_sm< Sys_eqs::dims, void > & gs
Domain Grid informations.
Definition: FDScheme.hpp:206
void impose_dit(const T &op, typename Sys_eqs::stype num, long int id, grid_dist_iterator_sub< Sys_eqs::dims, typename g_map_type::d_grid > it_d)
Impose an operator.
Definition: FDScheme.hpp:907
void copy_normal(Vct &v, Grid_dst &g_dst)
Copy a given solution vector in a normal grid.
Definition: FDScheme.hpp:401
T getHigh(int i) const
get the high interval of the box
Definition: Box.hpp:490
const Padding< Sys_eqs::dims > & getPadding()
Get the specified padding.
Definition: FDScheme.hpp:747
size_t size() const
Return the total number of points in the grid.
void setHigh(int i, T val)
set the high interval of the box
Definition: Box.hpp:467
void impose(const T &op, typename Sys_eqs::stype num, long int id, const long int(&start)[Sys_eqs::dims], const long int(&stop)[Sys_eqs::dims], bool skip_first=false)
Impose an operator.
Definition: FDScheme.hpp:817
size_t size()
Stub size.
Definition: map_vector.hpp:70
std::string to_string()
convert the information into a string
Definition: grid_key.hpp:349
Grid key for a distributed grid.
Sys_eqs::stype scal
scalar
Definition: FDScheme.hpp:140
void new_A()
In case we want to impose a new A re-using FDScheme we have to call This function.
Definition: FDScheme.hpp:864
grid_dist_iterator_sub< dim, device_grid > getSubDomainIterator(const grid_key_dx< dim > &start, const grid_key_dx< dim > &stop) const
It return an iterator that span the grid domain only in the specified part.
Sys_eqs::SparseMatrix_type & getA()
produce the Matrix
Definition: FDScheme.hpp:950
Encapsulation of the b term as constant.
Definition: FDScheme.hpp:137
void setStagPos(comb< Sys_eqs::dims >(&sp)[Sys_eqs::nvar])
set the staggered position for each property
Definition: FDScheme.hpp:718
Definition: Ghost.hpp:39
Implementation of VCluster class.
Definition: VCluster.hpp:36
void copy(Vct &v, Grid_dst &g_dst)
Copy the vector into the grid.
Definition: FDScheme.hpp:1040
void construct_gmap()
Construct the gmap structure.
Definition: FDScheme.hpp:644
grid_key_dx< dim > getGKey(const grid_dist_key_dx< dim > &k)
Convert a g_dist_key_dx into a global key.
This is a distributed grid.
void new_b()
In case we want to impose a new b re-using FDScheme we have to call This function.
Definition: FDScheme.hpp:857
Implementation of the staggered grid.
constant_b(typename Sys_eqs::stype scal)
Constrictor from a scalar.
Definition: FDScheme.hpp:147
this class is a functor for "for_each" algorithm
const g_map_type & getMap()
Return the map between the grid index position and the position in the distributed vector...
Definition: FDScheme.hpp:759
const size_t(& getSize() const)[N]
Return the size of the grid as an array.
Definition: grid_sm.hpp:677
is_grid_staggered analyse T if it has a property grid_type defined and indicate that the grid is stag...
Definition: common.hpp:42
key_and_eq from_row_to_key(size_t row)
From the row Matrix position to the spatial position.
Definition: FDScheme.hpp:251
Distributed grid iterator.
void setLow(int i, T val)
set the low interval of the box
Definition: Box.hpp:456
void impose_dit_b(bop num, long int id, const iterator &it_d)
Impose an operator.
Definition: FDScheme.hpp:458
void impose_dit(const T &op, b_term &b_t, const iterator &it_d, long int id=0)
Impose an operator.
Definition: FDScheme.hpp:932
Encapsulation of the b term as grid.
Definition: FDScheme.hpp:169
This class represent an N-dimensional box.
Definition: Box.hpp:56
Sys_eqs Sys_eqs_typ
Type that specify the properties of the system of equations.
Definition: FDScheme.hpp:132
This class is a trick to indicate the compiler a specific specialization pattern. ...
Definition: memory_c.hpp:201
void consistency()
Check if the Matrix is consistent.
Definition: FDScheme.hpp:295
void impose_dit(const T &op, bop num, long int id, const iterator &it_d)
Impose an operator.
Definition: FDScheme.hpp:499
openfpm::vector< size_t > pnt
Grid points that has each processor.
Definition: FDScheme.hpp:221
const std::vector< size_t > padding(const size_t(&sz)[Sys_eqs::dims], Padding< Sys_eqs::dims > &pd)
calculate the mapping grid size with padding
Definition: FDScheme.hpp:282
Sys_eqs::SparseMatrix_type::triplet_type triplet
Sparse matrix triplet type.
Definition: FDScheme.hpp:200
Sys_eqs::Vector_type & getB()
produce the B vector
Definition: FDScheme.hpp:968
size_t getLocalDomainSize() const
Get the total number of grid points for the calling processor.
void Initialize(const Box< Sys_eqs::dims, typename Sys_eqs::stype > &domain)
Definition: FDScheme.hpp:690
Equation id + space position.
Definition: FDScheme.hpp:235
void set_d(size_t i, mem_id id)
Set the i index.
Definition: grid_key.hpp:407
void impose_git(const T &op, bop num, long int id, const iterator &it_d)
Impose an operator.
Definition: FDScheme.hpp:579
Finite Differences.
Definition: FDScheme.hpp:124
grid_dist_id< Sys_eqs::dims, typename Sys_eqs::stype, aggregate< size_t >, typename Sys_eqs::b_grid::decomposition::extended_type > g_map_type
Distributed grid map.
Definition: FDScheme.hpp:129
this class is a functor for "for_each" algorithm
Definition: Vector_util.hpp:91
void copy(Vct &v, const long int(&start)[Sys_eqs_typ::dims], const long int(&stop)[Sys_eqs_typ::dims], Grid_dst &g_dst)
Copy the vector into the grid.
Definition: FDScheme.hpp:992
size_t getProcessingUnits()
Get the total number of processors.
Decomposition & getDecomposition()
Get the object that store the information about the decomposition.
size_t s_pnt
Definition: FDScheme.hpp:230
bool allGather(T &send, openfpm::vector< T, Mem, gr > &v)
Gather the data from all processors.
void impose_dit_b(b_term &b_t, const iterator &it_d, long int id=0)
Impose an operator.
Definition: FDScheme.hpp:883
size_t eq
equation id
Definition: FDScheme.hpp:241
Sys_eqs::SparseMatrix_type A
type of the sparse matrix
Definition: FDScheme.hpp:943
grid & gr
b term fo the grid
Definition: FDScheme.hpp:172
St spacing(size_t i) const
Get the spacing of the grid in direction i.