Flow123d  release_2.2.0-914-gf1a3a4f
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 <armadillo>
75 #include "system/system.hh"
76 
77 
78 /*
79  * Ordering of nodes and sides in reference elements
80  * =================================================
81  *
82  * TODO we want the following (22.10.):
83  *
84  * 1D element (line segment) 2D element (triangle) 3D element (tetrahedron)
85  *
86  * z
87  * .
88  * ,/
89  * /
90  * 3
91  * y ,/|`\
92  * ^ ,/ | `\
93  * | ,/ '. `\
94  * 2 ,/ | `\
95  * |`\ ,/ | `\
96  * | `\ 0-----------'.--------1 --> x
97  * | `\ `\. | ,/
98  * | `\ `\. | ,/
99  * | `\ `\. '. ,/
100  * 0----------1 --> x 0----------1 --> x `\. |/
101  * `2
102  * `\.
103  * `y
104  *
105  * side id node ids side id node ids side id node ids normal
106  * 0 0 0 0,1 0 0,1,2 OUT
107  * 1 1 1 0,2 1 0,1,3 IN
108  * 2 1,2 2 0,2,3 OUT
109  * 3 1,2,3 IN
110  *
111  *
112  * nodes coordinates:
113  * 0 [0] 0 [0,0] 0 [0,0,0]
114  * 1 [1] 1 [1,0] 1 [1,0,0]
115  * 2 [0,1] 2 [0,1,0]
116  * 3 [0,0,1]
117  *
118  * barycentric coordinates of nodes:
119  * 0 [1,0] 0 [1,0,0] 0 [1,0,0,0]
120  * 1 [0,1] 1 [0,1,0] 1 [0,1,0,0]
121  * 2 [0,0,1] 2 [0,0,1,0]
122  * 3 [0,0,0,1]
123  */
124 
125 /** Auxilliary class representing vector of indices (unsigned int).
126  * @tparam Size is the fixed size of the vector.
127  */
128 /*
129 template<unsigned int Size>
130 class IdxVector{
131  unsigned int data_[Size]; ///< Array with indices.
132 
133  public:
134  /// Constructor taking in array of indices.
135  IdxVector(std::array<unsigned int,Size> data_in);
136  /// Constructor enabling creating object with initializer list {...}.
137  IdxVector(std::initializer_list<unsigned int> data_in);
138  /// Getter for index @p idx.
139  unsigned int operator[](unsigned int idx) const;
140 };
141 */
142 
143 template<std::size_t Size>
144 using IdxVector = std::array<unsigned int, Size>;
145 
146 
147 /** Auxilliary structure that is used to pass template arguments into interact function of RefElement:
148  * RefElement<dim>::interact( Interaction<OutDim,InDim>(i) )
149  *
150  * This enables automatic deduction of dimensional template arguments.
151  * @see @p RefElement<dim>::interact
152  */
153 template <unsigned int OutDim, unsigned int InDim>
154 struct Interaction {
155  Interaction(unsigned int i) : i_(i) {}
156  unsigned int i_;
157 };
158 
159 template<unsigned int dim>
161 {
162 public:
163  typedef arma::vec::fixed<dim> LocalPoint;
164  /**
165  * Barycentric coordinates.
166  *
167  * e.g. coordinates (a,b,c) on triangle with vertices X, Y, Z
168  * represents a point: a*X+b*Y+c*Z
169  */
170  typedef arma::vec::fixed<dim+1> BaryPoint;
171  typedef arma::vec::fixed<dim> FaceBaryPoint;
172 
173  /**
174  * Return coordinates of given node.
175  * @see the class documentation @p RefElement
176  * @param nid Node number.
177  * NOTE: Implementation is dependent on current node and side numbering.
178  */
179  static LocalPoint node_coords(unsigned int nid);
180 
181  /**
182  * Compute normal vector to a given side.
183  * @param sid Side number.
184  */
185  static LocalPoint normal_vector(unsigned int sid);
186 
187 
188  /**
189  * If the given barycentric coordinate is in the ref. element, return unchanged.
190  * If the given barycentric coordinate is out of the ref. element,
191  * project it on the surface of the ref. element.
192  */
193  static BaryPoint clip(const BaryPoint &barycentric);
194 
195  /** Returns orientation of the normal of side @p sid. 0 -> OUT, 1 -> IN.
196  * NOTE: Implementation is dependent on current node and side numbering.
197  */
198  static unsigned int normal_orientation(unsigned int sid);
199 
200  static double side_measure(unsigned int sid);
201 
202  /**
203  * Returns index of the node that is oposite to side of given index @p sid.
204  * Note: It is dependent on current node and side numbering.
205  * @param sid Side number.
206  * NOTE: Implementation is dependent on current node and side numbering.
207  */
208  static unsigned int oposite_node(unsigned int sid);
209 
210  /**
211  * Return index of 1D line, shared by two faces @p f1 and @p f2 of the reference tetrahedron.
212  * Implemented only for @p dim == 3.
213  */
214  static unsigned int line_between_faces(unsigned int f1, unsigned int f2);
215 
216 
217  static const unsigned int n_sides = dim + 1; ///< Number of sides.
218  static const unsigned int n_nodes = dim + 1; ///< Number of nodes.
219  static const unsigned int n_nodes_per_side = dim; ///< Number of nodes on one side.
220  static const unsigned int n_lines_per_node = dim; ///< Number of lines with one common node.
221  static const unsigned int n_nodes_per_line = 2; ///< Number of nodes in one line.
222  static const unsigned int n_sides_per_line = 2; ///< Number of sides with one common line. @p dim == 3.
223  static const unsigned int n_sides_per_node = dim; ///< Number of sides with one common line.
224 
225  /// Number of lines on boundary of one side.
226  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
227 
228  /// Number of lines, i.e. @p object of dimension @p dim-2 on the boundary of the reference element.
229  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
230 
231 
232 // /**
233 // * Node numbers for each side.
234 // */
235 // static const unsigned int side_nodes[n_sides][n_nodes_per_side];
236 //
237 // /**
238 // * Indices of 1D lines of the 2D sides of an tetrahedron. Nonempty only for @p dim==3.
239 // */
240 // static const unsigned int side_lines[n_sides][n_lines_per_side];
241 //
242 // /**
243 // * Nodes of 1D lines of the tetrahedron.
244 // */
245 // static const unsigned int line_nodes[n_lines][2];
246 //
247 // /**
248 // * Indices of sides for each line. Nonempty only for @p dim==3 and @p dim==2.
249 // */
250 // static const unsigned int line_sides[n_lines][2];
251 
252 
254 
255  /**
256  * Number of permutations of nodes on sides.
257  * dim value
258  * -----------
259  * 1 1
260  * 2 2
261  * 3 6
262  */
263  static const unsigned int n_side_permutations = (dim+1)*(2*dim*dim-5*dim+6)/6;
264 
265  /**
266  * Permutations of nodes on sides.
267  */
268  static const unsigned int side_permutations[n_side_permutations][n_nodes_per_side];
269 
270  /**
271  * For a given permutation @p p of nodes finds its index within @p side_permutations.
272  * @param p Permutation of nodes.
273  */
274  static unsigned int permutation_index(unsigned int p[n_nodes_per_side]);
275 
276  /** @brief Converts from local to barycentric coordinates.
277  * @param lp point in local coordinates (x,y)
278  * @return point in barycentric coordinates (1-x-y, x, y)
279  */
280  static BaryPoint local_to_bary(const LocalPoint& lp);
281 
282  /** @brief Converts from barycentric to local coordinates.
283  * @param bp point in barycentric coordinates
284  * @return point in local coordinates
285  */
286  static LocalPoint bary_to_local(const BaryPoint& bp);
287 
289 
290  /**
291  * Used in the clip method.
292  */
293  static BarycentricUnitVec make_bary_unit_vec();
294 
295  /**
296  * For given barycentric coordinates on the ref element returns barycentric coordinates
297  * on the ref. element of given face. Assumes that the input point is on the face.
298  * Barycentric order: (complanatory, local_coords )
299  */
300  static FaceBaryPoint barycentric_on_face(const BaryPoint &barycentric, unsigned int i_face);
301 
302 
304  static CentersList centers_of_subelements(unsigned int sub_dim);
305 
306  /**
307  * Return (1) number of zeros and (2) positions of zeros in barycentric coordinates.
308  * @p tolerance serves for testing zero values of @p barycentric coordinates.
309  */
310  static std::pair<unsigned int, unsigned int> zeros_positions(const BaryPoint &barycentric,
312 
313  /**
314  * According to positions of zeros in barycentric coordinates, it gives the index of subdim-simplex
315  * in the reference element. Number of zeros must be equal to (3-subdim).
316  * e.g.:
317  * if 1 zeros, return index of side (subdim 2)
318  * if 2 zeros, return index of edge (subdim 1)
319  * if 3 zeros, return index of vertex (subdim 0)
320  */
321  template<unsigned int subdim> static unsigned int topology_idx(unsigned int zeros_positions);
322 
323  /** Function returns number of subdim-simplices inside dim-simplex.
324  * The aim is covering all the n_**** members with a single function.
325  * TODO: think of generalization for n_****_per_**** members, like function @p interact:
326  * template<unsigned int subdimA, unsigned int subdimB> static unsigned int count();
327  */
328  template<unsigned int subdim> static unsigned int count();
329 
330  /**
331  * @param sid - index of a sub-simplex in a simplex
332  * return an array of barycentric coordinates on <dim> simplex from <subdim> simplex
333  * for example: simplex<3> - ABCD and its subsubsimplex<1> AD (line index: 3)
334  * 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)
335  * this method creates array ((1,0,0,0),(0,0,0,1))
336  */
337  template<unsigned int subdim> static arma::mat::fixed<dim+1,subdim+1> bary_coords(unsigned int sid);
338 
339  /** Interpolate barycentric coords to a higher dimension of a simplex.
340  * @param coord - barycentric coords of a point on a sub-simplex
341  * @param sub_simplex_idx - id of sub-simplex on a simplex
342  */
343  template<unsigned int subdim> static BaryPoint interpolate(arma::vec::fixed<subdim+1> coord, int sub_simplex_idx);
344 
345 
346  /**
347  * Basic line interpolation.
348  */
349  static BaryPoint line_barycentric_interpolation(BaryPoint first_coords,
350  BaryPoint second_coords,
351  double first_theta, double second_theta, double theta);
352 
353  /**
354  * Usage:
355  * RefElement<3>::interact(Interaction<2,0>(1))
356  * (means: In tetrahedron <3>, give indices of sides <2>, connected by node <0> with index 1)
357  * RefElement<3>::interact(Interaction<2,0>(1))[1]
358  * (as above, but give only the side with index 1)
359  *
360  * Template usage: RefElement<dim>::interact(Interaction<OutDim, InDim>(i))[j]
361  * (means: on dim-dimensional reference element, go on InDim-dimensional subelement with index i,
362  * which connects OutDim-dimnesional subelements and select the one with index j)
363  *
364  * This method serves as an interface to topology information of the reference element.
365  * It returns indices of OutDim-dimensional object
366  * of InDim-dimnesional object of given index
367  * in dim-dimnesional reference element.
368  * @tparam interaction - auxilliary object carying the index and the template arguments OutDim and InDim
369  * @tparam OutDim - output dimension (give me node-0, line-1, side-2), <= dim
370  * @tparam InDim - input dimension (for node-0, line-1, side-2), <= dim
371  * @return vector of indices of OutDim-dimensional subelements represented by @p IdxVector object.
372  *
373  * possible calls:
374  * dim OutDim InDim return
375  * 1,2,3 0 1 InDim+1 - give me indices of nodes of line of given index
376  * 3 0 2 InDim+1 - give me indices of nodes of a side (triangle) of given index
377  * 3 1 2 InDim+1 - give me indices of lines of side (triangle) of given index
378  *
379  * 1,2,3 1 0 dim-InDim - give me indices of lines with common node of given index
380  * 3 2 0 dim-InDim - give me indices of sides (triangles) with common node of given index
381  * 3 2 1 dim-InDim - give me indices of sides (triangles) with common line of given index
382  *
383  */
384  template < template <unsigned int OutDim, unsigned int InDim> class TInteraction, unsigned int OutDim, unsigned int InDim>
385  static const IdxVector< (InDim>OutDim ? InDim+1 : dim-InDim) > interact( TInteraction<OutDim,InDim> interaction );
386 
387 
388 private:
389  /// Internal part of the interact function.
390  template<unsigned int OutDim, unsigned int InDim>
391  static const IdxVector< (InDim>OutDim ? InDim+1 : dim-InDim) > interact_(unsigned int index);
392 
393  static const IdxVector<n_nodes_per_line> line_nodes_[n_lines]; ///< For given line, returns its nodes indices.
394  static const IdxVector<n_lines_per_node> node_lines_[n_nodes]; ///< For given node, returns lines indices.
395  static const IdxVector<n_nodes_per_side> side_nodes_[n_sides]; ///< For given side, returns nodes indices. For @p dim == 3.
396  static const IdxVector<n_sides_per_node> node_sides_[n_nodes]; ///< For given node, returns sides indices. For @p dim == 3.
397  static const IdxVector<n_sides_per_line> line_sides_[n_lines]; ///< For given line, returns sides indices. For @p dim == 3.
398  static const IdxVector<n_lines_per_side> side_lines_[n_sides]; ///< For given side, returns lines indices. For @p dim == 3.
399 
400  //TODO: implement for 1d and 2d
401  static const IdxVector<(n_lines > n_nodes) ? n_lines : n_nodes> topology_zeros_[dim+1]; ///< Maps the zero mask of the barycentric coordinates to topology indices.
402 };
403 
404 
405 
406 
407 
408 
412 
416 
419 
421 
423 
424 template<> const unsigned int RefElement<1>::side_permutations[][n_nodes_per_side];
425 template<> const unsigned int RefElement<2>::side_permutations[][n_nodes_per_side];
426 template<> const unsigned int RefElement<3>::side_permutations[][n_nodes_per_side];
427 
430 
431 
432 
433 
434 
435 
436 /************************* template implementation ****************************/
437 
438 template<unsigned int dim>
439 template<unsigned int subdim> inline
440 arma::mat::fixed<dim+1,subdim+1> RefElement<dim>::bary_coords(unsigned int sid){
441  ASSERT_LT_DBG(subdim, dim).error("Dimension mismatch!");
442  arma::mat::fixed<dim+1,subdim+1> bary_c;
443 
444  for(unsigned int i = 0; i < subdim+1; i++){
445  unsigned int nid = interact_<0,subdim>(sid)[i];
446  bary_c.col(i).zeros();
447  bary_c.col(i)(nid) = 1;
448  }
449 
450  return bary_c;
451 };
452 
453 
454 template<unsigned int dim> inline
455 arma::vec::fixed<dim> RefElement<dim>::node_coords(unsigned int nid)
456 {
457  ASSERT_LT_DBG(nid, n_nodes).error("Node number is out of range!");
458 
459  arma::vec::fixed<dim> p;
460  p.zeros();
461 
462  if (nid > 0)
463  p(nid-1) = 1;
464 
465  return p;
466 }
467 
468 
469 template<unsigned int dim>
470 template<unsigned int subdim>
471 auto RefElement<dim>::interpolate(arma::vec::fixed<subdim+1> coord, int sub_simplex_idx) -> BaryPoint
472 {
473  return RefElement<dim>::bary_coords<subdim>(sub_simplex_idx)*coord;
474 };
475 /*
476 template <unsigned int Size>
477 IdxVector<Size>::IdxVector(std::array<unsigned int,Size> data_in)
478 : data_(data_in){}
479 
480 template <unsigned int Size>
481 IdxVector<Size>::IdxVector(std::initializer_list<unsigned int> data_in)
482 {
483  ASSERT_EQ_DBG(data_in.size(), Size).error("Incorrect data size.");
484  std::copy(data_in.begin(), data_in.end(), data_);
485 }
486 
487 template <unsigned int Size>
488 inline unsigned int IdxVector<Size>::operator[](unsigned int idx) const
489 { ASSERT_LT_DBG(idx, Size).error("Index out of bounds.");
490  return data_[idx]; }
491 
492 */
493 
494 template<> template<> inline unsigned int RefElement<3>::count<0>()
495 { return n_nodes; }
496 template<> template<> inline unsigned int RefElement<3>::count<1>()
497 { return n_lines; }
498 template<> template<> inline unsigned int RefElement<3>::count<2>()
499 { return n_sides; }
500 template<> template<> inline unsigned int RefElement<3>::count<3>()
501 { return 1; }
502 template<> template<> inline unsigned int RefElement<2>::count<0>()
503 { return n_nodes; }
504 template<> template<> inline unsigned int RefElement<2>::count<1>()
505 { return n_lines; }
506 template<> template<> inline unsigned int RefElement<2>::count<2>()
507 { return 1; }
508 template<> template<> inline unsigned int RefElement<2>::count<3>()
509 { return 0; }
510 template<> template<> inline unsigned int RefElement<1>::count<0>()
511 { return n_nodes; }
512 template<> template<> inline unsigned int RefElement<1>::count<1>()
513 { return 1; }
514 template<> template<> inline unsigned int RefElement<1>::count<2>()
515 { return 0; }
516 template<> template<> inline unsigned int RefElement<1>::count<3>()
517 { return 0; }
518 
519 template<unsigned int dim>
520 template<unsigned int subdim>
521 unsigned int RefElement<dim>::topology_idx(unsigned int zeros_positions)
522 {
523  for(unsigned int i=0; i < RefElement<dim>::count<subdim>(); i++){
524  if(zeros_positions == topology_zeros_[subdim][i]) return i;
525  }
526  ASSERT(0).error("Undefined zero pattern.");
527  return -1;
528 }
529 
530 
531 /// This function is for "side_nodes" - for given side, give me nodes (0->0, 1->1).
532 template<> template<> inline const IdxVector<1> RefElement<1>::interact_<0,0>(unsigned int i)
533 { ASSERT_LT_DBG(i, RefElement<1>::n_nodes).error("Index out of bounds.");
534  return IdxVector<1>({i});}
535 
536 /// For line i {0}, give me indices of its nodes.
537 template<> template<> inline const IdxVector<2> RefElement<1>::interact_<0,1>(unsigned int i)
538 { ASSERT_LT_DBG(i, RefElement<1>::n_lines).error("Index out of bounds.");
539  return line_nodes_[i];}
540 
541 /// For line i {0,1,2}, give me indices of its nodes.
542 template<> template<> inline const IdxVector<2> RefElement<2>::interact_<0,1>(unsigned int i)
543 { ASSERT_LT_DBG(i, RefElement<2>::n_lines).error("Index out of bounds.");
544  return line_nodes_[i];}
545 
546 /// For line i {0,1,2,3,4,5}, give me indices of its nodes.
547 template<> template<> inline const IdxVector<2> RefElement<3>::interact_<0,1>(unsigned int i)
548 { ASSERT_LT_DBG(i, RefElement<3>::n_lines).error("Index out of bounds.");
549  return line_nodes_[i];}
550 
551 /// For node i {0,1}, give me indices of lines.
552 template<> template<> inline const IdxVector<1> RefElement<1>::interact_<1,0>(unsigned int i)
553 { ASSERT_LT_DBG(i, RefElement<1>::n_nodes).error("Index out of bounds.");
554  return node_lines_[i];}
555 
556 /// For node i {0,1,2}, give me indices of lines.
557 template<> template<> inline const IdxVector<2> RefElement<2>::interact_<1,0>(unsigned int i)
558 { ASSERT_LT_DBG(i, RefElement<2>::n_nodes).error("Index out of bounds.");
559  return node_lines_[i];}
560 
561 /// For node i {0,1,2,3}, give me indices of lines.
562 template<> template<> inline const IdxVector<3> RefElement<3>::interact_<1,0>(unsigned int i)
563 { ASSERT_LT_DBG(i, RefElement<3>::n_nodes).error("Index out of bounds.");
564  return node_lines_[i];}
565 
566 /// For side i {0,1,2}, give me indices of its nodes.
567 template<> template<> inline const IdxVector<3> RefElement<3>::interact_<0,2>(unsigned int i)
568 { ASSERT_LT_DBG(i, RefElement<3>::n_sides).error("Index out of bounds.");
569  return side_nodes_[i];}
570 
571 /// For node i {0,1,2,3}, give me indices of sides.
572 template<> template<> inline const IdxVector<3> RefElement<3>::interact_<2,0>(unsigned int i)
573 { ASSERT_LT_DBG(i, RefElement<3>::n_sides).error("Index out of bounds.");
574  return node_sides_[i];}
575 
576 /// For line i {0,1,2,3}, give me indices of sides.
577 template<> template<> inline const IdxVector<2> RefElement<3>::interact_<2,1>(unsigned int i)
578 { ASSERT_LT_DBG(i, RefElement<3>::n_lines).error("Index out of bounds.");
579  return line_sides_[i];}
580 
581 /// For side i {0,1,2}, give me indices of its lines.
582 template<> template<> inline const IdxVector<3> RefElement<3>::interact_<1,2>(unsigned int i)
583 { ASSERT_LT_DBG(i, RefElement<3>::n_sides).error("Index out of bounds.");
584  return side_lines_[i];}
585 
586 template<unsigned int dim> template<unsigned int OutDim, unsigned int InDim>
587 inline const IdxVector< (InDim>OutDim ? InDim+1 : dim-InDim) > RefElement<dim>::interact_(unsigned int i)
588 {
589  ASSERT(false)(dim)(OutDim)(InDim)(i).error("Not implemented.");
590  //ASSERT_LT_DBG(OutDim, dim);
591  //ASSERT_LT_DBG(InDim, dim);
592 }
593 
594 
595 template<unsigned int dim>
596 template < template <unsigned int OutDim, unsigned int InDim> class TInteraction, unsigned int OutDim, unsigned int InDim>
597 inline const IdxVector< (InDim>OutDim ? InDim+1 : dim-InDim) > RefElement<dim>::interact( TInteraction<OutDim,InDim> interaction )
598 {
599  return interact_<OutDim,InDim>(interaction.i_);
600 }
601 
602 #endif /* REF_ELEMENT_HH_ */
static arma::mat::fixed< dim+1, subdim+1 > bary_coords(unsigned int sid)
Definition: ref_element.hh:440
static BaryPoint interpolate(arma::vec::fixed< subdim+1 > coord, int sub_simplex_idx)
std::vector< BaryPoint > BarycentricUnitVec
Definition: ref_element.hh:288
arma::vec::fixed< dim > LocalPoint
Definition: ref_element.hh:163
unsigned int i_
Definition: ref_element.hh:156
static const std::vector< std::vector< std::vector< unsigned int > > > nodes_of_subelements
Definition: ref_element.hh:253
#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:455
Interaction(unsigned int i)
Definition: ref_element.hh:155
const std::vector< LocalPoint > & CentersList
Definition: ref_element.hh:303
std::array< unsigned int, Size > IdxVector
Definition: ref_element.hh:144
const double epsilon
Definition: mathfce.h:23
arma::vec::fixed< dim > FaceBaryPoint
Definition: ref_element.hh:171
arma::vec::fixed< dim+1 > BaryPoint
Definition: ref_element.hh:170
IntersectionPoint * interpolate(const IntersectionPoint &A1, const IntersectionPoint &A2, double t)
static unsigned int topology_idx(unsigned int zeros_positions)
Definition: ref_element.hh:521
#define ASSERT_LT_DBG(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Definition: asserts.hh:299