openCARP
Doxygen code documentation for the open cardiac electrophysiology simulator openCARP
filament.cc
Go to the documentation of this file.
1 /*
2  * Filament detection routines
3  */
4 
5 #include "FMatrix.h"
6 #include "filament.h"
7 
8 #ifdef I
9 #undef I
10 #endif
11 
12 // this table is the mass matrix for the 1st subtriangle of the quadrilateral
13 // the rest of the mass matrices are just permutations
14 static const int TriFaces[1][3] = {
15  { 0, 1, 2 }
16 };
17 
18 static const int TetFaces[4][3] = {
19  { 0, 1, 2 },
20  { 0, 1, 3 },
21  { 1, 2, 3 },
22  { 0, 2, 3 }
23 };
24 
25 
26 #define FIL_PTS_FILE 0
27 #define FIL_ELEM_FILE 1
28 #define FIL_DATA_FILE 2
29 
30 const int MAX_ELEM_NODES = 8;
31 const int MAX_ELEM_FACES = 6;
32 
33 // prototypes
34 void output_filaments(filament *f, FILE *fh[]);
35 
36 bool fl_chkIso(float *d, float th, Elem *e);
37 bool fl_IsoTagFaces (float *d,float th, Elem *e, bool *faceTags );
38 bool fl_IsoTagTetFaces(float *d,float th, bool *faceTags);
39 bool fl_IsoTagTriFaces(float *d,float th, bool *faceTags);
40 bool fl_IsoTagTriangle(float *d, float th, const int idx[]);
41 bool fl_findIsoIntersect(float *d0, float *d1, float th, Point *p, Elem *e );
42 int fl_FindFilament(float *d0, float *d1, float th, Point *lp, Elem *e, singularity *sng);
43 bool fl_triFindPS(float *d0, float *d1, float th, Point *p, Point *PS,int meth);
44 bool fl_triFindPS_Lines(float *d0, float *d1, float th, Point *p, Point *PS);
45 bool fl_triFindPS_ShpFnc(float *d0,float *d1,float th,Point *p,Point *PS);
46 bool fl_findEdgeSection(float d_p0, float d_p1, float th, Point *p0, Point *p1, Point *pSec);
47 void interpoate_vm(int numNodes,float **d,int M,float **rd,float rDist);
48 void interpolateTBuffer(tbuffer *tbf);
49 float fl_triArea( Point *pts );
50 
51 // test point in triangle
52 bool SameSideOfLine (Point pTest, Point p0, Point p1, Point p2);
53 bool PointInTriangle(Point pTest, Point p0, Point p1, Point p2);
54 bool PointInTriangleBarycentric(Point pTest, Point p0, Point p1, Point p2);
55 bool PointInTriangleBarycentric(Point pTest, Point p0, Point p1, Point p2, Point &bcc);
56 bool PtInTriangleBarycentricSP(Point pTest, Point p1, Point p2);
57 bool PtInTriangleBarycentricSP(Point pTest, Point p1, Point p2, Point &bcc);
58 
59 
68 Point
69 UnTransformPoint( Point PTs, Point *opts, Point *spts)
70 {
71  Real whole = fabs(spts[1].x*spts[2].y/2.);
72  Point a[3] = { PTs, spts[1], spts[2] };
73  Real c0 = fl_triArea( a )/whole;
74  Real c2 = fabs(spts[1].x*PTs.y/2.)/whole;
75  if( c0+c2 > 1.02 ){
76  fprintf( stderr, "transformation back failed\n" );
77  }
78  return opts[0]*c0 + opts[2]*c2 + opts[1]*(1-c0-c2);
79 }
80 
81 
91 inline float
93 {
94  return 0.5*mag(cross(pts[1] - pts[0],pts[2] - pts[0]));
95 }
96 
107 Point
108 fl_triNormalLocal( Point *pts, bool u )
109 {
110  Point n=cross(pts[1] - pts[0],pts[2] - pts[0]);
111  return u?normalize(n):n;
112 }
113 
114 
123 void
125 {
126  int numNodes = tbf->nodes;
127  float **d = (float **)(tbf->d );
128  float **rd = (float **)(tbf->dr);
129  int M = tbf->wdth;
130  float rdist = tbf->rdist;
131 
132  interpolate_vm(numNodes,d,M,rd,rdist);
133 }
134 
135 
146 void
147 interpolate_vm(int numNodes,float **d,int M,float **rd,float rdist)
148 {
149  for(int j=0;j<numNodes;j++) {
150  rd[0][j] = d[0 ][j]+rdist*(d[1 ][j]-d[0 ][j]);
151  rd[1][j] = d[M-2][j]+rdist*(d[M-1][j]-d[M-2][j]);
152  }
153 }
154 
155 
170 bool
171 compute_filament(float *d0, float *d1, float vth, ElemList *elst,
172  NodeList *nlst, filament *f, SingularityDetector_t meth)
173 {
174  int PScnt = 0;
175 
176 #pragma omp parallel for
177  for(int i=0;i<elst->NElems; i++) {
178 
179  float d0_elem[MAX_ELEM_NODES];
180  float d1_elem[MAX_ELEM_NODES];
181 
182  Elem &e = elst->Ele[i];
183 
184  if((e.type!=Tri) && (e.type!=Tetra))
185  log_msg(NULL,2,0, "%s: Type in Element %d not supported yet.\n",
186  __func__, i );
187 
188  // copy data at element nodes to local array
189  for(int j=0;j<e.Nnode;j++) {
190  d0_elem[j] = d0[e.N[j]];
191  d1_elem[j] = d1[e.N[j]];
192  }
193 
194  // quick check, are there any isolines/surfaces on the element
195  bool iso0 = fl_chkIso(d0_elem,vth,&e);
196  bool iso1 = fl_chkIso(d1_elem,vth,&e);
197 
198  // if so, analyse in detail
199  if( iso0 && iso1 ) {
200 
201  // possibly a filament, at least we have two isosurfaces/isolines
202  // is this actually a filament?
203  bool ftags0[MAX_ELEM_FACES] = {false};
204  bool ftags1[MAX_ELEM_FACES] = {false};
205 
206  // all faces with with isolines tagged
207  bool tags0 = fl_IsoTagFaces(d0_elem, vth, &e, ftags0 );
208  bool tags1 = fl_IsoTagFaces(d1_elem, vth, &e, ftags1 );
209 
210  Point lp[MAX_ELEM_NODES];
211  for( int j=0; j<e.Nnode; j++ )
212  lp[j] = nlst->Node[e.N[j]];
213 
214  if(tags0 && tags1) {
215 
216  if((e.type==Tri) && (meth==vol_based)) {
217  lp[3] = lp[0] + fl_triNormalLocal(lp,false);
218  d0_elem[3] = d0_elem[0];
219  d1_elem[3] = d1_elem[0];
220  }
221 
222  singularity sng;
223  int FilamentType = fl_FindFilament(d0_elem,d1_elem,vth,lp,&e,&sng);
224 
225  if(FilamentType) {
226 #pragma omp critical
227  f->add_singularity(&sng,i);
228 #pragma omp atomic
229  PScnt++;
230  }
231 
232  //bool isFilament = fl_findIsoIntersect(d0_elem,d1_elem,vth,lp,e);
233  //if(isFilament)
234  // add_singularity2filament(f, elst->ElemGlobInds[i] );
235 
236  }
237  }
238  }
239  f->sort();
240  return PScnt>0;
241 }
242 
243 void
245 {
246  s.push_back(*sng);
247  numPts += sng->ps.type;
248  s.back().ps.eIndx = eIndx;
249 }
250 
251 
252 void
254 {
255  std::sort(s.begin(), s.end(),
256  [] (const singularity& lhs, const singularity& rhs)-> bool{return lhs.ps.eIndx<rhs.ps.eIndx;});
257 }
258 
259 
269 int
270 init_filament_IO(const char *meshname, int numT, FILE *fh[])
271 {
272  int ok = -1;
273 
274  char pname[512], ename[512], dname[512];
275 
276  snprintf(pname, sizeof pname, "%s.pts_t", meshname);
277  snprintf(ename, sizeof ename, "%s.elem_t", meshname);
278  snprintf(dname, sizeof dname, "%s.dat_t", meshname);
279 
280  static int init = 0;
281  if(!init) {
282 
283  FILE *fhp = fopen(pname,"wt"); fh[0] = fhp;
284  FILE *fhe = fopen(ename,"wt"); fh[1] = fhe;
285  FILE *fhd = fopen(dname,"wt"); fh[2] = fhd;
286 
287  if( (fhp==NULL) || (fhe==NULL) || (fhd==NULL) )
288  return ok;
289  else
290  ok = 0;
291 
292  // write comment lines
293  fprintf(fhp,"## Format x y z # singularity type (0=PS || 1=FilSeg) Global element index \n");
294  fprintf(fhp,"%d\n", numT);
295  fprintf(fhe,"## Format Ln p0 p1 || PS p0 \n");
296  fprintf(fhe,"%d\n", numT);
297  fprintf(fhd,"## Format data (matching with pts_t file) \n");
298  fprintf(fhd,"%d\n", numT);
299 
300  init = 1;
301  }
302 
303  return ok;
304 }
305 
306 
307 // temporary debug code
308 int
309 write_connections(ElemList *elst, NodeList *nlst, filament *f, FILE *filHdls[],
310  float t)
311 {
312  static int init = 0;
313 
314  if(!f->num_segs()) return 0;
315 
316  char fname[512];
317  char pname[512], ename[512], dname[512];
318 
319  snprintf(fname, sizeof fname, "filament_%.2f.cnnx", t);
320  FILE *fh = fopen(fname,"wt");
321  fprintf( fh,"%d\n", f->num_segs());
322  int i = 0;
323 
324  for( int i=0; i<f->num_segs(); i++ ) {
325  // pick the first edge of an element
326  Elem &e = elst->Ele[f->s[i].ps.eIndx];
327  //fl_getElem(elst,f->s[i].ps.eIndx,&e);
328  int vtx0 = e.N[0];
329  int vtx1 = e.N[1];
330 
331  fprintf( fh,"%d %d\n", vtx0, vtx1 );
332  }
333  fclose(fh);
334 
335  return f->num_segs();
336 }
337 
338 
347 void
348 write_filaments(filament *f,FILE *filHdls[],float t)
349 {
350  int msg1 = 20;
351  int msg2 = 30;
352  int N = f->num_segs();
353  int numPts = f->num_pts();
354 
355  for(int i=0; i<=FIL_DATA_FILE;i++)
356  fprintf(filHdls[i],"#%.6f\n",t);
357 
358  // preliminary, this needs to be fixed
359  fprintf(filHdls[FIL_PTS_FILE ],"%d\n",numPts);
360  //fprintf(filHdls[FIL_ELEM_FILE],"%d\n",N);
361  fprintf(filHdls[FIL_DATA_FILE],"%d\n",numPts);
362 
363  fprintf(stderr, "%d filament segments found\n",N);
364  output_filaments(f,filHdls);
365 }
366 
367 
374 void
375 output_filaments(filament *f, FILE *filHdls[])
376 {
377  int nele=0;
378  for( int i=0; i<f->num_segs(); i++ )
379  if( f->s[i].ps.type != PhaseSingularity ) nele++;
380  fprintf(filHdls[FIL_ELEM_FILE],"%d\n",nele);
381 
382  int ptCnt = 0;
383  for( int i=0; i<f->num_segs(); i++ ){
384  FILE *fhp = filHdls[0], *fhe = filHdls[1], *fhd = filHdls[2];
385 
386  FilSeg fs = f->s[i].fs;
387  fprintf(fhp,"%.3f %.3f %.3f # %d %d\n", fs.p0.x,fs.p0.y,fs.p0.z,
388  fs.type, fs.eIndx);
389  // for now we use the enclosing element index as data
390  // later we use a PS or filament ID
391  fprintf(fhd,"%d\n",fs.eIndx);
392 
393  if(f->s[i].ps.type==PhaseSingularity) {
394  //fprintf(fhe,"PS %d \n", ptCnt++);
395  ptCnt++;
396  }
397  else {
398  // filament segment
399  fprintf(fhp,"%.3f %.3f %.3f # %d %d\n", fs.p1.x,fs.p1.y,fs.p1.z,
400  fs.type, fs.eIndx);
401  fprintf(fhd,"%d\n",fs.eIndx);
402  fprintf(fhe,"Ln %d %d\n", ptCnt, ptCnt+1);
403  ptCnt += 2;
404  }
405  }
406 }
407 
408 
413 void
414 close_filament_IO(FILE *fh[])
415 {
416  for(int i=0;i<3;i++)
417  fclose(fh[i]);
418 }
419 
428 bool
429 fl_chkIso(float *d, float th, Elem *e)
430 {
431  if(d[0]>=th) {
432  for(int i=1;i<e->Nnode;i++)
433  if(d[i]<th) {
434  return true;
435  }
436  } else {
437  for(int i=1;i<e->Nnode;i++)
438  if(d[i]>=th) {
439  return true;
440  }
441  }
442  return false;
443 }
444 
445 
455 bool
456 fl_IsoTagFaces(float *d, float th, Elem *e, bool *faceTags )
457 {
458  bool flg = false;
459  switch(e->type) {
460  case Tri:
461  flg = fl_IsoTagTriFaces(d,th,faceTags);
462  break;
463  case Tetra:
464  flg = fl_IsoTagTetFaces(d,th,faceTags);
465  break;
466  default:
467  fprintf(stderr,"%s: Element type not implemented yet.\n", __func__);
468  }
469  return flg;
470 }
471 
472 
482 bool
483 fl_IsoTagTetFaces(float *d,float th, bool *faceTags)
484 {
485  int tagged = 0;
486  for(int i=0;i<4;i++)
487  if((faceTags[i]=fl_IsoTagTriangle(d,th,TetFaces[i])))
488  tagged++;
489 
490  return tagged;
491 }
492 
493 
503 bool
504 fl_IsoTagTriFaces(float *d,float th, bool *faceTags)
505 {
506  return faceTags[0]=fl_IsoTagTriangle(d,th,TriFaces[0]);
507 }
508 
509 
517 bool
518 fl_IsoTagTriangle(float *d, float th, const int idx[])
519 {
520  bool tag = false;
521 
522  if(d[idx[0]]>=th) {
523  if((d[idx[1]]<th) || (d[idx[2]]<th)) tag = true;
524  }
525  else
526  if((d[idx[1]]>=th) || (d[idx[2]]>=th)) tag = true;
527  return tag;
528 }
529 
530 
542 int
543 fl_FindFilament(float *d0,float *d1,float th,Point *lp,Elem *e,singularity *sng)
544 {
545  float d0_[3], d1_[3];
546  Point lp_[3];
547  int FilamentType = 0;
548  Point PhaseS[4];
549  int tagged[] = { 0, 0, 0, 0};
550  int tagIdx[] = {-1,-1,-1,-1};
551 
552  int meth = 1;
553  switch(e->type) {
554  case Tri:
555  if(fl_triFindPS(d0,d1,th,lp,PhaseS,meth)) {
556  int i = 0;
557  tagged[i]++;
558  tagIdx[FilamentType++] = i;
559  }
560  break;
561  case Tetra:
562  for(int i=0;i<4;i++) {
563  for(int j=0;j<3;j++) {
564  d0_[j] = d0[TetFaces[i][j]];
565  d1_[j] = d1[TetFaces[i][j]];
566  lp_[j] = lp[TetFaces[i][j]];
567  }
568  if(fl_triFindPS(d0_,d1_,th,lp_,PhaseS+i,meth)) {
569  tagged[i]++;
570  tagIdx[FilamentType++] = i;
571  }
572  }
573  break;
574  default:
575  fprintf(stderr,"Element type not supported by %s.\n", __func__ );
576  }
577 
578  if(FilamentType) {
579  sng->ps.type = FilamentType;
580  sng->ps.eIndx = 0;
581  sng->ps.p0 = PhaseS[tagIdx[0]];
582  if(FilamentType==2) // Filament segment
583  sng->fs.p1 = PhaseS[tagIdx[1]];
584  }
585 
586  return FilamentType;
587 }
588 
589 
601 bool
602 fl_triFindPS(float *d0, float *d1, float th, Point *lp, Point *PhaseS,
603  int meth)
604 {
605  bool found = false;
606  switch(meth) {
607  case 0:
608  found = fl_triFindPS_ShpFnc(d0,d1,th,lp,PhaseS);
609  break;
610  case 1:
611  default:
612  found = fl_triFindPS_Lines(d0,d1,th,lp,PhaseS);
613  }
614  return found;
615 }
616 
617 
628 bool
629 fl_triFindPS_ShpFnc(float *d0, float *d1, float th, Point *lp, Point *PhaseS)
630 {
631  fprintf(stderr, "Shape function method not implemented yet.\n");
632  return false;
633 }
634 
646 bool
647 fl_triFindPS_Lines(float *d0, float *d1, float th, Point *lp, Point *PhaseS)
648 {
649  bool isPS = false;
650  Point pt[3];
651 
652  // transform matrix into standard position first
653  //static FMatrix fwd = ZERO_FMATRIX, bwd = ZERO_FMATRIX, A = ZERO_FMATRIX;
654 
655  FMatrix A(2,2);
656 
657  /*
658  fl_getTransform3Dto2D(lp,&fwd,&bwd,&t);
659  int numPts = 3;
660  Point t;
661  fl_AffineTransformPoints(lp,pt,numPts,&fwd,&t,true);
662  */
663 
664  // transform to standard position
665  pt[0].set(0,0,0);
666  Point p01 = lp[1] - lp[0];
667  pt[1].set(mag(p01), 0., 0.);
668  Point p02 = lp[2] - lp[0];
669  double ptx = dot(p02,p01)/pt[1].x;
670  pt[2].set(ptx, sqrt(mag2(p02)-ptx*ptx),0.);
671 
672  int tEdg0[] = {0,0,0};
673  int tEdg1[] = {0,0,0};
674  Point tSec0[3], tSec1[3];
675 
676  // do our thingie
677  tEdg0[0] = fl_findEdgeSection(d0[0],d0[1],th,pt+0,pt+1,tSec0+0);
678  tEdg0[1] = fl_findEdgeSection(d0[1],d0[2],th,pt+1,pt+2,tSec0+1);
679  tEdg0[2] = fl_findEdgeSection(d0[0],d0[2],th,pt+0,pt+2,tSec0+2);
680 
681  tEdg1[0] = fl_findEdgeSection(d1[0],d1[1],th,pt+0,pt+1,tSec1+0);
682  tEdg1[1] = fl_findEdgeSection(d1[1],d1[2],th,pt+1,pt+2,tSec1+1);
683  tEdg1[2] = fl_findEdgeSection(d1[0],d1[2],th,pt+0,pt+2,tSec1+2);
684 
685  // singularity found?
686  int t0 = 0;
687  int t1 = 0;
688  int t0s[2], t1s[2];
689  int c0 = 0;
690  int c1 = 0;
691  for(int i=0;i<3;i++) {
692  t0 += tEdg0[i];
693  t1 += tEdg1[i];
694  if(tEdg0[i]) t0s[c0++] = i;
695  if(tEdg1[i]) t1s[c1++] = i;
696  }
697  if((t0==2) && (t1==2)) {
698  // equation in parametrized form
699  // a0, b0: intersection points of isoline 0 with tri
700  // a1, b1: intersection points of isoline 1 with tri
701  //
702  // There are two lines:
703  // l0 = a0 + lambda*(b0-a0) (1)
704  // l1 = a1 + gamma *(b1-a1) (2)
705  //
706  // The intersection is found by solving for l0=l1.
707  // this yields:
708  //
709  // (b0-a0)*lambda - (b1-a1)*gamma = (a1-a0)
710  //
711  // which gives a solution for lambda and gamma.
712  // Either (1) or (2) can be used to find the intersection coordinates
713 
714  Point &a0 = tSec0[t0s[0]];
715  Point &b0 = tSec0[t0s[1]];
716  Point &a1 = tSec1[t1s[0]];
717  Point &b1 = tSec1[t1s[1]];
718 
719  Point b0ma0 = b0 - a0;
720  Point b1ma1 = a1 - b1;
721  Point a1ma0 = a1 - a0;
722 
723  Real b[2];
724  A.Ent[0][0] = b0ma0.x; A.Ent[0][1] = b1ma1.x;
725  A.Ent[1][0] = b0ma0.y; A.Ent[1][1] = b1ma1.y;
726  b[0] = a1ma0.x; b[1] = a1ma0.y;
727 
728  // find params lambda and gamma
729  A.LUD();
730  if(A.decomposed) {
731  // A is not singular
732  A.LUsolve( b );
733 
734  // find intersection
735  Point PSt = a0 + scal_X(b0ma0,b[0]);
736 
737  // check whether PS is within triangle
738  Point bcc;
739 
740  if(PtInTriangleBarycentricSP(PSt,pt[1],pt[2],bcc)) {
741  //fprintf(stderr, "PS transformed found at x = %.2f, y = %.2f.\n", PSt.x, PSt.y );
742 
743  // transform back
744  // fl_AffineTransformPoints(&PSt,PhaseS,1,&bwd,&t,false);
745  // PhaseS[0] = UnTransformPoint(PSt, lp, pt);
746 
747  PhaseS[0] = lp[0]*bcc.x + lp[1]*bcc.y + lp[2]*bcc.z;
748 
749  // fprintf(stderr, "PS real found at x = %.2f, y = %.2f.\n", PhaseS->x, PhaseS->y );
750  isPS = true;
751 
752  bool debug = false;
753  if(debug) {
754  // check back transform - debug only!
755  Point ptb[3];
756  //fl_AffineTransformPoints(pt,ptb,3,&bwd,&t,false);
757  for( int j=0; j<3; j++ )
758  ptb[j] = dot(bcc,pt[j]);
759  }
760  }
761  }
762  }
763  return isPS;
764 }
765 
766 
779 bool
781 {
782  Point cp1 = cross(p1 - p0,pTest - p0);
783  Point cp2 = cross(p1 - p0,p2 - p0);
784  return dot(cp1,cp2)>=0;
785 }
786 
787 
798 bool
800 {
801  if(SameSideOfLine(pTest,p0,p1,p2) &&
802  SameSideOfLine(pTest,p1,p2,p0) &&
803  SameSideOfLine(pTest,p2,p0,p1) )
804  return true;
805  else
806  return false;
807 }
808 
810 bool
812 {
813  Point null;
814  return PtInTriangleBarycentricSP(pTest, p1, p2, null);
815 }
816 
829 bool
831 {
832  // Compute dot products
833  double dot00 = v0.x*v0.x+v0.y*v0.y;
834  double dot01 = v0.x*v1.x;
835  double dot02 = v0.x*v2.x+v0.y*v2.y;
836  double dot11 = v1.x*v1.x;
837  double dot12 = v1.x*v2.x;
838 
839  // Compute barycentric coordinates
840  double invDenom = 1/(dot00*dot11-dot01*dot01);
841  bcc.z = (dot11*dot02-dot01*dot12)*invDenom;
842  if( bcc.z < 0 ) return false;
843  bcc.y = (dot00*dot12-dot01*dot02)*invDenom;
844  bcc.x = 1. - bcc.z - bcc.y;
845 
846  // Check if point is in triangle
847  return (bcc.x>=0) && (bcc.y>=0);
848 }
849 
850 
851 bool
853 {
854  Point null;
855  return PointInTriangleBarycentric(pTest, p0, p1, p2, null);
856 }
857 
869 bool
871 {
872  // Compute vectors
873  Point v0 = p2 - p0;
874  Point v1 = p1 - p0;
875  Point v2 = pTest - p0;
876 
877  // Compute dot products
878  double dot00 = dot(v0,v0);
879  double dot01 = dot(v0,v1);
880  double dot02 = dot(v0,v2);
881  double dot11 = dot(v1,v1);
882  double dot12 = dot(v1,v2);
883 
884  // Compute barycentric coordinates
885  double invDenom = 1/(dot00*dot11-dot01*dot01);
886  bcc.z = (dot11*dot02-dot01*dot12)*invDenom;
887  bcc.y = (dot00*dot12-dot01*dot02)*invDenom;
888  bcc.x = 1. - bcc.z - bcc.y;
889 
890  // Check if point is in triangle
891  return (bcc.x>=0) && (bcc.y>=0) && (bcc.z>=0);
892 }
893 
894 
906 bool
907 fl_findEdgeSection(float d_p0, float d_p1, float th, Point *p0, Point *p1, Point *pSec)
908 {
909  bool hasIso = false;
910 
911  if(d_p0>=th) {
912  if (d_p1<th)
913  hasIso = true;
914  }
915  else {
916  if (d_p1>=th)
917  hasIso = true;
918  }
919 
920  if(hasIso) {
921  Point l = *p1 - *p0;
922  float len = mag(l);
923  float sec = (d_p0-th)/(d_p0-d_p1);
924  *pSec = *p0 + scal_X(l,sec);
925  }
926 
927  return hasIso;
928 }
929 
930 
940 bool
941 fl_findIsoIntersect(float *d0, float *d1, float th, Point *p, Elem *e )
942 {
943  bool isFilament = false;
944 
945  FMatrix M(4,4);
946  FMatrix Id(3,3);
947 
948  for (int i=0;i<4;i++) {
949  M.Ent[i][0] = 1.;
950  M.Ent[i][1] = p[i].x;
951  M.Ent[i][2] = p[i].y;
952  M.Ent[i][3] = p[i].z;
953  }
954 
955  double v0[4], v1[4], vf[4];
956  for(int i=0;i<4;i++) v0[i] = d0[i];
957  for(int i=0;i<4;i++) v1[i] = d1[i];
958 // for(int i=0;i<4;i++) vf[i] = phiIsoFace[i];
959 
960  M.LUD();
961 
962  // find base coefficients V(x,y,z) = a + bx +cy + dz
963  M.LUsolve( v0 );
964  M.LUsolve( v1 );
965 
966  // now we try to intesect the two hyperplanes with isosurface planes
967  // spanned by the faces of the tet, i.e. we assume a face is an isosurface hyperplane
968  int intersect = 0, NFaces = 4;
969  for(int i=0;i<NFaces;i++) {
970  // assign threshold values to nodes which span isosurface
971  for(int j=0;j<4;j++) vf[j] = th+1;
972  for(int j=0;j<3;j++) vf[TetFaces[i][j]] = th;
973 
974  // get base coeffs for this hyperplane
975  M.LUsolve( vf );
976 
977  // fill in rhs which is vth-a
978  double r[3];
979  r[0] = th - v0[0];
980  r[1] = th - v1[0];
981  r[2] = th - vf[0];
982 
983  for(int j=0;j<3;j++) {
984  Id.Ent[0][j] = v0[j+1];
985  Id.Ent[1][j] = v1[j+1];
986  Id.Ent[2][j] = vf[j+1];
987  }
988 
989  Id.LUD( );
990  if(Id.decomposed) {
991  Id.LUsolve( r );
992 
993  // this is primitive and not exactly correct,
994  // a proper test of whether the intersection is within the tet or not
995  // should go here instead
996  float xmn,xmx,ymn,ymx,zmn,zmx;
997  xmn = xmx = p[0].x;
998  ymn = ymx = p[0].y;
999  zmn = zmx = p[0].z;
1000 
1001  for(int j=1;j<4;j++) {
1002  if(p[j].x<xmn) xmn = p[j].x;
1003  if(p[j].x>xmx) xmx = p[j].x;
1004  if(p[j].y<ymn) ymn = p[j].y;
1005  if(p[j].y>ymx) ymx = p[j].y;
1006  if(p[j].z<zmn) zmn = p[j].z;
1007  if(p[j].z>zmx) zmx = p[j].z;
1008  }
1009 
1010  if((r[0]>=xmn) && (r[0]<=xmx) &&
1011  (r[1]>=ymn) && (r[1]<=ymx) &&
1012  (r[2]>=zmn) && (r[2]<=zmx) ) {
1013  intersect++;
1014  log_msg(NULL,0,0, "Intersection found at [x,y,z]=[%.2f,%.2f,%.2f]",
1015  r[0],r[1],r[2]);
1016  }
1017  }
1018 
1019  }
1020  if( ((e->type==Tri) && (intersect==1)) ||
1021  ((e->type==Tetra) && (intersect==2)) )
1022  isFilament = true;
1023 
1024  return isFilament;
1025 }
1026 
double Real
Definition: DataTypes.h:14
void dot(data_t *x, data_t *y, int n, int nc, FILE *out, bool pt)
perform dot product
Definition: IGBops.cc:550
void normalize(IGBheader *h, FILE *out, T *min, T *max)
normalize pixels over all time
Definition: IGBops.cc:637
bool decomposed
true if decomposed
Definition: FMatrix.h:28
Real ** Ent
matrix entries
Definition: FMatrix.h:27
void LUD()
Definition: FMatrix.cc:320
void LUsolve(Real *)
Definition: FMatrix.cc:355
int numPts
Definition: filament.h:43
int num_segs()
Definition: filament.h:41
int num_pts()
Definition: filament.h:42
void sort()
Definition: filament.cc:253
void add_singularity(singularity *s, int)
Definition: filament.cc:244
std::vector< singularity > s
Definition: filament.h:44
#define FIL_ELEM_FILE
Definition: filament.cc:27
float fl_triArea(Point *pts)
Definition: filament.cc:92
int init_filament_IO(const char *meshname, int numT, FILE *fh[])
Definition: filament.cc:270
void interpolate_vm(int numNodes, float **d, int M, float **rd, float rdist)
Definition: filament.cc:147
Point UnTransformPoint(Point PTs, Point *opts, Point *spts)
Definition: filament.cc:69
bool fl_chkIso(float *d, float th, Elem *e)
Definition: filament.cc:429
int write_connections(ElemList *elst, NodeList *nlst, filament *f, FILE *filHdls[], float t)
Definition: filament.cc:309
int fl_FindFilament(float *d0, float *d1, float th, Point *lp, Elem *e, singularity *sng)
Definition: filament.cc:543
bool PtInTriangleBarycentricSP(Point pTest, Point p1, Point p2)
Definition: filament.cc:811
void interpolateTBuffer(tbuffer *tbf)
Definition: filament.cc:124
#define FIL_DATA_FILE
Definition: filament.cc:28
bool SameSideOfLine(Point pTest, Point p0, Point p1, Point p2)
Definition: filament.cc:780
void interpoate_vm(int numNodes, float **d, int M, float **rd, float rDist)
bool fl_IsoTagFaces(float *d, float th, Elem *e, bool *faceTags)
Definition: filament.cc:456
bool fl_triFindPS_ShpFnc(float *d0, float *d1, float th, Point *p, Point *PS)
Definition: filament.cc:629
void write_filaments(filament *f, FILE *filHdls[], float t)
Definition: filament.cc:348
bool fl_findIsoIntersect(float *d0, float *d1, float th, Point *p, Elem *e)
Definition: filament.cc:941
const int MAX_ELEM_NODES
Definition: filament.cc:30
bool fl_IsoTagTriangle(float *d, float th, const int idx[])
Definition: filament.cc:518
void close_filament_IO(FILE *fh[])
Definition: filament.cc:414
bool PointInTriangle(Point pTest, Point p0, Point p1, Point p2)
Definition: filament.cc:799
const int MAX_ELEM_FACES
Definition: filament.cc:31
bool fl_findEdgeSection(float d_p0, float d_p1, float th, Point *p0, Point *p1, Point *pSec)
Definition: filament.cc:907
void output_filaments(filament *f, FILE *fh[])
Definition: filament.cc:375
bool PointInTriangleBarycentric(Point pTest, Point p0, Point p1, Point p2)
Definition: filament.cc:852
bool compute_filament(float *d0, float *d1, float vth, ElemList *elst, NodeList *nlst, filament *f, SingularityDetector_t meth)
Definition: filament.cc:171
#define FIL_PTS_FILE
Definition: filament.cc:26
bool fl_triFindPS_Lines(float *d0, float *d1, float th, Point *p, Point *PS)
Definition: filament.cc:647
bool fl_IsoTagTetFaces(float *d, float th, bool *faceTags)
Definition: filament.cc:483
Point fl_triNormalLocal(Point *pts, bool u)
Definition: filament.cc:108
bool fl_IsoTagTriFaces(float *d, float th, bool *faceTags)
Definition: filament.cc:504
bool fl_triFindPS(float *d0, float *d1, float th, Point *p, Point *PS, int meth)
Definition: filament.cc:602
SingularityDetector_t
Definition: filament.h:14
@ vol_based
Definition: filament.h:14
#define log_msg(F, L, O,...)
Definition: filament.h:8
@ PhaseSingularity
Definition: filament.h:15
@ Tri
Definition: filament.h:16
@ Tetra
Definition: filament.h:16
double mag(const Point &vect)
vector magnitude
Definition: SF_container.h:111
Point cross(const Point &a, const Point &b)
cross product
Definition: SF_container.h:84
vec3< V > scal_X(const vec3< V > &a, S k)
Definition: vect.h:150
vec3< POINT_REAL > Point
Definition: vect.h:93
V mag2(const vec3< V > &vect)
Definition: vect.h:138
Integer NElems
number of elements (owned+ghost in case of LOCAL_ELEM)
Definition: filament.h:60
Elem * Ele
Definition: filament.h:61
Definition: filament.h:53
Elem_t type
Definition: filament.h:54
int Nnode
Definition: filament.h:55
int * N
Definition: filament.h:56
Point * Node
list of nodes*
Definition: filament.h:49
Point p1
Definition: filament.h:28
Point p0
Definition: filament.h:27
int eIndx
Definition: filament.h:26
int type
Definition: filament.h:25
Definition: filament.h:18
int type
Definition: filament.h:19
Point p0
Definition: filament.h:21
int eIndx
Definition: filament.h:20
time slice buffer for post processing
Definition: filament.h:66
float ** dr
resampled data buffer
Definition: filament.h:76
int wdth
width of buffer, always 2*h_wdth+1
Definition: filament.h:74
float ** d
data buffer
Definition: filament.h:69
float rdist
rel dist of sample to center
Definition: filament.h:78
int nodes
number of nodes in time slice
Definition: filament.h:72
FilSeg fs
Definition: filament.h:33