Flow123d  JS_before_hm-2088-g1c00063e5
compute_intersection.cc
Go to the documentation of this file.
1 /*
2  * Author: viktor
3  */
4 
6 #include "mesh/ref_element.hh"
7 #include "mesh/mesh.h"
8 #include "mesh/accessors.hh"
9 #include "system/system.hh"
10 
11 #include "plucker.hh"
13 #include "intersection_aux.hh"
14 
15 using namespace std;
16 
17 /*************************************************************************************************************
18  * COMPUTE INTERSECTION FOR: 1D AND 2D
19  ************************************************************************************************************/
21 : computed_(false)
22 {
23  plucker_coordinates_abscissa_ = nullptr;
24  plucker_coordinates_triangle_.resize(3, nullptr);
25  plucker_products_.resize(3, nullptr);
26 }
27 
28 
30  ElementAccessor<3> triangle, FMT_UNUSED Mesh *mesh)
31 : computed_(false)
32 {
33  ASSERT_DBG(abscissa->dim() == 1);
34  ASSERT_DBG(triangle->dim() == 2);
35  // in this constructor, we suppose this is the final object -> we create all data members
36  plucker_coordinates_abscissa_ = new Plucker(*abscissa.node(0),
37  *abscissa.node(1), true);
38  scale_line_=plucker_coordinates_abscissa_->scale();
39 
40  plucker_coordinates_triangle_.resize(3);
41  plucker_products_.resize(3);
42  scale_triangle_=std::numeric_limits<double>::max();
43  for(unsigned int side = 0; side < 3; side++){
44  plucker_coordinates_triangle_[side] = new Plucker(*triangle.node(RefElement<2>::interact(Interaction<0,1>(side))[0]),
45  *triangle.node(RefElement<2>::interact(Interaction<0,1>(side))[1]),
46  true);
47  scale_triangle_ = std::min( scale_triangle_, plucker_coordinates_triangle_[side]->scale());
48 
49  // allocate and compute new Plucker products
50  plucker_products_[side] = new double((*plucker_coordinates_abscissa_)*(*plucker_coordinates_triangle_[side]));
51  }
52 }
53 
55 {
56  if(plucker_coordinates_abscissa_ != nullptr)
57  delete plucker_coordinates_abscissa_;
58 
59  for(unsigned int side = 0; side < RefElement<2>::n_sides; side++){
60  if(plucker_products_[side] != nullptr)
61  delete plucker_products_[side];
62  if(plucker_coordinates_triangle_[side] != nullptr)
63  delete plucker_coordinates_triangle_[side];
64  }
65 }
66 
67 
69  // unset all pointers
70  for(unsigned int side = 0; side < RefElement<2>::n_sides; side++)
71  {
72  plucker_products_[side] = nullptr;
73  plucker_coordinates_triangle_[side] = nullptr;
74  }
75  plucker_coordinates_abscissa_ = nullptr;
76 }
77 
79 
80  // if not already computed, compute plucker coordinates of abscissa
81  plucker_coordinates_abscissa_->compute();
82  scale_line_=plucker_coordinates_abscissa_->scale();
83 
84 // DBGMSG("Abscissa:\n");
85 // (*abscissa_)[0].point_coordinates().print();
86 // (*abscissa_)[1].point_coordinates().print();
87  scale_triangle_=std::numeric_limits<double>::max();
88  // if not already computed, compute plucker coordinates of triangle sides
89  for(unsigned int side = 0; side < RefElement<2>::n_sides; side++){
90  plucker_coordinates_triangle_[side]->compute();
91  scale_triangle_ = std::min( scale_triangle_, plucker_coordinates_triangle_[side]->scale());
92 
93  ASSERT_DBG(plucker_products_[side]).error("Undefined plucker product.");
94  if(*plucker_products_[side] == plucker_empty){
95  *plucker_products_[side] = (*plucker_coordinates_abscissa_)*(*plucker_coordinates_triangle_[side]);
96  }
97 // DBGMSG("Plucker product = %e\n", *(plucker_products_[side]));
98  }
99 }
100 
101 
103 {
104  double tol = geometry_epsilon * scale_line_;
105 
106  double t = RefElement<1>::bary_to_local(IP.local_bcoords_A())(0);
107 
108  int sign = (t < tol ? -2 : (t > 1+tol ? +2 : 0) ); // outside or inside
109  if (std::abs(t) <= tol){ // left end
110  sign = -1;
111  IP.set_topology_A(0,0);
112  }
113  else if (std::abs(1-t) <= tol){ // right end
114  sign = +1;
115  IP.set_topology_A(1,0);
116  }
117 
118  return sign;
119 }
120 
121 
123 {
124  // compute local barycentric coordinates of IP: see formula (3) on pg. 12 in BP VF
125  // local alfa = w2/sum; local beta = w1/sum; => local barycentric coordinates in the triangle
126  // where sum = w0+w1+w2
127 
128  // plucker products with correct signs according to ref_element, ordered by sides
129 // double w[3] = {signed_plucker_product(0),
130 // signed_plucker_product(1),
131 // signed_plucker_product(2)};
132 
133  // local barycentric coordinates of IP, depends on barycentric coordinates order !!!
134 // arma::vec3 local_triangle({w[2],w[1],w[0]});
135 // DBGMSG("Plucker product sum = %f\n",w[0]+w[1]+w[2]);
136 // local_triangle = local_triangle / (w[0]+w[1]+w[2]);
137 // double w_sum = local[0] + local[1] + local[2];
138 // DBGMSG("Plucker product sum = %e %e %e\n",w_sum, 1-rounding_epsilon, 1+rounding_epsilon);
139 
140  //assert inaccurate barycentric coordinates
141  ASSERT_DBG(fabs(1.0 - local[0] - local[1] - local[2]) < geometry_epsilon)(local[0]+local[1]+local[2])
142  (local[0])(local[1])(local[2]);
143 
144 
145 
146  arma::vec3 local_triangle({local[2],local[1],local[0]});
147 
148  // local coordinate T on the line
149  // for i-th coordinate it holds: (from formula (4) on pg. 12 in BP VF)
150  // T = localAbscissa= (- A(i) + ( 1 - alfa - beta ) * V0(i) + alfa * V1(i) + beta * V2 (i)) / U(i)
151  // let's choose [max,i] = max {U(i)}
152  arma::vec3 u = plucker_coordinates_abscissa_->get_u_vector();
153 
154  //find max in u in abs value:
155  unsigned int i = 0; // index of maximum in u
156  double max = u[0]; // maximum in u
157  if(fabs((double)u[1]) > fabs(max)){ max = u[1]; i = 1;}
158  if(fabs((double)u[2]) > fabs(max)){ max = u[2]; i = 2;}
159 
160  // global coordinates in triangle
161  double isect_coord_i =
162  local_triangle[0] * plucker_coordinates_triangle_[0]->point(RefElement<2>::normal_orientation(0))[i] +
163  local_triangle[1] * plucker_coordinates_triangle_[0]->point(1-RefElement<2>::normal_orientation(0))[i] +
164  local_triangle[2] * plucker_coordinates_triangle_[1]->point(RefElement<2>::normal_orientation(1))[i];
165 
166  //theta on abscissa
167  double t = (-plucker_coordinates_abscissa_->point(0)[i] + isect_coord_i)/max;
168  //DebugOut() << print_var(t) << print_var(isect_coord_i) << print_var(max);
169  arma::vec2 local_abscissa = {1-t, t};
170 
171  /*
172  DBGMSG("Coordinates: line; local and global triangle\n");
173  theta.print();
174  local_triangle.print();
175  global_triangle.print();
176  */
177 
178  IP.set_topology_A(0, 1);
179  IP.set_topology_B(0, 2);
180  IP.set_coordinates(local_abscissa,local_triangle);
181 
182  return true;
183 }
184 
185 
187 {
188  compute_plucker_products();
189  computed_ = true;
190 
191  //DebugOut() << "LINE: \n" << (*abscissa_)[0].point_coordinates()
192  // << (*abscissa_)[1].point_coordinates();
193 
194 
195  // convert plucker products to local coords
196  arma::vec3 w = {signed_plucker_product(0),
197  signed_plucker_product(1),
198  signed_plucker_product(2)};
199  double w_sum = w[0] + w[1] + w[2];
200 
201  unsigned int n_positive = 0;
202  unsigned int n_negative = 0;
203  unsigned int zero_idx_sum =0;
204  //DebugOut() << print_var(std::fabs(w_sum));
205  //DebugOut() << print_var(scale_line_);
206  //DebugOut() << print_var(scale_triangle_);
207 
208  double scaled_epsilon = geometry_epsilon*scale_line_*scale_triangle_*scale_triangle_;
209  if(std::fabs(w_sum) > scaled_epsilon) {
210  w = w / w_sum;
211  for (unsigned int i=0; i < 3; i++) {
212  //DebugOut() << print_var(i) << print_var(w[i]);
213  if (w[i] > geometry_epsilon) n_positive++;
214  else if ( w[i] > -geometry_epsilon) zero_idx_sum+=i;
215  else n_negative++;
216  }
217 
218  } else {
219  /* case 'w_sum == 0':
220  * 1] all products are zero => n_negative=0 and n_positive=0 => it is degenerate case (coplanar case)
221  * 2] at least two products are nonzero AND some of them must be negative => no intersection
222  * (it happens when line is parallel with the triangle but not coplanar; unit test line_triangle09.msh)
223  * See the IF conditions below.
224  */
225 
226  for (unsigned int i=0; i < 3; i++) {
227  //DebugOut().fmt("i: {} w[i]: {:g}", i, w[i]);
228  if (w[i] > scaled_epsilon || w[i] < -scaled_epsilon) n_negative++;
229  }
230  // n_positive == 0
231  //DebugOut() << print_var(n_negative);
232  }
233 
234 // w.print("w");
235 
236  // any negative barycentric coordinate means, no intersection
237  if (n_negative>0) return IntersectionResult::none;
238 
239  // test whether any plucker products is non-zero
240  if (n_positive > 0) {
241 
242  compute_plucker(IP, w);
243  // edge of triangle
244  unsigned int non_zero_idx=0;
245  if (n_positive == 2) {
246  // one zero product, intersection on the zero edge
247  // the zero edge index is equal to zero_idx_sum
248  IP.set_topology_B(zero_idx_sum, 1);
249  non_zero_idx = (zero_idx_sum + 1) % 3;
250  }
251  else if (n_positive == 1) {
252  // two zero products, intersection in vertex oposite to single non-zero edge
253  // index of the non-zero edge is 3-zero_idx_sum
254  IP.set_topology_B(RefElement<2>::oposite_node(3-zero_idx_sum), 0);
255  non_zero_idx = 3-zero_idx_sum;
256  }
257  //DebugOut() << print_var(non_zero_idx) << print_var(signed_plucker_product(non_zero_idx));
258 
259  IntersectionResult result = signed_plucker_product(non_zero_idx) > 0 ?
261  IP.set_result(result);
262 
263  return result;
264 
265  } else {
266  ASSERT_DBG(IP.topology_equal(IPAux12())); // check empty IP (none)
267  IP.set_result(IntersectionResult::degenerate); // set deg. result
269  }
270 }
271 
273 {
274  IPAux12 IP;
275  IntersectionResult result = compute(IP);
276 // DBGVAR((int)result);
277  // skip empty cases
278  if(result == IntersectionResult::none) return 0;
279 
280  // standard case with a single intersection corner
281  if(result < IntersectionResult::degenerate){
282 // DBGCOUT(<< "12d plucker case\n");
283  int sign = check_abscissa_topology(IP);
284 // DBGVAR(sign);
285  if(std::abs(sign) > 1) return 0;
286 
287  IP12s.push_back(IP);
288  return IP12s.size();
289  }
290  else{
292 // DBGCOUT(<< "12d degenerate case, all products are zero\n");
293  return compute_final_in_plane(IP12s);
294  }
295 }
296 
297 
299  IPAux12& IP)
300 {
301 // DBGMSG("PluckerProduct[%d]: %f\n",side, *plucker_products_[side]);
302 
303  // We solve following equation for parameters s,t:
304  /* A + sU = C + tV = intersection point
305  * sU - tV = C - A
306  *
307  * which is by components:
308  * (u1 -v1) (s) = (c1-a1)
309  * (u2 -v2) (t) = (c2-a2)
310  * (u3 -v3) = (c3-a3)
311  *
312  * these are 3 equations for variables s,t
313  * see (4.3) on pg. 19 in DP VF
314  *
315  * We will solve this using Crammer's rule for the maximal subdeterminant det_ij of the matrix.
316  * s = detX_ij / det_ij
317  * t = detY_ij / det_ij
318  */
319 
320  // starting point of abscissa
321  arma::vec3 A = plucker_coordinates_abscissa_->point(0);
322  // direction vector of abscissa
323  arma::vec3 U = plucker_coordinates_abscissa_->get_u_vector();
324  // vertex of triangle side
325  arma::vec3 C = plucker_coordinates_triangle_[side]->point(RefElement<2>::normal_orientation(side));
326  // direction vector of triangle side
327  arma::vec3 V = plucker_coordinates_triangle_[side]->get_u_vector();
328  // right hand side
329  arma::vec3 K = C - A;
330  // subdeterminants det_ij of the system equal minus normal vector to common plane of U and V
331  // det_12 = (-UxV)[1]
332  // det_13 = -(-UxV)[2]
333  // det_23 = (-UxV)[3]
334  arma::vec3 Det = -arma::cross(U,V);
335 
336 
337 // A.print("A");
338 // U.print("U");
339 // C.print("C");
340 // V.print("V");
341 // K.print("K");
342 // Det.print();
343 // cout <<endl;
344 
345  unsigned int max_index = 0;
346  double maximum = fabs(Det[0]);
347  if(fabs((double)Det[1]) > maximum){
348  maximum = fabs(Det[1]);
349  max_index = 1;
350  }
351  if(fabs((double)Det[2]) > maximum){
352  maximum = fabs(Det[2]);
353  max_index = 2;
354  }
355 // DBGVAR(maximum);
356  //abscissa is parallel to triangle side
357  //TODO What tolerance should we use here and why?
358 // if(std::abs(maximum) <= std::sqrt(geometry_epsilon)) return false;
359  if(std::abs(maximum) <= geometry_epsilon) return false;
360 
361  // map maximum index in {-UxV} to i,j of subdeterminants
362  // i j
363  // max_index 0: 1 2
364  // 1: 2 0 (switch due to sign change)
365  // 2: 0 1
366  unsigned int i = (max_index+1)%3,
367  j = (max_index+2)%3;
368 
369  double DetX = -K[i]*V[j] + K[j]*V[i];
370  double DetY = -K[i]*U[j] + K[j]*U[i];
371 
372  double s = DetX/Det[max_index]; //parameter on abscissa
373  double t = DetY/Det[max_index]; //parameter on triangle side
374 
375  // change sign according to side orientation
376  if(RefElement<2>::normal_orientation(side)) t=-t;
377 
378 // DBGVAR(s);
379 // DBGVAR(t);
380 
381  //TODO correct tolerance with scale; compute scale without plucker coords
382  double tol = geometry_epsilon;
383  // if IP is inside of triangle side
384  if(t >= -tol && t <= 1+tol){
385 
386  IP.set_result(IntersectionResult::degenerate); // set orientation as a pathologic case ( > 1)
387 
388  // possibly set triangle vertex {0,1,2}
390  else if(fabs(1-t) <= tol) { IP.set_topology_B(RefElement<2>::interact(Interaction<0,1>(side))[1-RefElement<2>::normal_orientation(side)],0);}
391  else { IP.set_topology_B(side,1);} // no vertex, side side, dim = 1
392 
393  arma::vec2 local_abscissa({1-s, s});
394  arma::vec3 local_triangle({0,0,0});
395 
396  // set local triangle barycentric coordinates according to nodes of the triangle side:
397  local_triangle[RefElement<2>::interact(Interaction<0,1>(side))[RefElement<2>::normal_orientation(side)]] = 1 - t;
399 // local_abscissa.print();
400 // local_triangle.print();
401 
402  IP.set_coordinates(local_abscissa,local_triangle);
403  return true; // IP found
404  }
405 
406  return false; // IP NOT found
407 }
408 
409 
411 {
412  std::vector<int> sign;
413 
414  for(unsigned int i = 0; i < 3;i++){
415 // DBGCOUT( << "side [" << i << "]\n");
416  if(IP12s.size() == 2) break;
417 
418  IPAux12 IP;
419  if (compute_degenerate(i,IP))
420  {
421  uint s = check_abscissa_topology(IP);
422 // cout << IP;
423  // check whether we have found the same IP (e.g. vertex case)
424  if(IP12s.size() > 0 && IP.topology_equal(IP12s.back())) continue;
425 
426  sign.push_back(s);
427  IP12s.push_back(IP);
428  }
429  }
430  if(IP12s.size() == 0) return IP12s.size();
431 
432  // in the case, that line goes through vertex, but outside tetrahedron (touching vertex)
433  if(IP12s.size() == 1){
434  if(std::abs(sign[0]) > 1) IP12s.pop_back(); // outside abscissa
435  return IP12s.size();
436  }
437 
438  // 2 IPs CASE:
439  ASSERT_EQ_DBG(2, IP12s.size());
440 
441 // DBGVAR(sign[0]);
442 // DebugOut() << IP12s[0];
443 // DBGVAR(sign[1]);
444 // DebugOut() << IP12s[1];
445 
446  // intersection outside of abscissa => NO intersection
447  if (sign[0] == sign[1] && std::abs(sign[0]) > 1) {
448  IP12s.clear();
449  return 0;
450  }
451 
452  // order IPs according to the abscissa parameter
453  if(IP12s[0].local_bcoords_A()[1] > IP12s[1].local_bcoords_A()[1]){
454  std::swap(IP12s[0], IP12s[1]);
455  std::swap(sign[0], sign[1]);
456  }
457 
458  // possibly cut IPs to abscissa ends and interpolate tetrahedron bcoords
459  const int ip_sign[] = {-2, +2}; // states to cut
460  for( unsigned int ip=0; ip<2; ip++) {
461  // cut every IP to its end of the abscissa
462  if (sign[ip] == ip_sign[ip]) {
463  sign[ip] /=2; // -2 to -1; +2 to +1
464  correct_triangle_ip_topology(double(ip), ip, IP12s);
465  IP12s[ip].set_topology_A(ip, 0);
466  }
467  }
468 
469  // if IPs are the same, then throw the second one away
470  if(IP12s[0].topology_equal(IP12s[1])){
471 // DebugOut() << "IP same, pop up\n";
472  IP12s.pop_back();
473  }
474 
475  return IP12s.size();
476 }
477 
479  double t, unsigned int ip, std::vector<IPAux12> &ips)
480 {
482  ips[0].local_bcoords_B(), ips[1].local_bcoords_B(),
483  ips[0].local_bcoords_A()[1], ips[1].local_bcoords_A()[1], t);
484  arma::vec2 local_abscissa({1 - t, t}); // abscissa local barycentric coords
485  ips[ip].set_coordinates(local_abscissa, local_tria);
486 
487  // create mask for zeros in barycentric coordinates
488  // coords (*,*,*,*) -> byte bitwise xxxx
489  // only least significant one byte used from the integer
490  std::pair<unsigned int, unsigned int> zeros =
491  RefElement<2>::zeros_positions(ips[ip].local_bcoords_B(), geometry_epsilon);
492 
493  /**
494  * TODO:
495  * 1. Try to avoid setting topology from coords. Try to use just topology information.
496  * 2. If can not be done. Use interact method to setup a map mapping 16 possible zeros positions to appropriate topology,
497  * remove topology_idx from RefElement.
498  */
499 
500 // DebugOut() << "zeros: <" << zeros.first << ", " << zeros.second << ">\n";
501  switch(zeros.first)
502  {
503  default: ips[ip].set_topology_B(0,2); //inside tetrahedon
504  break;
505  case 1: ips[ip].set_topology_B(RefElement<2>::topology_idx<1>(zeros.second),1);
506  break;
507  case 2: ips[ip].set_topology_B(RefElement<2>::topology_idx<0>(zeros.second),0);
508  break;
509  }
510 }
511 
512 
514  os << "\tPluckerCoordinates Abscissa[0]";
515  if(plucker_coordinates_abscissa_ == nullptr){
516  os << "NULL" << endl;
517  }else{
518  os << plucker_coordinates_abscissa_;
519  }
520  for(unsigned int i = 0; i < 3;i++){
521  os << "\tPluckerCoordinates Triangle[" << i << "]";
522  if(plucker_coordinates_triangle_[i] == nullptr){
523  os << "NULL" << endl;
524  }else{
525  os << plucker_coordinates_triangle_[i];
526  }
527  }
528 }
529 
530 
531 
532 
533 /*************************************************************************************************************
534  * COMPUTE INTERSECTION FOR: 2D AND 2D
535  ************************************************************************************************************/
537 {
538  plucker_coordinates_.resize(2*RefElement<2>::n_sides, nullptr);
539  plucker_products_.resize(3*RefElement<2>::n_sides, nullptr);
540 }
541 
543  ElementAccessor<3> triaB,
544  FMT_UNUSED Mesh *mesh)
545 {
546  ASSERT_DBG(triaA->dim() == 2);
547  ASSERT_DBG(triaB->dim() == 2);
548  plucker_coordinates_.resize(2*RefElement<2>::n_sides);
549  plucker_products_.resize(3*RefElement<2>::n_sides);
550 
551  for(unsigned int side = 0; side < RefElement<2>::n_sides; side++){
552  plucker_coordinates_[side] = new Plucker(*triaA.node(RefElement<2>::interact(Interaction<0,1>(side))[0]),
553  *triaA.node(RefElement<2>::interact(Interaction<0,1>(side))[1]));
554  plucker_coordinates_[RefElement<2>::n_sides+side]
555  = new Plucker(*triaB.node(RefElement<2>::interact(Interaction<0,1>(side))[0]),
556  *triaB.node(RefElement<2>::interact(Interaction<0,1>(side))[1]));
557  }
558 
559  // compute Plucker products for each pair triangle A side and triangle B side
560  for(unsigned int p = 0; p < 3*RefElement<2>::n_sides; p++){
561  plucker_products_[p] = new double(plucker_empty);
562  }
563 }
564 
566 {
567  // unset pointers:
568  for(unsigned int side = 0; side < 2*RefElement<2>::n_sides; side++)
569  CI12[side].clear_all();
570 
571  // then delete objects:
572  for(unsigned int side = 0; side < 2*RefElement<2>::n_sides; side++){
573  if(plucker_coordinates_[side] != nullptr)
574  delete plucker_coordinates_[side];
575  }
576 
577  for(unsigned int p = 0; p < 3*RefElement<2>::n_sides; p++){
578  if(plucker_products_[p] != nullptr)
579  delete plucker_products_[p];
580  }
581 }
582 
584 {
585  // unset all pointers
586  for(unsigned int side = 0; side < 2*RefElement<2>::n_sides; side++)
587  {
588  plucker_coordinates_[side] = nullptr;
589  }
590  for(unsigned int p = 0; p < 3*RefElement<2>::n_sides; p++){
591  plucker_products_[p] = nullptr;
592  }
593 }
594 
596 
597 // DBGMSG("init\n");
598  for(unsigned int i = 0; i < RefElement<2>::n_sides; i++){
599  // set side A vs triangle B
600  for(unsigned int j = 0; j < RefElement<2>::n_sides; j++)
601  CI12[i].set_pc_triangle(plucker_coordinates_[3+j], j); // set triangle B
602  // set side of triangle A
603  CI12[i].set_pc_abscissa(plucker_coordinates_[i]);
604 
605  // set side B vs triangle A
606  for(unsigned int j = 0; j < RefElement<2>::n_sides; j++)
607  CI12[RefElement<2>::n_sides + i].set_pc_triangle(plucker_coordinates_[j], j); // set triangle A
608  // set side of triangle B
609  CI12[RefElement<2>::n_sides + i].set_pc_abscissa(plucker_coordinates_[RefElement<2>::n_sides + i]);
610 
611  // set plucker products
612  for(unsigned int j = 0; j < RefElement<2>::n_sides; j++)
613  {
614  //for A[i]_B set pp. A[i] x B[j]
615  CI12[i].set_plucker_product(get_plucker_product(i,j),j);
616  //for B[i]_A set pp. A[j] x B[i]
617  CI12[RefElement<2>::n_sides + i].set_plucker_product(get_plucker_product(j,i),j);
618  }
619 
620  }
621 }
622 
623 
625 {
626  // final intersection points
627  std::vector<IPAux22> &IP22s = intersection.points();
628  // temporary vector for lower dimensional IPs
629  std::vector<IPAux12> IP12s;
630  IP12s.reserve(2);
631  unsigned int ip_coutner = 0;
632 
633  // loop over CIs (side vs triangle): [A0_B, A1_B, A2_B, B0_A, B1_A, B2_A].
634  for(unsigned int i = 0; i < 2*RefElement<2>::n_sides && ip_coutner < 2; i++){
635  if(!CI12[i].is_computed()) // if not computed yet
636  {
637 // DBGVAR(i);
638  if(CI12[i].compute_final(IP12s) == 0) continue;
639 
640  unsigned int triangle_side = i%3; //i goes from 0 to 5 -> i%3 = 0,1,2,0,1,2
641 
642  for(IPAux12 &IP : IP12s)
643  {
644  IPAux22 IP22(IP.switch_objects(), triangle_side); // swicth dim 12 -> 21; interpolate dim 21 -> 22
645 
646  if(i < 3){
647  //switch back to keep order of triangles [A,B]
648  IP22 = IP22.switch_objects();
649 
650  if( IP22.dim_A() == 0 ) // if IP is vertex of triangle
651  {
652 // DBGCOUT("set_node A\n");
653  // we are on line of the triangle A, and IP.idx_A contains local node of the line
654  // we know vertex index
655 
656  // set flag on all sides of tetrahedron connected by the node
657  for(unsigned int s=0; s < RefElement<2>::n_sides_per_node; s++)
658  CI12[RefElement<2>::interact(Interaction<1,0>(IP22.idx_A()))[s]].set_computed();
659  }
660  if( IP22.dim_B() == 0 ) // if IP is vertex of triangle
661  {
662 // DBGCOUT("set_node B\n");
663  // set flag on both sides of triangle connected by the node
664  for(unsigned int s=0; s < RefElement<2>::n_sides_per_node; s++)
665  CI12[3 + RefElement<2>::interact(Interaction<1,0>(IP22.idx_B()))[s]].set_computed();
666  }
667  else if( IP22.dim_B() == 1 ) // if IP is vertex of triangle
668  {
669 // DBGCOUT("set line B\n");
670  // set flag on both sides of triangle connected by the node
671  CI12[3 + IP22.idx_B()].set_computed();
672  }
673  }
674  else if( IP22.dim_B() == 0 ) // if IP is vertex of triangle B (triangles switched!!! A <-> B)
675  {
676  //we do not need to look back to triangle A, because if IP was at vertex, we would know already
677 // DBGCOUT("set_node B\n");
678  // we are on line of the triangle, and IP.idx_A contains local node of the line
679  // we know vertex index
680 
681  // set flag on both sides of triangle connected by the node
682  for(unsigned int s=0; s < RefElement<2>::n_sides_per_node; s++)
683  CI12[3 + RefElement<2>::interact(Interaction<1,0>(IP22.idx_B()))[s]].set_computed();
684  }
685 // DBGCOUT( << IP22);
686  ip_coutner++;
687  IP22s.push_back(IP22);
688  }
689  IP12s.clear();
690  }
691  }
692 
693  return ip_coutner;
694 }
695 
696 // void ComputeIntersection< 2, 2>::correct_triangle_ip_topology(IntersectionPointAux<2,2>& ip)
697 // {
698 // // create mask for zeros in barycentric coordinates
699 // // coords (*,*,*,*) -> byte bitwise xxxx
700 // // only least significant one byte used from the integer
701 // unsigned int zeros = 0;
702 // unsigned int n_zeros = 0;
703 // for(char i=0; i < 3; i++){
704 // if(std::fabs(ip.local_bcoords_B()[i]) < geometry_epsilon)
705 // {
706 // zeros = zeros | (1 << (2-i));
707 // n_zeros++;
708 // }
709 // }
710 //
711 // switch(n_zeros)
712 // {
713 // default: ip.set_topology_B(0,2); //inside triangle
714 // break;
715 // case 1: ip.set_topology_B(RefElement<2>::topology_idx<1>(zeros),2);
716 // break;
717 // case 2: ip.set_topology_B(RefElement<2>::topology_idx<0>(zeros),1);
718 // break;
719 // }
720 // }
721 
722 
724 
725  for(unsigned int i = 0; i < RefElement<2>::n_lines; i++){
726  os << "\tPluckerCoordinates Triangle A[" << i << "]";
727  if(plucker_coordinates_[i] == nullptr){
728  os << "NULL" << endl;
729  }else{
730  os << plucker_coordinates_[i];
731  }
732  }
733  for(unsigned int i = 0; i < RefElement<2>::n_lines; i++){
734  os << "\tPluckerCoordinates Triangle B[" << i << "]";
735  if(plucker_coordinates_[RefElement<2>::n_lines + i] == nullptr){
736  os << "NULL" << endl;
737  }else{
738  os << plucker_coordinates_[RefElement<2>::n_lines + i];
739  }
740  }
741 }
742 
744  os << "ComputeIntersection<2,2> Plucker Coordinates Tree:" << endl;
745  print_plucker_coordinates(os);
746  for(unsigned int i = 0; i < 6;i++){
747  os << "ComputeIntersection<1,2>["<< i <<"] Plucker Coordinates:" << endl;
748  CI12[i].print_plucker_coordinates(os);
749  }
750 }
751 
752 
753 /*************************************************************************************************************
754  * COMPUTE INTERSECTION FOR: 1D AND 3D
755  ************************************************************************************************************/
757 {
758  plucker_coordinates_abscissa_ = nullptr;
759  plucker_coordinates_tetrahedron.resize(6, nullptr);
760  plucker_products_.resize(6, nullptr);
761 }
762 
764  ElementAccessor<3> tetrahedron,
765  FMT_UNUSED Mesh *mesh)
766 {
767  ASSERT_DBG(abscissa->dim() == 1);
768  ASSERT_DBG(tetrahedron->dim() == 3);
769  ASSERT_DBG(tetrahedron.sign() * tetrahedron.jacobian_S3() > 0).add_value(tetrahedron.input_id(),"element index").error(
770  "Tetrahedron element (%d) has wrong numbering or is degenerated (negative Jacobian).");
771 
772 
773  plucker_coordinates_abscissa_ = new Plucker(*abscissa.node(0), *abscissa.node(1));
774  plucker_coordinates_tetrahedron.resize(6);
775  plucker_products_.resize(6);
776 
777  for(unsigned int line = 0; line < RefElement<3>::n_lines; line++){
778  plucker_coordinates_tetrahedron[line] = new Plucker(*tetrahedron.node(RefElement<3>::interact(Interaction<0,1>(line))[0]),
779  *tetrahedron.node(RefElement<3>::interact(Interaction<0,1>(line))[1]));
780  // compute Plucker products (abscissa X tetrahedron line)
781  plucker_products_[line] = new double(plucker_empty);
782  }
783 }
784 
786 {
787  // unset pointers:
788  for(unsigned int side = 0; side < RefElement<3>::n_sides; side++)
789  CI12[side].clear_all();
790 
791  // then delete objects:
792  if(plucker_coordinates_abscissa_ != nullptr)
793  delete plucker_coordinates_abscissa_;
794 
795  for(unsigned int line = 0; line < RefElement<3>::n_lines; line++){
796  if(plucker_products_[line] != nullptr)
797  delete plucker_products_[line];
798  if(plucker_coordinates_tetrahedron[line] != nullptr)
799  delete plucker_coordinates_tetrahedron[line];
800  }
801 }
802 
804 {
805  // unset all pointers
806  for(unsigned int side = 0; side < RefElement<3>::n_lines; side++)
807  {
808  plucker_products_[side] = nullptr;
809  plucker_coordinates_tetrahedron[side] = nullptr;
810  }
811  plucker_coordinates_abscissa_ = nullptr;
812 }
813 
815 
816 
817  for(unsigned int side = 0; side < RefElement<3>::n_sides; side++){
818  for(unsigned int line = 0; line < RefElement<3>::n_lines_per_side; line++){
819  CI12[side].set_pc_triangle(plucker_coordinates_tetrahedron[
820  RefElement<3>::interact(Interaction<1,2>(side))[line]], line);
821 
822  CI12[side].set_plucker_product(
823  plucker_products_[RefElement<3>::interact(Interaction<1,2>(side))[line]],
824  line);
825  }
826  CI12[side].set_pc_abscissa(plucker_coordinates_abscissa_);
827  }
828 }
829 
830 
832 {
833  return compute(intersection.i_points_);
834 }
835 
837  ASSERT_EQ_DBG(0, IP13s.size());
838 
839  // Topological position of abscissa -2, -1, 0, 1, 2; see 'check_abscissa_topology' method
840  std::vector<int> position_on_abscissa;
841 
842 
843  // loop over faces of tetrahedron
844  for(unsigned int face = 0; face < RefElement<3>::n_sides && IP13s.size() < 2; face++){
845 
846  if (CI12[face].is_computed()) continue;
847 
849  IntersectionResult result = CI12[face].compute(IP12);
850  //DebugOut() << print_var(face) << print_var(int(result)) << "1d-3d";
851 
852  if (int(result) < int(IntersectionResult::degenerate) ) {
853  // fix topology of A (abscissa) in IP12 and save sign
854  position_on_abscissa.push_back(CI12[face].check_abscissa_topology(IP12));
855 
856  IPAux IP13(IP12, face);
857  //DebugOut() << IP13;
858  IP13s.push_back(IP13);
859 
860  // set the 'computed' flag on the connected sides by IP
861  if(IP13.dim_B() == 0) // IP is vertex of triangle
862  {
863  // set flag on all sides of tetrahedron connected by the node
864  for(unsigned int node_face : RefElement<3>::interact(Interaction<2,0>(IP13.idx_B())))
865  CI12[node_face].set_computed();
866  }
867  else if(IP13.dim_B() == 1) // IP is on edge of triangle
868  {
869  for(unsigned int edge_face : RefElement<3>::interact(Interaction<2,1>(IP13.idx_B())))
870  CI12[edge_face].set_computed();
871  }
872  }
873  }
874  if (IP13s.size() == 0) return IP13s.size();
875 
876  // in the case, that line goes through vertex, but outside tetrahedron (touching vertex)
877  if(IP13s.size() == 1){
878  if(std::abs(position_on_abscissa[0]) > 1) IP13s.pop_back(); // outside abscissa
879  return IP13s.size();
880  }
881 
882  // 2 IPs CASE:
883  ASSERT_EQ_DBG(2, IP13s.size());
884 
885  // intersection outside of abscissa => NO intersection
886  if (position_on_abscissa[0] == position_on_abscissa[1] && std::abs(position_on_abscissa[0]) > 1) {
887  IP13s.clear();
888  return 0;
889  }
890 
891  // order IPs according to the abscissa parameter
892  if(IP13s[0].local_bcoords_A()[1] > IP13s[1].local_bcoords_A()[1]){
893  std::swap(IP13s[0], IP13s[1]);
894  std::swap(position_on_abscissa[0], position_on_abscissa[1]);
895  }
896 
897  // possibly cut IPs to abscissa ends and interpolate tetrahedron bcoords
898  const int ip_sign[] = {-2, +2}; // states to cut
899  for( unsigned int ip=0; ip<2; ip++) {
900  // cut every IP to its end of the abscissa
901  if (position_on_abscissa[ip] == ip_sign[ip]) {
902  position_on_abscissa[ip] /=2; // -2 to -1; +2 to +1
903  correct_tetrahedron_ip_topology(double(ip), ip, IP13s);
904  IP13s[ip].set_topology_A(ip, 0);
905  }
906  }
907 
908  // if IPs are the same, then throw the second one away
909  if(IP13s[0].topology_equal(IP13s[1])){
910  IP13s.pop_back();
911  }
912 
913  return IP13s.size();
914 }
915 
917  double t, unsigned int ip, std::vector<IPAux> &ips)
918 {
919  arma::vec4 local_tetra = RefElement<3>::line_barycentric_interpolation(
920  ips[0].local_bcoords_B(), ips[1].local_bcoords_B(),
921  ips[0].local_bcoords_A()[1], ips[1].local_bcoords_A()[1], t);
922  arma::vec2 local_abscissa({1 - t, t}); // abscissa local barycentric coords
923  ips[ip].set_coordinates(local_abscissa, local_tetra);
924 
925  // create mask for zeros in barycentric coordinates
926  // coords (*,*,*,*) -> byte bitwise xxxx
927  // only least significant one byte used from the integer
928  unsigned int zeros = 0;
929  unsigned int n_zeros = 0;
930  for(char i=0; i < 4; i++){
931  if(std::fabs(ips[ip].local_bcoords_B()[i]) < geometry_epsilon)
932  {
933  zeros = zeros | (1 << i);
934  n_zeros++;
935  }
936  }
937 
938  /**
939  * TODO:
940  * 1. Try to avoid setting topology from coords. Try to use just topology information.
941  * 2. If can not be done. Use interact method to setup a map mapping 16 possible zeros positions to appropriate topology,
942  * remove topology_idx from RefElement.
943  */
944 
945  switch(n_zeros)
946  {
947  default: ips[ip].set_topology_B(0,3); //inside tetrahedon
948  break;
949  case 1: ips[ip].set_topology_B(RefElement<3>::topology_idx<2>(zeros),2);
950  break;
951  case 2: ips[ip].set_topology_B(RefElement<3>::topology_idx<1>(zeros),1);
952  break;
953  case 3: ips[ip].set_topology_B(RefElement<3>::topology_idx<0>(zeros),0);
954  break;
955  }
956 }
957 
958 
960  os << "\tPluckerCoordinates Abscissa[0]";
961  if(plucker_coordinates_abscissa_ == nullptr){
962  os << "NULL" << endl;
963  }else{
964  os << plucker_coordinates_abscissa_;
965  }
966 
967  for(unsigned int i = 0; i < 6;i++){
968  os << "\tPluckerCoordinates Tetrahedron[" << i << "]";
969  if(plucker_coordinates_tetrahedron[i] == nullptr){
970  os << "NULL" << endl;
971  }else{
972  os << plucker_coordinates_tetrahedron[i];
973  }
974  }
975 }
976 
978  os << "ComputeIntersection<1,3> Plucker Coordinates Tree:" << endl;
979  print_plucker_coordinates(os);
980  for(unsigned int i = 0; i < 4;i++){
981  os << "ComputeIntersection<1,2>["<< i <<"] Plucker Coordinates:" << endl;
982  CI12[i].print_plucker_coordinates(os);
983  }
984 }
985 
986 
987 
988 /*************************************************************************************************************
989  * COMPUTE INTERSECTION FOR: 2D AND 3D
990  ************************************************************************************************************/
992 : no_idx(100),
993 s3_dim_starts({0, 4, 10, 14}), // vertices, edges, faces, volume
994 s2_dim_starts({15, 18, 21}), // vertices, sides, surface
995 object_next(22, no_idx), // 4 vertices, 6 edges, 4 faces, 1 volume, 3 corners, 3 sides, 1 surface; total 22
996 on_faces(_on_faces())
997  {
998 
999  plucker_coordinates_triangle_.resize(3, nullptr);
1000  plucker_coordinates_tetrahedron.resize(6, nullptr);
1001  plucker_products_.resize(3*6, nullptr);
1002 }
1003 
1004 
1006  ElementAccessor<3> tetrahedron,
1007  Mesh *mesh)
1009 {
1010  ASSERT_DBG(tetrahedron.sign() * tetrahedron.jacobian_S3() > 0)
1011  (tetrahedron.sign())(tetrahedron.jacobian_S3()).add_value(tetrahedron.input_id(),"element index").error(
1012  "Tetrahedron element (%d) has wrong numbering or is degenerated (negative Jacobian).");
1013 
1014  mesh_ = mesh;
1015  plucker_coordinates_triangle_.resize(3);
1016  plucker_coordinates_tetrahedron.resize(6);
1017 
1018  // set CI object for 1D-2D intersection 'tetrahedron edge - triangle'
1019  for(unsigned int i = 0; i < RefElement<3>::n_lines; i++){
1020  plucker_coordinates_tetrahedron[i] = new Plucker(*tetrahedron.node(RefElement<3>::interact(Interaction<0,1>(i))[0]),
1021  *tetrahedron.node(RefElement<3>::interact(Interaction<0,1>(i))[1]));
1022  }
1023  // set CI object for 1D-3D intersection 'triangle side - tetrahedron'
1024  for(unsigned int i = 0; i < RefElement<2>::n_lines;i++){
1025  plucker_coordinates_triangle_[i] = new Plucker(*triangle.node(RefElement<2>::interact(Interaction<0,1>(i))[0]),
1026  *triangle.node(RefElement<2>::interact(Interaction<0,1>(i))[1]));
1027  }
1028 
1029  // compute Plucker products (triangle side X tetrahedron line)
1030  // order: triangle sides X tetrahedron lines:
1031  // TS[0] X TL[0..6]; TS[1] X TL[0..6]; TS[1] X TL[0..6]
1032  unsigned int np = RefElement<2>::n_sides * RefElement<3>::n_lines;
1033  plucker_products_.resize(np, nullptr);
1034  for(unsigned int line = 0; line < np; line++){
1035  plucker_products_[line] = new double(plucker_empty);
1036 
1037  }
1038 }
1039 
1041 {
1042  // unset pointers:
1043  for(unsigned int triangle_side = 0; triangle_side < RefElement<2>::n_sides; triangle_side++)
1044  CI13[triangle_side].clear_all();
1045 
1046  for(unsigned int line = 0; line < RefElement<3>::n_lines; line++)
1047  CI12[line].clear_all();
1048 
1049  // then delete objects:
1050  unsigned int np = RefElement<2>::n_sides * RefElement<3>::n_lines;
1051  for(unsigned int line = 0; line < np; line++){
1052  if(plucker_products_[line] != nullptr)
1053  delete plucker_products_[line];
1054  }
1055 
1056  for(unsigned int i = 0; i < RefElement<3>::n_lines;i++){
1057  if(plucker_coordinates_tetrahedron[i] != nullptr)
1058  delete plucker_coordinates_tetrahedron[i];
1059  }
1060  for(unsigned int i = 0; i < RefElement<2>::n_sides;i++){
1061  if(plucker_coordinates_triangle_[i] != nullptr)
1062  delete plucker_coordinates_triangle_[i];
1063  }
1064 }
1065 
1067 
1068  // set pointers to Plucker coordinates for 1D-2D
1069  // set pointers to Plucker coordinates for 1D-3D
1070  // distribute Plucker products - CI13
1071  for(unsigned int triangle_side = 0; triangle_side < RefElement<2>::n_sides; triangle_side++){
1072  for(unsigned int line = 0; line < RefElement<3>::n_lines; line++){
1073  CI13[triangle_side].set_plucker_product(
1074  plucker_products_[triangle_side * RefElement<3>::n_lines + line],
1075  line);
1076  CI12[line].set_plucker_product(
1077  plucker_products_[triangle_side * RefElement<3>::n_lines + line],
1078  triangle_side);
1079 
1080  CI13[triangle_side].set_pc_tetrahedron(plucker_coordinates_tetrahedron[line],line);
1081  CI12[line].set_pc_triangle(plucker_coordinates_triangle_[triangle_side],triangle_side);
1082  }
1083  CI13[triangle_side].set_pc_abscissa(plucker_coordinates_triangle_[triangle_side]);
1084  CI13[triangle_side].init();
1085  }
1086 
1087  // set pointers to Plucker coordinates for 1D-2D
1088  for(unsigned int line = 0; line < RefElement<3>::n_lines; line++)
1089  CI12[line].set_pc_abscissa(plucker_coordinates_tetrahedron[line]);
1090 
1091 }
1092 
1093 
1094 
1096  ASSERT_LT_DBG(i_obj, object_next.size());
1097  unsigned int ip = object_next[i_obj];
1098  if (ip == no_idx) return false;
1099  ASSERT_LT_DBG(ip, IP_next.size());
1100  return IP_next[ip] == i_obj;
1101 }
1102 
1103 /**
1104  * Set links: obj_before -> IP -> obj_after
1105  * if obj_after have null successor, set obj_after -> IP (backlink)
1106  */
1107 void ComputeIntersection<2,3>::set_links(uint obj_before_ip, uint ip_idx, uint obj_after_ip) {
1108  if (have_backlink(obj_after_ip)) {
1109  // target object is already target of other IP, so it must be source object
1110  std::swap(obj_before_ip, obj_after_ip);
1111  }
1112 
1113  unsigned int ip = object_next[obj_after_ip];
1114  ASSERT_DBG( ! have_backlink(obj_after_ip) )
1115  (mesh_->element_accessor(intersection_->component_ele_idx()).idx())
1116  (mesh_->element_accessor(intersection_->bulk_ele_idx()).idx())
1117  (obj_before_ip)
1118  (ip_idx)
1119  (obj_after_ip)
1120  (IP_next[ip])
1121  (object_next[obj_after_ip]); // at least one could be target object
1122 
1123  //DebugOut().fmt("SET LINK: {} -> ip: {} -> {}\n", obj_before_ip, ip_idx, obj_after_ip );
1124  object_next[obj_before_ip] = ip_idx;
1125  IP_next.push_back( obj_after_ip);
1126  if (object_next[obj_after_ip] == no_idx) {
1127  object_next[obj_after_ip] = ip_idx;
1128  }
1129 }
1130 
1131 
1132 
1134 {
1135  S3_inverted = mesh_->element_accessor(intersection.bulk_ele_idx()).inverted();
1136  intersection_= &intersection;
1137  //DebugOut().fmt("2d ele: {} 3d ele: {}\n",
1138  // intersection.component_ele_idx(),
1139  // intersection.bulk_ele_idx());
1140 
1141  IP23_list.clear();
1142  IP_next.clear();
1143  std::fill(object_next.begin(), object_next.end(), no_idx);
1144  std::vector<IPAux13> IP13s;
1145 
1146  std::array<bool, 6> edge_touch={false,false, false, false, false, false};
1147  unsigned int object_before_ip, object_after_ip;
1148 
1149  //unsigned int last_triangle_vertex=30; // no vertex at last IP
1150 
1151  // pass through the ccwise oriented sides in ccwise oriented order
1152  // How to make this in independent way?
1153  // Move this into RefElement?
1154  std::vector<unsigned int> side_cycle_orientation = { 0, 0, 1};
1155  std::vector<unsigned int> cycle_sides = {0, 2, 1};
1156 
1157  // TODO:
1158  // better mechanism for detecting vertex duplicities, do no depend on cyclic order of sides
1159  // still need cyclic orientation
1160 
1161  /**
1162  * Compute triangle side - tetrahedron intersections
1163  * IPs order is given by orientation of triangle sides, sign of IP results influenced by S3 inversion
1164  * is irrelevant.
1165  */
1166  for(unsigned int _i_side = 0; _i_side < RefElement<2>::n_lines; _i_side++) { // go through triangle lines
1167  unsigned int i_side = cycle_sides[_i_side];
1168  IP13s.clear();
1169  CI13[ i_side ].compute(IP13s);
1170  ASSERT_DBG(IP13s.size() < 3);
1171  if (IP13s.size() == 0) continue;
1172  for(unsigned int _ip=0; _ip < IP13s.size(); _ip++) {
1173  //int ip_topo_position = _ip*2-1; // -1 (incoming ip), +1 (outcoming ip), 0 both
1174 
1175  // fix order of IPs
1176  unsigned int ip = (side_cycle_orientation[_i_side] + _ip) % IP13s.size();
1177 
1178 // DebugOut().fmt("rside: {} cside: {} rip: {} cip: {}", _i_side, i_side, _ip, ip);
1179 
1180  // convert from 13 to 23 IP
1181  IntersectionPointAux<3,1> IP31 = IP13s[ip].switch_objects(); // switch idx_A and idx_B and coords
1182  IntersectionPointAux<3,2> IP32(IP31, i_side); // interpolation uses local_bcoords_B and given idx_B
1183  IPAux23 IP23 = IP32.switch_objects(); // switch idx_A and idx_B and coords back
1184  //DebugOut() << IP;
1185 
1186  // Tracking info
1187  unsigned int tetra_object = s3_dim_starts[IP23.dim_B()] + IP23.idx_B();
1188  unsigned int side_object = s2_dim_starts[1] + i_side;
1189 
1190  object_before_ip = tetra_object;
1191  object_after_ip = side_object;
1192 
1193 
1194  // IP is vertex of triangle,
1195  if( IP23.dim_A() == 0 && IP23.dim_B() == 3)
1196  {
1197  // we are on line of the triangle, and IP.idx_A contains local node of the line
1198  // E-E, we know vertex index
1199  object_before_ip = s2_dim_starts[0] + IP23.idx_A();
1200  }// else current_triangle_vertex=3+IP23_list.size(); // no vertex, and unique
1201 
1202  // side of triangle touching S3, in vertex or in edge
1203  if (IP13s.size() == 1 ) {
1204  if (IP23.dim_B() == 0) {
1205  continue; // skip, S3 vertices are better detected in phase 2
1206  }
1207  if (IP23.dim_A() == 0) { // vertex of triangle
1208  object_before_ip = tetra_object;
1209  object_after_ip = s2_dim_starts[0]+IP23.idx_A();
1210 
1211  // source vertex of the side vector (oriented CCwise)
1212  //if ( (IP.idx_A()+side_cycle_orientation[_i_side])%2 == 0)
1213  // std::swap(object_before_ip, object_after_ip);
1214  } else {
1215  // touch in edge
1216 
1217  //continue;
1218  ASSERT_EQ_DBG(IP23.dim_B(), 1);
1219  edge_touch[IP23.idx_B()]=true;
1220  std::swap(object_before_ip, object_after_ip);
1221  }
1222  }
1223 
1224  unsigned int ip_idx = add_ip(IP23);
1225 
1226 // DebugOut().fmt("before: {} after: {} ip: {}\n", object_before_ip, object_after_ip, ip_idx);
1227  ASSERT_EQ_DBG(IP23_list.size(), IP_next.size()+1);
1228  set_links(object_before_ip, ip_idx, object_after_ip);
1229  }
1230 
1231  }
1232 
1233  // now we have at most single true degenerate IP in IP23
1234  // TODO:
1235  // - deal with degenerate IPs in the edge-trinagle phase
1236  // - remove degenerate_list (just count degen points)
1237  // - remove check for duplicities in final list copy
1238  // - add more comment to individual cases in order to be sure that any case in particular branch is
1239  // treated right
1240 
1241 
1242  //TODO: cannot the two cycles be merged in one now?
1243  // and instead of processed_edge use CI12[tetra_edge].set_computed();
1244  // and I think we don't have to generate dummy IP12s
1245  // so we don't have to have IP12s vector here
1246  IP12s_.clear();
1247  // S3 Edge - S2 intersections; collect all signs, make dummy intersections
1248  for(unsigned int tetra_edge = 0; tetra_edge < 6; tetra_edge++) {
1249  IPAux12 IP12;
1250  CI12[tetra_edge].compute(IP12);
1251  // IntersectionResult result = CI12[tetra_edge].compute(IP12);
1252  // DebugOut() << print_var(tetra_edge) << print_var(int(result));
1253  // in degenerate case: IP12 is empty with degenerate result
1254  IP12s_.push_back(IP12);
1255  }
1256  vector<uint> processed_edge(6, 0);
1257  FacePair face_pair;
1258  for(unsigned int tetra_edge = 0; tetra_edge < 6;tetra_edge++) {
1259  if (! processed_edge[tetra_edge]) {
1260 // DBGVAR(tetra_edge);
1261  IPAux12 &IP12 = IP12s_[tetra_edge];
1262  if(IP12.result() >= IntersectionResult::degenerate) continue;
1263 
1264  int sign = CI12[tetra_edge].check_abscissa_topology(IP12);
1265 // DBGVAR(sign);
1266  if(std::abs(sign) > 1) continue;
1267 
1268  IPAux23 IP23(IP12s_[tetra_edge].switch_objects(), tetra_edge);
1269 
1270  const uint edge_dim = IP23.dim_B();
1271  const uint i_edge = IP23.idx_B();
1272  ASSERT_LT_DBG(edge_dim, 2);
1273 
1274  unsigned int ip_idx = add_ip(IP23);
1275 // DebugOut() << print_var(edge_dim);
1276 
1277  if ( edge_dim == 0) {
1278  face_pair = vertex_faces(i_edge);
1279  // mark edges coincident with the vertex
1280  for( uint ie : RefElement<3>::interact(Interaction<1,0>(i_edge), S3_inverted) )
1281  processed_edge[ie] = 1;
1282  }
1283  else
1284  face_pair = edge_faces(i_edge);
1285 
1286 // DebugOut() << print_var(face_pair[0])<< print_var(face_pair[1]);
1287 
1288 
1289 
1290  unsigned int s3_object = s3_dim_starts[edge_dim] + i_edge;
1291 
1292 // DebugOut() << print_var(edge_touch[i_edge]) << print_var(have_backlink(s3_object));
1293  if (IP23.dim_A() < 2
1294  && (! edge_touch[i_edge])
1295  && object_next[s3_object] != no_idx) { // boundary of S2, these ICs are duplicate
1296 
1297  if ( have_backlink(s3_object) ) {
1298  set_links(s3_object, ip_idx, face_pair[1]);
1299  } else {
1300  set_links(face_pair[0], ip_idx, s3_object);
1301  }
1302 
1303  } else { // interior of S2, just use the face pair
1304 // DebugOut() << print_var(face_pair[0])<< print_var(face_pair[1]);
1305  set_links(face_pair[0], ip_idx, face_pair[1]);
1306 
1307  if ( have_backlink(s3_object) ) {
1308  object_next[s3_object]=ip_idx;
1309  }
1310  }
1311  }
1312  }
1313 
1314 
1315  // Return IPs in correct order and remove duplicates
1316  ASSERT_EQ(0, intersection.size());
1317 
1318  if (IP23_list.size() == 0) return; // empty intersection
1319 
1320  // detect first IP, this needs to be done only in the case of
1321  // point or line intersections, where IPs links do not form closed cycle
1322  // Possibly we do this only if we detect such case through previous phases.
1323  vector<char> have_predecessor(IP23_list.size(), 0);
1324  for(auto obj : IP_next) {
1325  ASSERT_LT_DBG(obj, object_next.size());
1326  unsigned int ip = object_next[obj];
1327  if (ip < IP_next.size()) have_predecessor[ip]=1;
1328  }
1329  unsigned int ip_init=0;
1330  for(unsigned int i=0; i< IP23_list.size(); i++) if (! have_predecessor[i]) ip_init=i;
1331 
1332  // count number of IPs on the S3 faces.
1333  arma::uvec::fixed<RefElement<3>::n_sides> ips_face_counter;
1334  ips_face_counter.zeros();
1335 
1336  // regular case, merge duplicit IPs
1337  unsigned int ip=ip_init;
1338  ASSERT_EQ_DBG(IP_next.size(), IP23_list.size());
1339  intersection.points().push_back(IP23_list[ip]);
1340  //DebugOut() << print_var(ip) << IP23_list[ip];
1341  while (1) {
1342  //DebugOut() << print_var(ip) << IP23_list[ip];
1343 
1344  unsigned int object = IP_next[ip];
1345  //IP_next[ip]=no_idx;
1346  ASSERT_LT_DBG(object, object_next.size());
1347  ip = object_next[object];
1348  object_next[object]=no_idx;
1349  if ((ip == no_idx)) break;
1350  ASSERT_LT_DBG(ip, IP_next.size());
1351 
1352  if ( ! IP23_list[ip].topology_equal(intersection.points().back()) ) {
1353  IPAux23 &IP = IP23_list[ip];
1354  //DebugOut() << print_var(ip) << IP23_list[ip];
1355  intersection.points().push_back(IP);
1356  if(IP.dim_B() < 3)
1357  ips_face_counter += on_faces[IP.dim_B()][IP.idx_B()];
1358  }
1359  }
1360 
1361  if (intersection.points().size() == 1) return;
1362 
1363  if (IP23_list[ip_init].topology_equal(intersection.points().back()) )
1364  intersection.points().pop_back();
1365 
1366  if (intersection.points().size() > 2){
1367  for(uint i=0; i< ips_face_counter.n_elem; i++){
1368  if(ips_face_counter(i) >= intersection.points().size()){ // all IPs lie in a single face
1369 // DBGCOUT(<< "all IPs in face: " << i <<"\n");
1370 // ips_face_counter.print("face_counter");
1371  intersection.set_ips_in_face(i);
1372  }
1373  }
1374  }
1375 }
1376 
1378 {
1379  auto &line_faces=RefElement<3>::interact(Interaction<2,1>(i_edge), S3_inverted);
1380  unsigned int ip_ori = (unsigned int)(IP12s_[i_edge].result());
1381  ASSERT_DBG(ip_ori < 2); // no degenerate case
1382  //ip_ori ^= (int)(mesh_->element_accessor(intersection_->bulk_ele_idx).inverted());
1383 
1384  // RefElement returns edge faces in clockwise order (edge pointing to us)
1385  // negative ip sign (ori 0) = faces counter-clockwise
1386  // positive ip sign (ori 1) = faces clockwise
1387  return { s3_dim_starts[2] + line_faces[1-ip_ori], s3_dim_starts[2] + line_faces[ip_ori] };
1388 }
1389 
1391 {
1392  // vertex edges clockwise when looking from outside of tetrahedron
1393  const IdxVector<3> &vtx_edges = RefElement<3>::interact(Interaction<1,0>(i_vertex), S3_inverted);
1394  std::array<unsigned int, 3> n_ori, sum_idx;
1395  n_ori.fill(0);
1396  sum_idx.fill(0);
1397  for(unsigned int ie=0; ie <3; ie++) {
1398  unsigned int edge_ip_ori = (unsigned int)(IP12s_[ vtx_edges[ie]].result());
1399  if (RefElement<3>::interact(Interaction<0,1>(vtx_edges[ie]))[0] != i_vertex
1400  && edge_ip_ori!= int(IntersectionResult::degenerate) )
1401  edge_ip_ori = (edge_ip_ori +1)%2;
1402  //ASSERT_LT_DBG(edge_ip_ori, 3)(ie);
1403  if (edge_ip_ori == 3) edge_ip_ori=2; // none as degenerate
1404  n_ori[edge_ip_ori]++;
1405  sum_idx[edge_ip_ori]+=ie;
1406  }
1407  unsigned int n_degen = n_ori[ int(IntersectionResult::degenerate) ];
1408  unsigned int sum_degen = sum_idx[ int(IntersectionResult::degenerate) ];
1409  unsigned int n_positive = n_ori[ int(IntersectionResult::positive) ];
1410  unsigned int n_negative= n_ori[ int(IntersectionResult::negative) ];
1411 // DebugOut().fmt("nd: {} sd: {} np: {} nn: {}", n_degen, sum_degen, n_positive, n_negative);
1412  if ( n_degen == 2 ) {
1413  // S2 plane match a face of S3, we treat degenerated edges as the faces
1414  // incident with the single regualr edge.
1415 
1416  unsigned int i_edge = 3 - sum_degen; // regular edge index
1417  FacePair pair = edge_faces(vtx_edges[i_edge]);
1418  auto &vtx_faces = RefElement<3>::interact(Interaction<2,0>(i_vertex), S3_inverted);
1419  // replace faces by edges
1420  if (pair[0] == s3_dim_starts[2] + vtx_faces[(i_edge+1)%3])
1421  return { s3_dim_starts[1] + (i_edge+2)%3, s3_dim_starts[1] + (i_edge+1)%3 };
1422  else
1423  return { s3_dim_starts[1] + (i_edge+1)%3, s3_dim_starts[1] + (i_edge+2)%3 };
1424 
1425  } else if (n_degen == 1) {
1426  // One edge in S2 plane.
1427  unsigned int i_edge = sum_degen;
1428  ASSERT_DBG( n_positive + n_negative == 2);
1429  if ( n_positive == 1) {
1430  // opposite signs, S2 plane cuts S3
1431 
1432  FacePair pair = edge_faces(vtx_edges[(i_edge+1)%3]);
1433  // Get the face opposite to the degenerated edge.
1434  unsigned int face = RefElement<3>::interact(Interaction<2,0>(i_vertex), S3_inverted)[i_edge];
1435  // assign edges to faces
1436  //DebugOut().fmt("vtx: {} edg: {} face: {}", i_vertex, i_edge, face);
1437  if (pair[0] == s3_dim_starts[2] + face)
1438  return { s3_dim_starts[2] + face, s3_dim_starts[1] + vtx_edges[i_edge]};
1439  else
1440  return { s3_dim_starts[1] + vtx_edges[i_edge], s3_dim_starts[2] + face };
1441  } else {
1442  // same signs; S2 plane touch S3 vertex and a single edge
1443  //DebugOut() << "Touch in edge.";
1444  // same signs S2 plane touchs S3
1445  ASSERT_DBG(n_positive == 0 || n_positive== 2);
1446  return { s3_dim_starts[0]+i_vertex, s3_dim_starts[1] + vtx_edges[i_edge]};
1447  }
1448 
1449 
1450  } else {
1451  ASSERT_DBG(n_degen == 0);
1452  ASSERT_DBG( n_positive + n_negative == 3);
1453 
1454  if (n_positive == 1) {
1455  unsigned int i_edge = sum_idx[ int(IntersectionResult::positive) ];
1456  return edge_faces(vtx_edges[i_edge]);
1457  } else if (n_negative == 1) {
1458  unsigned int i_edge = sum_idx[ int(IntersectionResult::negative) ];
1459  return edge_faces(vtx_edges[i_edge]);
1460  } else {
1461  // S2 touch vertex of S3 in
1462  ASSERT_DBG( n_positive == 0 || n_positive == 3);
1463  return { s3_dim_starts[0]+i_vertex, s3_dim_starts[0]+i_vertex};
1464  }
1465  }
1466 }
1467 
1468 
1470 {
1472 
1473  on_faces.resize(3);
1474  arma::uvec::fixed<RefElement<3>::n_sides> v; v.zeros();
1475 
1476  on_faces[0].resize(RefElement<3>::n_nodes, v);
1477  for(uint i=0; i<RefElement<3>::n_nodes; i++) {
1478  auto faces = RefElement<3>::interact(Interaction<2,0>(i), false); // order doesn't matter
1479  for(uint j=0; j<RefElement<3>::n_sides_per_node; j++)
1480  on_faces[0][i](faces[j]) = 1;
1481  }
1482 
1483  on_faces[1].resize(RefElement<3>::n_lines, v);
1484  for(uint i=0; i<RefElement<3>::n_lines; i++) {
1485  auto faces = RefElement<3>::interact(Interaction<2,1>(i), false); // order doesn't matter
1486  for(uint j=0; j<RefElement<3>::n_sides_per_line; j++)
1487  on_faces[1][i](faces[j]) = 1;
1488  }
1489 
1490  on_faces[2].resize(RefElement<3>::n_sides, v);
1491  for(uint i=0; i<RefElement<3>::n_sides; i++) {
1492  on_faces[2][i](i) = 1;
1493  }
1494 
1495 // DBGCOUT("Print on_faces:\n");
1496 // for(uint d=0; d<on_faces.size(); d++)
1497 // for(uint i=0; i<on_faces[d].size(); i++){
1498 // for(uint j=0; j<on_faces[d][i].n_elem; j++)
1499 // cout << on_faces[d][i](j) << " ";
1500 // cout << endl;
1501 // }
1502  return on_faces;
1503 }
1504 
1505 
1507  for(unsigned int i = 0; i < 3;i++){
1508  os << "\tPluckerCoordinates Triangle[" << i << "]";
1509  if(plucker_coordinates_triangle_[i] == nullptr){
1510  os << "NULL" << endl;
1511  }else{
1512  os << plucker_coordinates_triangle_[i];
1513  }
1514  }
1515  for(unsigned int i = 0; i < 6;i++){
1516  os << "\tPluckerCoordinates Tetrahedron[" << i << "]";
1517  if(plucker_coordinates_tetrahedron[i] == nullptr){
1518  os << "NULL" << endl;
1519  }else{
1520  os << plucker_coordinates_tetrahedron[i];
1521  }
1522  }
1523 }
1524 
1526  os << "ComputeIntersection<2,3> Plucker Coordinates Tree:" << endl;
1527  print_plucker_coordinates(os);
1528  for(unsigned int i = 0; i < 6;i++){
1529  os << "ComputeIntersection<1,2>["<< i <<"] Plucker Coordinates:" << endl;
1530  CI12[i].print_plucker_coordinates(os);
1531  }
1532  for(unsigned int i = 0; i < 3;i++){
1533  CI13[i].print_plucker_coordinates_tree(os);
1534  }
1535 }
1536 
IntersectionPointAux::set_topology_B
void set_topology_B(unsigned int idx, unsigned int dim_B)
Sets the topology of object B in Simplex<M>.
Definition: intersection_point_aux.hh:234
ref_element.hh
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc.
Plucker
Plucker coordinates representing line given by points A,B.
Definition: plucker.hh:45
RefElement::line_barycentric_interpolation
static BaryPoint line_barycentric_interpolation(BaryPoint first_coords, BaryPoint second_coords, double first_theta, double second_theta, double theta)
Definition: ref_element.cc:437
RefElement
Definition: ref_element.hh:339
RefElement::normal_orientation
static unsigned int normal_orientation(unsigned int sid)
Definition: ref_element.cc:216
std::swap
void swap(nlohmann::json &j1, nlohmann::json &j2) noexcept(is_nothrow_move_constructible< nlohmann::json >::value and is_nothrow_move_assignable< nlohmann::json >::value)
exchanges the values of two JSON objects
Definition: json.hpp:8688
IntersectionAux::i_points_
std::vector< IntersectionPointAux< dimA, dimB > > i_points_
Vector of internal intersection points.
Definition: intersection_aux.hh:39
RefElement::interact
static const IdxVector<(InDim >OutDim ? InDim+1 :dim-InDim) > interact(TInteraction< OutDim, InDim > interaction, bool inv=false)
Element::dim
unsigned int dim() const
Definition: elements.h:120
IntersectionPointAux::set_result
void set_result(IntersectionResult result)
Setter orientation flag.
Definition: intersection_point_aux.hh:209
ASSERT_DBG
#define ASSERT_DBG(expr)
Definition: include_fadbad.hh:28
geometry_epsilon
static const double geometry_epsilon
Definition: intersection_point_aux.hh:28
std::vector
Definition: doxy_dummy_defs.hh:7
ElementAccessor< 3 >
RefElement::zeros_positions
static std::pair< unsigned int, unsigned int > zeros_positions(const BaryPoint &barycentric, double tolerance=std::numeric_limits< double >::epsilon() *2)
Definition: ref_element.cc:294
system.hh
arma::vec3
Definition: doxy_dummy_defs.hh:17
ASSERT_EQ_DBG
#define ASSERT_EQ_DBG(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Definition: asserts.hh:332
uint
unsigned int uint
Definition: mh_dofhandler.hh:101
IntersectionAux::set_ips_in_face
void set_ips_in_face(unsigned int face_idx)
Definition: intersection_aux.hh:138
IntersectionPointAux::topology_equal
bool topology_equal(const IntersectionPointAux< N, M > &other) const
Returns true, if other intersection point has the same topology.
Definition: intersection_point_aux.cc:133
IntersectionAux::points
std::vector< IntersectionPointAux< dimA, dimB > > & points()
Returns intersection points by a reference.
Definition: intersection_aux.hh:100
plucker.hh
Plucker coordinates class.
IntersectionPointAux::local_bcoords_A
const arma::vec::fixed< N+1 > & local_bcoords_A() const
Returns barycentric coordinates in the Simplex<N>.
Definition: intersection_point_aux.hh:221
accessors.hh
IntersectionAux::bulk_ele_idx
unsigned int bulk_ele_idx() const
Returns index of bulk element.
Definition: intersection_aux.hh:122
RefElement::bary_to_local
static LocalPoint bary_to_local(const BaryPoint &bp)
Converts from barycentric to local coordinates.
Definition: ref_element.cc:201
IntersectionPointAux::set_topology_A
void set_topology_A(unsigned int idx, unsigned int dim_A)
Sets the topology of object A in Simplex<N>.
Definition: intersection_point_aux.hh:229
IntersectionPointAux::idx_A
unsigned int idx_A() const
Returns the index of Simplex<N>.
Definition: intersection_point_aux.hh:239
IdxVector
std::array< unsigned int, Size > IdxVector
Definition: ref_element.hh:276
ElementAccessor::sign
double sign() const
Definition: accessors.hh:124
ElementAccessor::jacobian_S3
double jacobian_S3() const
Definition: accessors.hh:132
ASSERT_EQ
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual)
Definition: asserts.hh:328
IntersectionPointAux
Internal auxiliary class represents an intersection point of simplex<N> and simplex<M>.
Definition: compute_intersection.hh:50
IntersectionPointAux::dim_A
unsigned int dim_A() const
Returns dimension of object A.
Definition: intersection_point_aux.hh:213
IntersectionResult
IntersectionResult
Definition: intersection_point_aux.hh:38
mesh.h
ElementAccessor::node
NodeAccessor< 3 > node(unsigned int ni) const
Definition: accessors.hh:233
ElementAccessor::input_id
unsigned int input_id() const
Return the element ID in the input mesh. Should be only used for special output.
Definition: accessors.hh:223
Interaction
Definition: ref_element.hh:286
intersection_aux.hh
Internal class representing intersection object.
IntersectionPointAux::switch_objects
IntersectionPointAux< M, N > switch_objects() const
Switches the object A and B.
Definition: intersection_point_aux.cc:101
Mesh
Definition: mesh.h:355
IntersectionAux::size
unsigned int size() const
Returns number of intersection points.
Definition: intersection_aux.hh:114
IntersectionAux< 2, 2 >
std
Definition: doxy_dummy_defs.hh:5
compute_intersection.hh
Fundamental simplicial intersections.
IntersectionPointAux::set_coordinates
void set_coordinates(const arma::vec::fixed< N+1 > &lcA, const arma::vec::fixed< M+1 > &lcB)
Setter for coordinates.
Definition: intersection_point_aux.hh:196
IntersectionResult::positive
@ positive
ASSERT_LT_DBG
#define ASSERT_LT_DBG(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Definition: asserts.hh:300
IntersectionResult::negative
@ negative
IntersectionResult::degenerate
@ degenerate
IntersectionPointAux::idx_B
unsigned int idx_B() const
Returns the index of Simplex<M>.
Definition: intersection_point_aux.hh:243
ComputeIntersection
Definition: compute_intersection.hh:47
IntersectionPointAux::result
IntersectionResult result() const
Result: 0 - negative sign, 1 - positive sign, 2 - degenerate (zero for all sides),...
Definition: intersection_point_aux.hh:247
plucker_empty
static const double plucker_empty
Auxiliary value for Plucker product. If equal this value, it is supposed not to be set yet.
Definition: compute_intersection.hh:53
IntersectionPointAux::dim_B
unsigned int dim_B() const
Returns dimension of object B.
Definition: intersection_point_aux.hh:217
IntersectionResult::none
@ none
ComputeIntersection< 2, 3 >::FacePair
std::array< uint, 2 > FacePair
Definition: compute_intersection.hh:514
intersection_point_aux.hh
Internal class representing intersection point.
FMT_UNUSED
#define FMT_UNUSED
Definition: posix.h:75