Flow123d  release_3.0.0-1008-gca43bb7
ref_element.hh
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file ref_element.hh
15  * @brief Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc.
16  * @author Jan Stebel
17  * @todo
18  *
19  * TODO: reconsider following whether it is actual...
20  * - design interface in such a way, that we can change numbering
21  * - design numbering and orientations on ref element that is consistent (orientation and numbering od 2d el. match sides of 3d),
22  * and possibly allows permute vertices of elements so that sides sharing an edge match numbering and orientation (we dos'nt need permute faces)
23  *
24  * Proposal(prefers combinatoric order) :
25  * 1D - orientation V0 -> V1
26  *
27  * 2D - edges: E0: V0 -> V1,
28  * E1: V0 -> V2
29  * E2: V1 -> V2
30  * This maximize number of edge orientations achievable by edge permutations
31  * edge numbering edge orientation( in original numbering)
32  * 0 1 2 + + +
33  * 0 2 1 - + +
34  * 1 0 2 + + -
35  * 1 2 0 - - +
36  * 2 0 1 + - -
37  * 2 1 0 - - -
38  *
39  * vertices edges normal (out = +)
40  * 3D - sides: S0: 0 1 2 E0 E1 E3 -
41  * S1: 0 1 3 E0 E2 E4 +
42  * S2: 0 2 3 E1 E2 E5 -
43  * S3: 1 2 3 E3 E4 E5 -
44  *
45  * edges: E0: V0 -> V1 x direction
46  * E1: V0 -> V2 y direction
47  * E2: V0 -> V3 z direction
48  * E3: V1 -> V2
49  * E4: V1 -> V3
50  * E5: V2 -> V3
51  *
52  * - functions from DEAL.ii:
53  * bool is_inside_unit_cell( point )
54  * line_to_cell_vertices(line, vertex) vertex index on line to index on whole element
55  * face_to_cell_vertices(face, vertex, Orientation), Orientation should be some class describing permutation of shared face to element's face/side
56  * face_to_cel_lines
57  * standard_to_real_face_vertex(vertex, Orientation), maps vertex to permuted face
58  * real_to_standard_face_vertex - inverse
59  * ... same for line; we should need also something like standard_to_real_line_vertex; Deal dosn;t support Orientation changes for 2D element faces
60  * Point unit_cell_vertex(vertex) - coordinates
61  * project_to_unit_cell
62  * distance_to_unit_cell
63  * d_linear_shape_function
64  *
65  * - can not change numbering of element sides due to DarcyFlow, which use hardwired side numbering in construction of basis functions
66  * - any change of side numbering requires also change in flow/old_bcd.cc
67  *
68  *
69  */
70 
71 #ifndef REF_ELEMENT_HH_
72 #define REF_ELEMENT_HH_
73 
74 #include <vector> // for vector
75 #include <array>
76 #include <armadillo>
77 #include "system/asserts.hh"
78 
79 
80 /*
81  * Ordering of nodes and sides in reference elements
82  * =================================================
83  *
84  * TODO we want the following (22.10.):
85  *
86  * 1D element (line segment) 2D element (triangle) 3D element (tetrahedron)
87  *
88  * z
89  * .
90  * ,/
91  * /
92  * 3
93  * y ,/|`\
94  * ^ ,/ | `\
95  * | ,/ '. `\
96  * 2 ,/ | `\
97  * |`\ ,/ | `\
98  * | `\ 0-----------'.--------1 --> x
99  * | `\ `\. | ,/
100  * | `\ `\. | ,/
101  * | `\ `\. '. ,/
102  * 0----------1 --> x 0----------1 --> x `\. |/
103  * `2
104  * `\.
105  * `y
106  *
107  * side id node ids side id node ids side id node ids normal
108  * 0 0 0 0,1 0 0,1,2 OUT
109  * 1 1 1 0,2 1 0,1,3 IN
110  * 2 1,2 2 0,2,3 OUT
111  * 3 1,2,3 IN
112  *
113  *
114  * nodes coordinates:
115  * 0 [0] 0 [0,0] 0 [0,0,0]
116  * 1 [1] 1 [1,0] 1 [1,0,0]
117  * 2 [0,1] 2 [0,1,0]
118  * 3 [0,0,1]
119  *
120  * barycentric coordinates of nodes:
121  * 0 [1,0] 0 [1,0,0] 0 [1,0,0,0]
122  * 1 [0,1] 1 [0,1,0] 1 [0,1,0,0]
123  * 2 [0,0,1] 2 [0,0,1,0]
124  * 3 [0,0,0,1]
125  */
126 
127 /** Auxilliary class representing vector of indices (unsigned int).
128  * @tparam Size is the fixed size of the vector.
129  */
130 /*
131 template<unsigned int Size>
132 class IdxVector{
133  unsigned int data_[Size]; ///< Array with indices.
134 
135  public:
136  /// Constructor taking in array of indices.
137  IdxVector(std::array<unsigned int,Size> data_in);
138  /// Constructor enabling creating object with initializer list {...}.
139  IdxVector(std::initializer_list<unsigned int> data_in);
140  /// Getter for index @p idx.
141  unsigned int operator[](unsigned int idx) const;
142 };
143 */
144 
145 template<std::size_t Size>
146 using IdxVector = std::array<unsigned int, Size>;
147 
148 
149 /** Auxilliary structure that is used to pass template arguments into interact function of RefElement:
150  * RefElement<dim>::interact( Interaction<OutDim,InDim>(i) )
151  *
152  * This enables automatic deduction of dimensional template arguments.
153  * @see @p RefElement<dim>::interact
154  */
155 template <unsigned int OutDim, unsigned int InDim>
156 struct Interaction {
157  Interaction(unsigned int i) : i_(i) {}
158  unsigned int i_;
159 };
160 
161 template<unsigned int dim>
163 {
164 public:
165  typedef arma::vec::fixed<dim> LocalPoint;
166  /**
167  * Barycentric coordinates.
168  *
169  * e.g. coordinates (a,b,c) on triangle with vertices X, Y, Z
170  * represents a point: a*X+b*Y+c*Z
171  */
172  typedef arma::vec::fixed<dim+1> BaryPoint;
173  typedef arma::vec::fixed<dim> FaceBaryPoint;
174 
175  /**
176  * Return coordinates of given node.
177  * @see the class documentation @p RefElement
178  * @param nid Node number.
179  * NOTE: Implementation is dependent on current node and side numbering.
180  */
181  static LocalPoint node_coords(unsigned int nid);
182 
183  /**
184  * Compute normal vector to a given side.
185  * @param sid Side number.
186  */
187  static LocalPoint normal_vector(unsigned int sid);
188 
189 
190  /**
191  * If the given barycentric coordinate is in the ref. element, return unchanged.
192  * If the given barycentric coordinate is out of the ref. element,
193  * project it on the surface of the ref. element.
194  */
195  static BaryPoint clip(const BaryPoint &barycentric);
196 
197  /** Returns orientation of the normal of side @p sid. 0 -> OUT, 1 -> IN.
198  * NOTE: Implementation is dependent on current node and side numbering.
199  */
200  static unsigned int normal_orientation(unsigned int sid);
201 
202  static double side_measure(unsigned int sid);
203 
204  /**
205  * Returns index of the node that is oposite to side of given index @p sid.
206  * Note: It is dependent on current node and side numbering.
207  * @param sid Side number.
208  * NOTE: Implementation is dependent on current node and side numbering.
209  */
210  static unsigned int oposite_node(unsigned int sid);
211 
212  /**
213  * Return index of 1D line, shared by two faces @p f1 and @p f2 of the reference tetrahedron.
214  * Implemented only for @p dim == 3.
215  */
216  static unsigned int line_between_faces(unsigned int f1, unsigned int f2);
217 
218 
219  static const unsigned int n_sides = dim + 1; ///< Number of sides.
220  static const unsigned int n_nodes = dim + 1; ///< Number of nodes.
221  static const unsigned int n_nodes_per_side = dim; ///< Number of nodes on one side.
222  static const unsigned int n_lines_per_node = dim; ///< Number of lines with one common node.
223  static const unsigned int n_nodes_per_line = 2; ///< Number of nodes in one line.
224  static const unsigned int n_sides_per_line = 2; ///< Number of sides with one common line. @p dim == 3.
225  static const unsigned int n_sides_per_node = dim; ///< Number of sides with one common line.
226 
227  /// Number of lines on boundary of one side.
228  static const unsigned int n_lines_per_side = (unsigned int)((dim * (dim - 1)) / 2);//( dim == 3 ? 3 : 0);// Kombinační číslo dim nad dvěma
229 
230  /// Number of lines, i.e. @p object of dimension @p dim-2 on the boundary of the reference element.
231  static const unsigned int n_lines = (unsigned int)((dim * (dim + 1)) / 2); //( dim == 3 ? 6 : dim == 2 ? 3 : dim == 1 ? 1 : 0); součet posloupnosti
232 
233 
234 // /**
235 // * Node numbers for each side.
236 // */
237 // static const unsigned int side_nodes[n_sides][n_nodes_per_side];
238 //
239 // /**
240 // * Indices of 1D lines of the 2D sides of an tetrahedron. Nonempty only for @p dim==3.
241 // */
242 // static const unsigned int side_lines[n_sides][n_lines_per_side];
243 //
244 // /**
245 // * Nodes of 1D lines of the tetrahedron.
246 // */
247 // static const unsigned int line_nodes[n_lines][2];
248 //
249 // /**
250 // * Indices of sides for each line. Nonempty only for @p dim==3 and @p dim==2.
251 // */
252 // static const unsigned int line_sides[n_lines][2];
253 
254 
256 
257  /**
258  * Number of permutations of nodes on sides.
259  * dim value
260  * -----------
261  * 1 1
262  * 2 2
263  * 3 6
264  */
265  static const unsigned int n_side_permutations = (dim+1)*(2*dim*dim-5*dim+6)/6;
266 
267  /**
268  * Permutations of nodes on sides.
269  */
270  static const unsigned int side_permutations[n_side_permutations][n_nodes_per_side];
271 
272  /**
273  * For a given permutation @p p of nodes finds its index within @p side_permutations.
274  * @param p Permutation of nodes.
275  */
276  static unsigned int permutation_index(unsigned int p[n_nodes_per_side]);
277 
278  /** @brief Converts from local to barycentric coordinates.
279  * @param lp point in local coordinates (x,y)
280  * @return point in barycentric coordinates (1-x-y, x, y)
281  */
282  static BaryPoint local_to_bary(const LocalPoint& lp);
283 
284  /** @brief Converts from barycentric to local coordinates.
285  * @param bp point in barycentric coordinates
286  * @return point in local coordinates
287  */
288  static LocalPoint bary_to_local(const BaryPoint& bp);
289 
291 
292  /**
293  * Used in the clip method.
294  */
295  static BarycentricUnitVec make_bary_unit_vec();
296 
297  /**
298  * For given barycentric coordinates on the ref element returns barycentric coordinates
299  * on the ref. element of given face. Assumes that the input point is on the face.
300  * Barycentric order: (complanatory, local_coords )
301  */
302  static FaceBaryPoint barycentric_on_face(const BaryPoint &barycentric, unsigned int i_face);
303 
304 
306  static CentersList centers_of_subelements(unsigned int sub_dim);
307 
308  /**
309  * Return (1) number of zeros and (2) positions of zeros in barycentric coordinates.
310  * @p tolerance serves for testing zero values of @p barycentric coordinates.
311  */
312  static std::pair<unsigned int, unsigned int> zeros_positions(const BaryPoint &barycentric,
313  double tolerance = std::numeric_limits<double>::epsilon()*2);
314 
315  /**
316  * According to positions of zeros in barycentric coordinates, it gives the index of subdim-simplex
317  * in the reference element. Number of zeros must be equal to (3-subdim).
318  * e.g.:
319  * if 1 zeros, return index of side (subdim 2)
320  * if 2 zeros, return index of edge (subdim 1)
321  * if 3 zeros, return index of vertex (subdim 0)
322  */
323  template<unsigned int subdim> static unsigned int topology_idx(unsigned int zeros_positions);
324 
325  /** Function returns number of subdim-simplices inside dim-simplex.
326  * The aim is covering all the n_**** members with a single function.
327  * TODO: think of generalization for n_****_per_**** members, like function @p interact:
328  * template<unsigned int subdimA, unsigned int subdimB> static unsigned int count();
329  */
330  template<unsigned int subdim> static unsigned int count();
331 
332  /**
333  * @param sid - index of a sub-simplex in a simplex
334  * return an array of barycentric coordinates on <dim> simplex from <subdim> simplex
335  * for example: simplex<3> - ABCD and its subsubsimplex<1> AD (line index: 3)
336  * AD has barycoords for A (1,0), for D (0,1), but A in ABCD is (1,0,0,0) and D is (0,0,0,1)
337  * this method creates array ((1,0,0,0),(0,0,0,1))
338  */
339  template<unsigned int subdim> static arma::mat::fixed<dim+1,subdim+1> bary_coords(unsigned int sid);
340 
341  /** Interpolate barycentric coords to a higher dimension of a simplex.
342  * @param coord - barycentric coords of a point on a sub-simplex
343  * @param sub_simplex_idx - id of sub-simplex on a simplex
344  */
345  template<unsigned int subdim> static BaryPoint interpolate(arma::vec::fixed<subdim+1> coord, int sub_simplex_idx);
346 
347 
348  /**
349  * Basic line interpolation.
350  */
351  static BaryPoint line_barycentric_interpolation(BaryPoint first_coords,
352  BaryPoint second_coords,
353  double first_theta, double second_theta, double theta);
354 
355  /**
356  * Usage:
357  * RefElement<3>::interact(Interaction<2,0>(1))
358  * (means: In tetrahedron <3>, give indices of sides <2>, connected by node <0> with index 1)
359  * RefElement<3>::interact(Interaction<2,0>(1))[1]
360  * (as above, but give only the side with index 1)
361  *
362  * Template usage: RefElement<dim>::interact(Interaction<OutDim, InDim>(i))[j]
363  * (means: on dim-dimensional reference element, go on InDim-dimensional subelement with index i,
364  * which connects OutDim-dimnesional subelements and select the one with index j)
365  *
366  * This method serves as an interface to topology information of the reference element.
367  * It returns indices of OutDim-dimensional object
368  * of InDim-dimnesional object of given index
369  * in dim-dimnesional reference element.
370  * @tparam interaction - auxilliary object carying the index and the template arguments OutDim and InDim
371  * @tparam OutDim - output dimension (give me node-0, line-1, side-2), <= dim
372  * @tparam InDim - input dimension (for node-0, line-1, side-2), <= dim
373  * @return vector of indices of OutDim-dimensional subelements represented by @p IdxVector object.
374  *
375  * possible calls:
376  * dim OutDim InDim return
377  * 1,2,3 0 1 InDim+1 - give me indices of nodes of line of given index
378  * 3 0 2 InDim+1 - give me indices of nodes of a side (triangle) of given index
379  * 3 1 2 InDim+1 - give me indices of lines of side (triangle) of given index
380  *
381  * 1,2,3 1 0 dim-InDim - give me indices of lines with common node of given index
382  * 3 2 0 dim-InDim - give me indices of sides (triangles) with common node of given index
383  * 3 2 1 dim-InDim - give me indices of sides (triangles) with common line of given index
384  *
385  */
386  template < template <unsigned int OutDim, unsigned int InDim> class TInteraction, unsigned int OutDim, unsigned int InDim>
387  static const IdxVector< (InDim>OutDim ? InDim+1 : dim-InDim) > interact( TInteraction<OutDim,InDim> interaction );
388 
389 
390 private:
391  /// Internal part of the interact function.
392  template<unsigned int OutDim, unsigned int InDim>
393  static const IdxVector< (InDim>OutDim ? InDim+1 : dim-InDim) > interact_(unsigned int index);
394 
395  static const IdxVector<n_nodes_per_line> line_nodes_[n_lines]; ///< For given line, returns its nodes indices.
396  static const IdxVector<n_lines_per_node> node_lines_[n_nodes]; ///< For given node, returns lines indices.
397  static const IdxVector<n_nodes_per_side> side_nodes_[n_sides]; ///< For given side, returns nodes indices. For @p dim == 3.
398  static const IdxVector<n_sides_per_node> node_sides_[n_nodes]; ///< For given node, returns sides indices. For @p dim == 3.
399  static const IdxVector<n_sides_per_line> line_sides_[n_lines]; ///< For given line, returns sides indices. For @p dim == 3.
400  static const IdxVector<n_lines_per_side> side_lines_[n_sides]; ///< For given side, returns lines indices. For @p dim == 3.
401 
402  //TODO: implement for 1d and 2d
403  /**
404  * Consider an n-face (node, edge, face, bulk) with dimension `subdim` and
405  * index within subdimension `idx`. Barycentric coordinates of all points
406  * on the n-face have unique pattern of zero coordinates.
407  *
408  * topology_zeros_[subdim][idx] is a bitfield with '1' where the pattern have zeros.
409  */
410  static const IdxVector<(n_lines > n_nodes) ? n_lines : n_nodes> topology_zeros_[dim+1];
411 };
412 
413 
414 
415 
416 
417 
421 
425 
428 
430 
432 
433 template<> const unsigned int RefElement<1>::side_permutations[][n_nodes_per_side];
434 template<> const unsigned int RefElement<2>::side_permutations[][n_nodes_per_side];
435 template<> const unsigned int RefElement<3>::side_permutations[][n_nodes_per_side];
436 
440 
441 
442 
443 
444 
445 
446 /************************* template implementation ****************************/
447 
448 template<unsigned int dim>
449 template<unsigned int subdim> inline
450 arma::mat::fixed<dim+1,subdim+1> RefElement<dim>::bary_coords(unsigned int sid){
451  ASSERT_LT_DBG(subdim, dim).error("Dimension mismatch!");
452  arma::mat::fixed<dim+1,subdim+1> bary_c;
453 
454  for(unsigned int i = 0; i < subdim+1; i++){
455  unsigned int nid = interact_<0,subdim>(sid)[i];
456  bary_c.col(i).zeros();
457  bary_c.col(i)(nid) = 1;
458  }
459 
460  return bary_c;
461 };
462 
463 
464 template<unsigned int dim> inline
465 arma::vec::fixed<dim> RefElement<dim>::node_coords(unsigned int nid)
466 {
467  ASSERT_LT_DBG(nid, n_nodes).error("Node number is out of range!");
468 
469  arma::vec::fixed<dim> p;
470  p.zeros();
471 
472  if (nid > 0)
473  p(nid-1) = 1;
474 
475  return p;
476 }
477 
478 
479 template<unsigned int dim>
480 template<unsigned int subdim>
481 auto RefElement<dim>::interpolate(arma::vec::fixed<subdim+1> coord, int sub_simplex_idx) -> BaryPoint
482 {
483  return RefElement<dim>::bary_coords<subdim>(sub_simplex_idx)*coord;
484 };
485 /*
486 template <unsigned int Size>
487 IdxVector<Size>::IdxVector(std::array<unsigned int,Size> data_in)
488 : data_(data_in){}
489 
490 template <unsigned int Size>
491 IdxVector<Size>::IdxVector(std::initializer_list<unsigned int> data_in)
492 {
493  ASSERT_EQ_DBG(data_in.size(), Size).error("Incorrect data size.");
494  std::copy(data_in.begin(), data_in.end(), data_);
495 }
496 
497 template <unsigned int Size>
498 inline unsigned int IdxVector<Size>::operator[](unsigned int idx) const
499 { ASSERT_LT_DBG(idx, Size).error("Index out of bounds.");
500  return data_[idx]; }
501 
502 */
503 
504 template<> template<> inline unsigned int RefElement<3>::count<0>()
505 { return n_nodes; }
506 template<> template<> inline unsigned int RefElement<3>::count<1>()
507 { return n_lines; }
508 template<> template<> inline unsigned int RefElement<3>::count<2>()
509 { return n_sides; }
510 template<> template<> inline unsigned int RefElement<3>::count<3>()
511 { return 1; }
512 template<> template<> inline unsigned int RefElement<2>::count<0>()
513 { return n_nodes; }
514 template<> template<> inline unsigned int RefElement<2>::count<1>()
515 { return n_lines; }
516 template<> template<> inline unsigned int RefElement<2>::count<2>()
517 { return 1; }
518 template<> template<> inline unsigned int RefElement<2>::count<3>()
519 { return 0; }
520 template<> template<> inline unsigned int RefElement<1>::count<0>()
521 { return n_nodes; }
522 template<> template<> inline unsigned int RefElement<1>::count<1>()
523 { return 1; }
524 template<> template<> inline unsigned int RefElement<1>::count<2>()
525 { return 0; }
526 template<> template<> inline unsigned int RefElement<1>::count<3>()
527 { return 0; }
528 
529 template<unsigned int dim>
530 template<unsigned int subdim>
531 unsigned int RefElement<dim>::topology_idx(unsigned int zeros_positions)
532 {
533  for(unsigned int i=0; i < RefElement<dim>::count<subdim>(); i++){
534  if(zeros_positions == topology_zeros_[subdim][i]) return i;
535  }
536  ASSERT(0).error("Undefined zero pattern.");
537  return -1;
538 }
539 
540 
541 /// This function is for "side_nodes" - for given side, give me nodes (0->0, 1->1).
542 template<> template<> inline const IdxVector<1> RefElement<1>::interact_<0,0>(unsigned int i)
543 { ASSERT_LT_DBG(i, RefElement<1>::n_nodes).error("Index out of bounds.");
544  return IdxVector<1>({i});}
545 
546 /// For line i {0}, give me indices of its nodes.
547 template<> template<> inline const IdxVector<2> RefElement<1>::interact_<0,1>(unsigned int i)
548 { ASSERT_LT_DBG(i, RefElement<1>::n_lines).error("Index out of bounds.");
549  return line_nodes_[i];}
550 
551 /// For line i {0,1,2}, give me indices of its nodes.
552 template<> template<> inline const IdxVector<2> RefElement<2>::interact_<0,1>(unsigned int i)
553 { ASSERT_LT_DBG(i, RefElement<2>::n_lines).error("Index out of bounds.");
554  return line_nodes_[i];}
555 
556 /// For line i {0,1,2,3,4,5}, give me indices of its nodes.
557 template<> template<> inline const IdxVector<2> RefElement<3>::interact_<0,1>(unsigned int i)
558 { ASSERT_LT_DBG(i, RefElement<3>::n_lines).error("Index out of bounds.");
559  return line_nodes_[i];}
560 
561 /// For node i {0,1}, give me indices of lines.
562 template<> template<> inline const IdxVector<1> RefElement<1>::interact_<1,0>(unsigned int i)
563 { ASSERT_LT_DBG(i, RefElement<1>::n_nodes).error("Index out of bounds.");
564  return node_lines_[i];}
565 
566 /// For node i {0,1,2}, give me indices of lines.
567 template<> template<> inline const IdxVector<2> RefElement<2>::interact_<1,0>(unsigned int i)
568 { ASSERT_LT_DBG(i, RefElement<2>::n_nodes).error("Index out of bounds.");
569  return node_lines_[i];}
570 
571 /// For node i {0,1,2,3}, give me indices of lines.
572 template<> template<> inline const IdxVector<3> RefElement<3>::interact_<1,0>(unsigned int i)
573 { ASSERT_LT_DBG(i, RefElement<3>::n_nodes).error("Index out of bounds.");
574  return node_lines_[i];}
575 
576 /// For side i {0,1,2}, give me indices of its nodes.
577 template<> template<> inline const IdxVector<3> RefElement<3>::interact_<0,2>(unsigned int i)
578 { ASSERT_LT_DBG(i, RefElement<3>::n_sides).error("Index out of bounds.");
579  return side_nodes_[i];}
580 
581 /// For node i {0,1,2,3}, give me indices of sides.
582 template<> template<> inline const IdxVector<3> RefElement<3>::interact_<2,0>(unsigned int i)
583 { ASSERT_LT_DBG(i, RefElement<3>::n_sides).error("Index out of bounds.");
584  return node_sides_[i];}
585 
586 /// For line i {0,1,2,3}, give me indices of sides.
587 template<> template<> inline const IdxVector<2> RefElement<3>::interact_<2,1>(unsigned int i)
588 { ASSERT_LT_DBG(i, RefElement<3>::n_lines).error("Index out of bounds.");
589  return line_sides_[i];}
590 
591 /// For side i {0,1,2}, give me indices of its lines.
592 template<> template<> inline const IdxVector<3> RefElement<3>::interact_<1,2>(unsigned int i)
593 { ASSERT_LT_DBG(i, RefElement<3>::n_sides).error("Index out of bounds.");
594  return side_lines_[i];}
595 
596 template<unsigned int dim> template<unsigned int OutDim, unsigned int InDim>
597 inline const IdxVector< (InDim>OutDim ? InDim+1 : dim-InDim) > RefElement<dim>::interact_(unsigned int i)
598 {
599  ASSERT(false)(dim)(OutDim)(InDim)(i).error("Not implemented.");
600  //ASSERT_LT_DBG(OutDim, dim);
601  //ASSERT_LT_DBG(InDim, dim);
602 }
603 
604 
605 template<unsigned int dim>
606 template < template <unsigned int OutDim, unsigned int InDim> class TInteraction, unsigned int OutDim, unsigned int InDim>
607 inline const IdxVector< (InDim>OutDim ? InDim+1 : dim-InDim) > RefElement<dim>::interact( TInteraction<OutDim,InDim> interaction )
608 {
609  return interact_<OutDim,InDim>(interaction.i_);
610 }
611 
612 #endif /* REF_ELEMENT_HH_ */
static arma::mat::fixed< dim+1, subdim+1 > bary_coords(unsigned int sid)
Definition: ref_element.hh:450
static BaryPoint interpolate(arma::vec::fixed< subdim+1 > coord, int sub_simplex_idx)
std::vector< BaryPoint > BarycentricUnitVec
Definition: ref_element.hh:290
arma::vec::fixed< dim > LocalPoint
Definition: ref_element.hh:165
Definitions of ASSERTS.
unsigned int i_
Definition: ref_element.hh:158
static const std::vector< std::vector< std::vector< unsigned int > > > nodes_of_subelements
Definition: ref_element.hh:255
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Definition: asserts.hh:346
static LocalPoint node_coords(unsigned int nid)
Definition: ref_element.hh:465
Interaction(unsigned int i)
Definition: ref_element.hh:157
const std::vector< LocalPoint > & CentersList
Definition: ref_element.hh:305
std::array< unsigned int, Size > IdxVector
Definition: ref_element.hh:146
arma::vec::fixed< dim > FaceBaryPoint
Definition: ref_element.hh:173
arma::vec::fixed< dim+1 > BaryPoint
Definition: ref_element.hh:172
static unsigned int topology_idx(unsigned int zeros_positions)
Definition: ref_element.hh:531
#define ASSERT_LT_DBG(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Definition: asserts.hh:299