openCARP
Doxygen code documentation for the open cardiac electrophysiology simulator openCARP
SF_container.h
Go to the documentation of this file.
1 // ----------------------------------------------------------------------------
2 // openCARP is an open cardiac electrophysiology simulator.
3 //
4 // Copyright (C) 2020 openCARP project
5 //
6 // This program is licensed under the openCARP Academic Public License (APL)
7 // v1.0: You can use and redistribute it and/or modify it in non-commercial
8 // academic environments under the terms of APL as published by the openCARP
9 // project v1.0, or (at your option) any later version. Commercial use requires
10 // a commercial license (info@opencarp.org).
11 //
12 // This program is distributed without any warranty; see the openCARP APL for
13 // more details.
14 //
15 // You should have received a copy of the openCARP APL along with this program
16 // and can find it online: http://www.opencarp.org/license
17 // ----------------------------------------------------------------------------
18 
27 #ifndef _SF_CONTAINER_H
28 #define _SF_CONTAINER_H
29 
30 #include <cstring>
31 #include <mpi.h>
32 #include <map>
33 #include <math.h>
34 #include <set>
35 
36 #include "dense_mat.hpp"
37 #include "SF_globals.h"
38 #include "SF_vector.h"
39 #include "SF_linalg_utils.h"
40 #include "hashmap.hpp"
41 
42 // This datatypes will be used in conjunction with mesh data. We might want to
43 // have increased / reduced precision when storing meshes, mesh related numberings, etc.
44 // Therefore, we will use (mesh_int_t,mesh_real_t) for data fields that scale with the mesh size,
45 // to be able to tweak the used datatype later. -Aurel
46 typedef int mesh_int_t;
47 //typedef float mesh_real_t;
48 typedef double mesh_real_t;
49 
50 namespace SF {
51 
53 enum elem_t
54 {
55  Tetra = 0,
61  Tri,
62  Line
63 };
64 
66 struct Point {
67  double x;
68  double y;
69  double z;
70 };
71 
72 template<class V>
73 inline Point arr_to_point(V* arr)
74 {
75  return {arr[0], arr[1], arr[2]};
76 }
77 
78 template<class V>
79 inline void point_to_arr(Point & p, V* arr)
80 {
81  arr[0] = p.x, arr[1] = p.y, arr[2] = p.z;
82 }
83 
85 inline Point cross(const Point & a, const Point & b)
86 {
87  Point c = {a.y*b.z - b.y*a.z, b.x*a.z - a.x*b.z, a.x*b.y - a.y*b.x};
88  return c;
89 }
90 
91 inline double inner_prod(const Point & a, const Point & b)
92 {
93  return a.x*b.x + a.y*b.y + a.z*b.z;
94 }
95 
96 inline void outer_prod(const Point & a, const Point & b, const double s, double* buff, const bool add = false)
97 {
98  if(add == false)
99  memset(buff, 0, 9*sizeof(double));
100 
101  buff[0] += a.x*b.x*s, buff[1] += a.x*b.y*s, buff[2] += a.x*b.z*s;
102  buff[3] += a.y*b.x*s, buff[4] += a.y*b.y*s, buff[5] += a.y*b.z*s;
103  buff[6] += a.z*b.x*s, buff[7] += a.z*b.y*s, buff[8] += a.z*b.z*s;
104 }
105 
106 inline void outer_prod(const Point & a, const Point & b, double* buff)
107 {
108  outer_prod(a, b, 1.0, buff);
109 }
110 
112 inline double mag(const Point & vect)
113 {
114  return sqrt(inner_prod(vect, vect));
115 }
116 
117 inline Point normalize(const Point & vect)
118 {
119  double norm = mag(vect);
120  return {vect.x / norm, vect.y / norm, vect.z / norm};
121 }
122 
123 inline Point operator-(const Point & a, const Point & b)
124 {
125  return {a.x - b.x, a.y - b.y, a.z - b.z};
126 }
127 
128 inline Point operator+(const Point & a, const Point & b)
129 {
130  return {a.x + b.x, a.y + b.y, a.z + b.z};
131 }
132 
133 inline Point operator*(const Point & a, const double & s)
134 {
135  return {a.x * s, a.y * s, a.z * s};
136 }
137 
139 inline Point project(const Point & a, const Point & b)
140 {
141  return a * inner_prod(normalize(a), normalize(b));
142 }
143 
144 inline Point orthogonalize(const Point & a, const Point & b)
145 {
146  return a - project(a, b);
147 }
148 
149 
150 inline double distance(const Point & a, const Point & b)
151 {
152  return mag(a - b);
153 }
154 
156 inline elem_t getElemTypeID(char *eletype)
157 {
158  elem_t ret;
159 
160  if ( !strcmp( eletype, "Tt" ) ) {
161  ret = Tetra;
162  } else if ( !strcmp( eletype, "Hx" ) ) {
163  ret = Hexa;
164  } else if ( !strcmp( eletype, "Oc" ) ) {
165  ret = Octa;
166  } else if ( !strcmp( eletype, "Py" ) ) {
167  ret = Pyramid;
168  } else if ( !strcmp( eletype, "Pr" ) ) {
169  ret = Prism;
170  } else if ( !strcmp( eletype, "Qd" ) ) {
171  ret = Quad;
172  } else if ( !strcmp( eletype, "Tr" ) ) {
173  ret = Tri;
174  } else {
175  ret = Line;
176  }
177  return ret;
178 }
179 
180 template<class S, class POINT> inline
181 void array_to_points(const vector<S> & arr, vector<POINT> & pts)
182 {
183  pts.resize(arr.size() / 3);
184  for(size_t i=0; i<pts.size(); i++)
185  pts[i] = {arr[i*3+0], arr[i*3+1], arr[i*3+2]};
186 }
187 
189 typedef enum {
196 
197 // forward declaration of parallel layouts since they are used in mesh.
198 template<class T> class overlapping_layout;
199 template<class T> class non_overlapping_layout;
200 
207 template<class T>
209 {
210  private:
213 
214  public:
217  {}
218 
225  index_mapping(const vector<T> & a, const vector<T> & b)
226  {
227  this->assign(a, b);
228  }
229 
236  inline void assign(const vector<T> & a, const vector<T> & b)
237  {
238  // a 1:1 mapping between old and new indices is required
239  assert( a.size() == b.size() );
240 
241  _fwd.clear();
242  _bwd.clear();
243 
244  for(size_t i=0; i<a.size(); i++)
245  {
246  T aidx = a[i], bidx = b[i];
247  _fwd[aidx] = bidx;
248  _bwd[bidx] = aidx;
249  }
250  }
251 
253  inline T forward_map(T idx) const
254  {
255  typename hashmap::unordered_map<T, T>::const_iterator it = _fwd.find(idx);
256  if(it != _fwd.end())
257  return it->second;
258  else
259  return T(-1);
260  }
261 
263  inline T backward_map(T idx) const
264  {
265  typename hashmap::unordered_map<T, T>::const_iterator it = _bwd.find(idx);
266  if(it != _bwd.end())
267  return it->second;
268  else
269  return T(-1);
270  }
271 
272 
281  inline void forward_map(vector<T> & idx) const
282  {
283  size_t lsize = idx.size();
284  vector<T> perm(lsize);
285  interval(perm, 0, lsize);
286 
287  binary_sort_copy(idx, perm);
288 
289  size_t c = 0;
291 
292  while(c < lsize)
293  {
294  T cc = idx[c];
295  it = _fwd.find(cc);
296 
297  T val = it != _fwd.end() ? it->second : T(-1);
298 
299  while( (c < lsize) && (cc == idx[c]) )
300  {
301  idx[c] = val;
302  c++;
303  }
304  }
305 
306  binary_sort_copy(perm, idx);
307  }
308 
317  inline void backward_map(vector<T> & idx) const
318  {
319  size_t lsize = idx.size();
320  vector<T> perm(lsize);
321  interval(perm, 0, lsize);
322 
323  binary_sort_copy(idx, perm);
324 
325  size_t c = 0;
327 
328  while(c < lsize)
329  {
330  T cc = idx[c];
331  it = _bwd.find(cc);
332 
333  T val = it != _bwd.end() ? *it : T(-1);
334 
335  while( (c < lsize) && (cc == idx[c]) )
336  {
337  idx[c] = val;
338  c++;
339  }
340  }
341 
342  binary_sort_copy(perm, idx);
343  }
344 
346  size_t size() const
347  {
348  return _fwd.size();
349  }
351  bool in_a(const T idx) {
352  return _fwd.count(idx) == 1;
353  }
355  bool in_b(const T idx) {
356  return _bwd.count(idx) == 1;
357  }
358 
360  {
361  return _fwd;
362  }
363 
365  {
366  return _bwd;
367  }
368 
370  {
371  return _fwd;
372  }
373 
375  {
376  return _bwd;
377  }
378 };
379 
383 template<class T, class S>
384 class meshdata
385 {
386 public:
387  size_t g_numelem;
388  size_t l_numelem;
389  size_t g_numpts;
390  size_t l_numpts;
391 
393  MPI_Comm comm;
394 
396  std::string name;
397 
398  // element connectivity : the nodes forming the elements
402 
403  // element data
404  // one value per element
408  // multiple values per element
411 
414 
415  std::map< SF_nbr, vector<T> > nbr;
418 
421 
424  comm(SF_COMM), name("unnamed")
425  {}
426 
435  {
436  if(nbr.count(nbr_type) == 0)
437  nbr[nbr_type] = vector<T>();
438 
439  return nbr[nbr_type];
440  }
441 
454  inline vector<T> & get_numbering(SF_nbr nbr_type)
455  {
456  typename std::map< SF_nbr, vector<T> >::iterator it = nbr.find(nbr_type);
457  assert(it != nbr.end());
458  return it->second;
459  }
472  inline const vector<T> & get_numbering(SF_nbr nbr_type) const
473  {
474  typename std::map< SF_nbr, vector<T> >::const_iterator it = nbr.find(nbr_type);
475  assert(it != nbr.end());
476  return it->second;
477  }
478 
486  inline void localize(SF_nbr nbr_type)
487  {
488  // register numbering
489  vector<T> & nod = this->register_numbering(nbr_type);
490 
492  for(const T & c : this->con) g2l[c] = 0;
493  g2l.sort();
494 
495  // compute numbering
496  nod.resize(g2l.size());
497 
498  size_t widx = 0;
499  for(auto it = g2l.begin(); it != g2l.end(); ++it, widx++) {
500  nod[widx] = it->first;
501  it->second = widx;
502  }
503 
504  global_to_local(g2l, this->con, true);
505  this->l_numpts = nod.size();
506  }
507 
515  inline void globalize(SF_nbr nbr_type)
516  {
517  // register numbering
518  vector<T> & nod = this->get_numbering(nbr_type);
519  for(T & c : con) c = nod[c];
520  }
521 
522 
528  inline void generate_par_layout()
529  {
530  const vector<T> & rnod = this->get_numbering(NBR_REF);
531  const vector<T> & reidx = this->get_numbering(NBR_ELEM_REF);
532 
533  pl.assign(rnod, comm);
534  g_numpts = pl.num_global_idx();
535 
536  epl.assign(reidx, comm);
537  }
538 
542  inline void clear_data()
543  {
544  g_numelem = 0;
545  g_numpts = 0;
546  l_numelem = 0;
547  l_numpts = 0;
548 
549  con.resize(0); con.reallocate();
550  dsp.resize(0); dsp.reallocate();
551  tag.resize(0); tag.reallocate();
552  type.resize(0); type.reallocate();
553 
554  fib.resize(0); fib.reallocate();
555  she.resize(0); she.reallocate();
556  xyz.resize(0); xyz.reallocate();
558 
559  extr_tag.clear();
560  nbr.clear();
561  }
562 
563 
564 };
565 
573 template<class T, class S>
574 inline void nodal_connectivity_graph(const meshdata<T, S> & mesh,
575  vector<T> & n2n_cnt,
576  vector<T> & n2n_con)
577 {
578  vector<T> e2n_cnt, n2e_cnt, n2e_con;
579  const vector<T> & e2n_con = mesh.con;
580 
581  cnt_from_dsp(mesh.dsp, e2n_cnt);
582 
583  transpose_connectivity(e2n_cnt, e2n_con, n2e_cnt, n2e_con);
584  multiply_connectivities(n2e_cnt, n2e_con,
585  e2n_cnt, e2n_con,
586  n2n_cnt, n2n_con);
587 }
588 
598 template<class T, class S>
600 {
601  vector<T> n2n_cnt, n2n_con;
602  nodal_connectivity_graph(mesh, n2n_cnt, n2n_con);
603 
604  int nmax = 0;
605  for(size_t nidx=0; nidx < n2n_cnt.size(); nidx++)
606  if(nmax < n2n_cnt[nidx]) nmax = n2n_cnt[nidx];
607 
608  MPI_Allreduce(MPI_IN_PLACE, &nmax, 1, MPI_INT, MPI_MAX, mesh.comm);
609 
610  return nmax;
611 }
612 
613 template<class T, class S>
614 void get_alg_mask(const meshdata<T,S> & mesh, vector<bool> & alg_mask)
615 {
616  const vector<T> & alg_nod = mesh.pl.algebraic_nodes();
617 
618  alg_mask.assign(mesh.l_numpts, false);
619  for(const T & n : alg_nod) alg_mask[n] = true;
620 }
621 
622 
627 template<class T>
629 {
630  public:
635 
637  inline void resize(size_t size)
638  {
639  scnt.assign(size, T()), rcnt.assign(size, T()), sdsp.assign(size+1, T()), rdsp.assign(size+1, T());
640  }
641 
643  template<class V>
644  inline void scale(V fac)
645  {
646  for(size_t i=0; i<scnt.size(); i++)
647  scnt[i] *= fac, rcnt[i] *= fac;
648 
649  for(size_t i=0; i<sdsp.size(); i++)
650  sdsp[i] *= fac, rdsp[i] *= fac;
651  }
652 
654  inline void transpose()
655  {
656  vector<T> tcnt(scnt), tdsp(sdsp);
657  scnt = rcnt, sdsp = rdsp;
658  rcnt = tcnt, rdsp = tdsp;
659  }
660 
669  template<class V>
670  inline void configure(const vector<V> & dest, MPI_Comm comm)
671  {
672  int size, rank;
673  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
674 
675  this->resize(size);
676 
677  count(dest, scnt);
678  MPI_Alltoall(scnt.data(), sizeof(T), MPI_BYTE, rcnt.data(), sizeof(T), MPI_BYTE, comm);
681  }
691  template<class V>
692  inline void configure(const vector<V> & dest, const vector<T> & cnt, MPI_Comm comm)
693  {
694  int size, rank;
695  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
696 
697  this->resize(size);
698 
699  scnt.zero();
700  for(size_t i=0; i<dest.size(); i++) scnt[dest[i]] += cnt[i];
701 
702  MPI_Alltoall(scnt.data(), sizeof(T), MPI_BYTE, rcnt.data(), sizeof(T), MPI_BYTE, comm);
703 
706  }
716  template<class V>
717  inline void source_ranks(vector<V> & source)
718  {
719  size_t size = rcnt.size();
720  source.resize(sum(rcnt));
721 
722  for(size_t i=0, widx=0; i<size; i++)
723  for(T j=0; j<rcnt[i]; j++, widx++)
724  source[widx] = i;
725  }
726 };
727 
728 template<class T>
729 struct tuple {
730  T v1;
731  T v2;
732 };
733 
734 template<class T, class S>
736 {
737  T v1;
738  S v2;
739 };
740 
741 template<class T>
742 struct triple {
743  T v1;
744  T v2;
745  T v3;
746 };
747 
748 template<class T, class S, class V>
749 struct mixed_triple {
750  T v1;
751  S v2;
752  V v3;
753 };
754 
755 template<class T>
756 struct quadruple {
757  T v1;
758  T v2;
759  T v3;
760  T v4;
761 };
762 
763 template<class T>
764 bool operator<(const struct tuple<T> & lhs, const struct tuple<T> & rhs)
765 {
766  if(lhs.v1 != rhs.v1)
767  return lhs.v1 < rhs.v1;
768  else
769  return lhs.v2 < rhs.v2;
770 }
771 template<class T, class S>
772 bool operator<(const struct mixed_tuple<T, S> & lhs, const struct mixed_tuple<T, S> & rhs)
773 {
774  if(lhs.v1 != rhs.v1)
775  return lhs.v1 < rhs.v1;
776  else
777  return lhs.v2 < rhs.v2;
778 }
779 
780 template<class T>
781 bool operator<(const struct triple<T> & lhs, const struct triple<T> & rhs)
782 {
783  if(lhs.v1 != rhs.v1)
784  return lhs.v1 < rhs.v1;
785  else if(lhs.v2 != rhs.v2)
786  return lhs.v2 < rhs.v2;
787  else
788  return lhs.v3 < rhs.v3;
789 }
790 
791 template<class T>
792 bool operator<(const struct quadruple<T> & lhs, const struct quadruple<T> & rhs)
793 {
794  if(lhs.v1 != rhs.v1)
795  return lhs.v1 < rhs.v1;
796  else if(lhs.v2 != rhs.v2)
797  return lhs.v2 < rhs.v2;
798  else if(lhs.v3 != rhs.v3)
799  return lhs.v3 < rhs.v3;
800  else
801  return lhs.v4 < rhs.v4;
802 }
803 
804 // expand hashing for custom triple struct
805 
806 }
807 
808 namespace hashmap {
809 
810 template<typename T>
811 struct hash_ops< SF::tuple<T> >
812 {
813  static inline
815  return a.v1 == b.v1 && a.v2 == b.v2;
816  }
817 
818  static inline
823  return h;
824  }
825 };
826 
827 template<typename T>
828 struct hash_ops< SF::triple<T> >
829 {
830  static inline
832  return a.v1 == b.v1 && a.v2 == b.v2 && a.v3 == b.v3;
833  }
834 
835  static inline
837  hm_uint h = mkhash_init;
841  return h;
842  }
843 };
844 
845 template<typename T>
846 struct hash_ops< SF::quadruple<T> >
847 {
848  static inline
850  return a.v1 == b.v1 && a.v2 == b.v2 && a.v3 == b.v3 && a.v4 == b.v4;
851  }
852 
853  static inline
855  hm_uint h = mkhash_init;
860  return h;
861  }
862 };
863 }
864 
865 template<class T>
866 struct tri_sele {
867  T v1;
868  T v2;
869  T v3;
870  T eidx;
871 };
872 
873 template<class T>
874 struct quad_sele {
875  T v1;
876  T v2;
877  T v3;
878  T v4;
879  T eidx;
880 };
881 
883 template<class T>
884 void sort_triple(const T in1, const T in2, const T in3, T & out1, T & out2, T & out3)
885 {
886  bool t12 = in1 < in2;
887  bool t13 = in1 < in3;
888  bool t23 = in2 < in3;
889 
890  if(t12) {
891  //(123),(312),(132)
892  if(t23) {
893  //123
894  out1=in1; out2=in2; out3=in3;
895  }
896  else {
897  //(312),(132)
898  if(t13) {
899  //132
900  out1=in1; out2=in3; out3=in2;
901  }
902  else {
903  //312
904  out1=in3; out2=in1; out3=in2;
905  }
906  }
907  }
908  else {
909  //(213),(231),(321)
910  if(t23) {
911  //(213),(231)
912  if(t13) {
913  //213
914  out1=in2; out2=in1; out3=in3;
915  }
916  else {
917  //231
918  out1=in2; out2=in3; out3=in1;
919  }
920  }
921  else {
922  //321
923  out1=in3; out2=in2; out3=in1;
924  }
925  }
926 }
927 
928 #endif
929 
void sort_triple(const T in1, const T in2, const T in3, T &out1, T &out2, T &out3)
sort the "in" triple into the "out" triple
Definition: SF_container.h:884
int mesh_int_t
Definition: SF_container.h:46
double mesh_real_t
Definition: SF_container.h:48
#define SF_COMM
the default SlimFem MPI communicator
Definition: SF_globals.h:27
Sequential linear algebra kernels.
The vector class and related algorithms.
Definition: mesher.cc:237
The class holds the communication graph for a MPI_Exchange() call.
Definition: SF_container.h:629
vector< T > rcnt
Number of elements received from each rank.
Definition: SF_container.h:633
vector< T > scnt
Number of elements sent to each rank.
Definition: SF_container.h:631
void resize(size_t size)
Resize all vectors to size.
Definition: SF_container.h:637
void configure(const vector< V > &dest, MPI_Comm comm)
Set up the communication graph.
Definition: SF_container.h:670
void configure(const vector< V > &dest, const vector< T > &cnt, MPI_Comm comm)
Set up the communication graph.
Definition: SF_container.h:692
vector< T > sdsp
Displacements w.r.t. scnt.
Definition: SF_container.h:632
vector< T > rdsp
Displacements w.r.t. rcnt.
Definition: SF_container.h:634
void source_ranks(vector< V > &source)
For every received data element, get the rank indices it was receive from.
Definition: SF_container.h:717
void transpose()
transpose comm graph (receive becomes send, and vice versa)
Definition: SF_container.h:654
void scale(V fac)
scale comm graph layout data
Definition: SF_container.h:644
Index mapping class. This is a bijective mapping.
Definition: SF_container.h:209
void assign(const vector< T > &a, const vector< T > &b)
Set up the index mapping between a and b.
Definition: SF_container.h:236
const hashmap::unordered_map< T, T > & get_fwd_map() const
Definition: SF_container.h:359
void backward_map(vector< T > &idx) const
Map a whole array of indices in b to indices in a.
Definition: SF_container.h:317
hashmap::unordered_map< T, T > & get_bwd_map()
Definition: SF_container.h:374
index_mapping(const vector< T > &a, const vector< T > &b)
Constructor that uses assign() to set up the index mapping.
Definition: SF_container.h:225
size_t size() const
number of entries in bijective map
Definition: SF_container.h:346
void forward_map(vector< T > &idx) const
Map a whole array of indices in a to indices in b.
Definition: SF_container.h:281
hashmap::unordered_map< T, T > & get_fwd_map()
Definition: SF_container.h:369
const hashmap::unordered_map< T, T > & get_bwd_map() const
Definition: SF_container.h:364
T forward_map(T idx) const
Map one index from a to b.
Definition: SF_container.h:253
bool in_a(const T idx)
return whether idx is in set A
Definition: SF_container.h:351
bool in_b(const T idx)
return whether idx is in set B
Definition: SF_container.h:355
T backward_map(T idx) const
Map one index from b to a.
Definition: SF_container.h:263
index_mapping()
empty constructor.
Definition: SF_container.h:216
The mesh storage class. It contains both element and vertex data.
Definition: SF_container.h:385
void globalize(SF_nbr nbr_type)
Localize the connectivity data w.r.t. a given numbering.
Definition: SF_container.h:515
void clear_data()
Clear the mesh data from memory.
Definition: SF_container.h:542
overlapping_layout< T > pl
nodal parallel layout
Definition: SF_container.h:419
size_t g_numpts
global number of points
Definition: SF_container.h:389
meshdata()
construct empty mesh
Definition: SF_container.h:423
vector< T > dsp
connectivity starting index of each element
Definition: SF_container.h:405
vector< S > she
sheet direction
Definition: SF_container.h:410
vector< S > fib
fiber direction
Definition: SF_container.h:409
size_t l_numelem
local number of elements
Definition: SF_container.h:388
vector< elem_t > type
element type
Definition: SF_container.h:407
std::map< SF_nbr, vector< T > > nbr
container for different numberings
Definition: SF_container.h:415
std::string name
the mesh name
Definition: SF_container.h:396
void generate_par_layout()
Set up the parallel layout.
Definition: SF_container.h:528
vector< T > & register_numbering(SF_nbr nbr_type)
Register a new numbering to the mesh and return the associated index vector.
Definition: SF_container.h:434
vector< T > con
Definition: SF_container.h:401
vector< S > xyz_Euler
node coordinates in Euler formulation
Definition: SF_container.h:417
vector< S > xyz
node coordinates in Lagrange formulation
Definition: SF_container.h:416
const vector< T > & get_numbering(SF_nbr nbr_type) const
Get the vector defining a certain numbering.
Definition: SF_container.h:472
size_t l_numpts
local number of points
Definition: SF_container.h:390
size_t g_numelem
global number of elements
Definition: SF_container.h:387
void localize(SF_nbr nbr_type)
Localize the connectivity data w.r.t. a given numbering.
Definition: SF_container.h:486
MPI_Comm comm
the parallel mesh is defined on a MPI world
Definition: SF_container.h:393
vector< T > & get_numbering(SF_nbr nbr_type)
Get the vector defining a certain numbering.
Definition: SF_container.h:454
vector< T > tag
element tag
Definition: SF_container.h:406
hashmap::unordered_set< int > extr_tag
the element tags based on which the mesh has been extracted
Definition: SF_container.h:413
non_overlapping_layout< T > epl
element parallel layout
Definition: SF_container.h:420
The parallel layout of non overlapping indices.
The overlapping_layout class contains the algorithms related to managing overlapping parallel index s...
void assign(const vector< T > &idx, MPI_Comm comm)
Initialization function.
size_t size() const
The current size of the vector.
Definition: SF_vector.h:104
void resize(size_t n)
Resize a vector.
Definition: SF_vector.h:209
void assign(InputIterator s, InputIterator e)
Assign a memory range.
Definition: SF_vector.h:161
void reallocate()
Definition: SF_vector.h:257
iterator find(const K &key)
Search for key. Return iterator.
Definition: hashmap.hpp:593
hm_int count(const K &key) const
Check if key exists.
Definition: hashmap.hpp:579
void sort(Compare comp=Compare())
Sort data entries.
Definition: hashmap.hpp:659
size_t size() const
Definition: hashmap.hpp:687
Dense matrix class and associated funcs.
Classes similar to unordered_set and unordered_map, but with better performance.
unsigned long int hm_uint
Definition: hashmap.hpp:43
Definition: dense_mat.hpp:34
double mag(const Point &vect)
vector magnitude
Definition: SF_container.h:112
void cnt_from_dsp(const vector< T > &dsp, vector< T > &cnt)
Compute counts from displacements.
Definition: SF_vector.h:319
Point arr_to_point(V *arr)
Definition: SF_container.h:73
void dsp_from_cnt(const vector< T > &cnt, vector< T > &dsp)
Compute displacements from counts.
Definition: SF_vector.h:310
void transpose_connectivity(const vector< T > &a_cnt, const vector< T > &a_con, vector< T > &b_cnt, vector< T > &b_con)
Transpose CRS matrix graph A into B.
bool operator<(const struct tuple< T > &lhs, const struct tuple< T > &rhs)
Definition: SF_container.h:764
dmat< S > operator-(const dmat< S > &a, const dmat< S > &b)
Definition: dense_mat.hpp:421
void interval(vector< T > &vec, size_t start, size_t end)
Create an integer interval between start and end.
Definition: SF_vector.h:350
void get_alg_mask(const meshdata< T, S > &mesh, vector< bool > &alg_mask)
Definition: SF_container.h:614
double inner_prod(const Point &a, const Point &b)
Definition: SF_container.h:91
void binary_sort_copy(vector< T > &_V, vector< S > &_W)
Definition: SF_sort.h:296
T sum(const vector< T > &vec)
Compute sum of a vector's entries.
Definition: SF_vector.h:340
void count(const vector< T > &data, vector< S > &cnt)
Count number of occurrences of indices.
Definition: SF_vector.h:332
void point_to_arr(Point &p, V *arr)
Definition: SF_container.h:79
Point project(const Point &a, const Point &b)
project b onto a
Definition: SF_container.h:139
void outer_prod(const Point &a, const Point &b, const double s, double *buff, const bool add=false)
Definition: SF_container.h:96
void multiply_connectivities(const vector< T > &a_cnt, const vector< T > &a_con, const vector< T > &b_cnt, const vector< T > &b_con, vector< T > &c_cnt, vector< T > &c_con)
int max_nodal_edgecount(const meshdata< T, S > &mesh)
Compute the maximum number of node-to-node edges for a mesh.
Definition: SF_container.h:599
dmat< S > operator*(const dmat< S > &a, const dmat< S > &b)
Definition: dense_mat.hpp:373
Point normalize(const Point &vect)
Definition: SF_container.h:117
void array_to_points(const vector< S > &arr, vector< POINT > &pts)
Definition: SF_container.h:181
Point orthogonalize(const Point &a, const Point &b)
Definition: SF_container.h:144
void global_to_local(const vector< T > &glob, vector< T > &data, bool sortedData, bool doWarn)
Definition: SF_sort.h:540
elem_t getElemTypeID(char *eletype)
Generate element type enum from string.
Definition: SF_container.h:156
double distance(const Point &a, const Point &b)
Definition: SF_container.h:150
dmat< S > operator+(const dmat< S > &a, const dmat< S > &b)
Definition: dense_mat.hpp:413
Point cross(const Point &a, const Point &b)
cross product
Definition: SF_container.h:85
elem_t
element type enum
Definition: SF_container.h:54
@ Line
Definition: SF_container.h:62
@ Tri
Definition: SF_container.h:61
@ Prism
Definition: SF_container.h:59
@ Pyramid
Definition: SF_container.h:58
@ Tetra
Definition: SF_container.h:55
@ Quad
Definition: SF_container.h:60
@ Octa
Definition: SF_container.h:57
@ Hexa
Definition: SF_container.h:56
void nodal_connectivity_graph(const meshdata< T, S > &mesh, vector< T > &n2n_cnt, vector< T > &n2n_con)
Compute the node-to-node connectivity.
Definition: SF_container.h:574
SF_nbr
Enumeration encoding the different supported numberings.
Definition: SF_container.h:189
@ NBR_PETSC
PETSc numbering of nodes.
Definition: SF_container.h:192
@ NBR_ELEM_REF
The element numbering of the reference mesh (the one stored on HD).
Definition: SF_container.h:193
@ NBR_REF
The nodal numbering of the reference mesh (the one stored on HD).
Definition: SF_container.h:190
@ NBR_SUBMESH
Submesh nodal numbering: The globally ascending sorted reference indices are reindexed.
Definition: SF_container.h:191
@ NBR_ELEM_SUBMESH
Submesh element numbering: The globally ascending sorted reference indices are reindexed.
Definition: SF_container.h:194
const hm_uint mkhash_init
Definition: hashmap.hpp:57
hm_uint mkhash(hm_uint a, hm_uint b)
Definition: hashmap.hpp:52
Point and vector struct.
Definition: SF_container.h:66
double y
Definition: SF_container.h:68
double z
Definition: SF_container.h:69
double x
Definition: SF_container.h:67
static hm_uint hash(SF::quadruple< T > a)
Definition: SF_container.h:854
static bool cmp(SF::quadruple< T > a, SF::quadruple< T > b)
Definition: SF_container.h:849
static bool cmp(SF::triple< T > a, SF::triple< T > b)
Definition: SF_container.h:831
static hm_uint hash(SF::triple< T > a)
Definition: SF_container.h:836
static hm_uint hash(SF::tuple< T > a)
Definition: SF_container.h:819
static bool cmp(SF::tuple< T > a, SF::tuple< T > b)
Definition: SF_container.h:814
Base hashing class.
Definition: hashmap.hpp:88