1 /*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkQuadraticPyramid.cxx
5
6 Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7 All rights reserved.
8 See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9
10 This software is distributed WITHOUT ANY WARRANTY; without even
11 the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12 PURPOSE. See the above copyright notice for more information.
13
14 =========================================================================*/
15 #include "vtkQuadraticPyramid.h"
16
17 #include "vtkObjectFactory.h"
18 #include "vtkCellData.h"
19 #include "vtkDoubleArray.h"
20 #include "vtkTetra.h"
21 #include "vtkPyramid.h"
22 #include "vtkMath.h"
23 #include "vtkPointData.h"
24 #include "vtkQuadraticEdge.h"
25 #include "vtkQuadraticQuad.h"
26 #include "vtkQuadraticTriangle.h"
27 #include "vtkPoints.h"
28
29 vtkStandardNewMacro(vtkQuadraticPyramid);
30
31 //----------------------------------------------------------------------------
32 // Construct the pyramid with 13 points + 1 extra point for internal
33 // computation.
34 //
vtkQuadraticPyramid()35 vtkQuadraticPyramid::vtkQuadraticPyramid()
36 {
37 // At creation time the cell looks like it has 14 points (during interpolation)
38 // We initially allocate for 14.
39 this->Points->SetNumberOfPoints(14);
40 this->PointIds->SetNumberOfIds(14);
41 for (int i = 0; i < 14; i++)
42 {
43 this->Points->SetPoint(i, 0.0, 0.0, 0.0);
44 this->PointIds->SetId(i,0);
45 }
46 this->Points->SetNumberOfPoints(13);
47 this->PointIds->SetNumberOfIds(13);
48
49 this->Edge = vtkQuadraticEdge::New();
50 this->Face = vtkQuadraticQuad::New();
51 this->TriangleFace = vtkQuadraticTriangle::New();
52 this->Tetra = vtkTetra::New();
53 this->Pyramid = vtkPyramid::New();
54
55 this->PointData = vtkPointData::New();
56 this->CellData = vtkCellData::New();
57 this->CellScalars = vtkDoubleArray::New();
58 this->CellScalars->SetNumberOfTuples(14);
59 this->Scalars = vtkDoubleArray::New();
60 this->Scalars->SetNumberOfTuples(5); //num of vertices
61 }
62
63 //----------------------------------------------------------------------------
~vtkQuadraticPyramid()64 vtkQuadraticPyramid::~vtkQuadraticPyramid()
65 {
66 this->Edge->Delete();
67 this->Face->Delete();
68 this->TriangleFace->Delete();
69 this->Tetra->Delete();
70 this->Pyramid->Delete();
71
72 this->PointData->Delete();
73 this->CellData->Delete();
74 this->Scalars->Delete();
75 this->CellScalars->Delete();
76 }
77
78 //----------------------------------------------------------------------------
79 static int LinearPyramids[10][5] = { {0,5,13,8,9},
80 {5,1,6,13,10},
81 {8,13,7,3,12},
82 {13,6,2,7,11},
83 {9,10,11,12,4},
84 {9,12,11,10,13},
85 {5,10,9,13,0},
86 {6,11,10,13,0},
87 {7,12,11,13,0},
88 {8,9,12,13,0} };
89
90 static int PyramidFaces[5][8] = { {0,3,2,1,8,7,6,5},
91 {0,1,4,5,10,9,0,0},
92 {1,2,4,6,11,10,0,0},
93 {2,3,4,7,12,11,0,0},
94 {3,0,4,8,9,12,0,0}};
95
96 static int PyramidEdges[8][3] = { {0,1,5}, {1,2,6}, {2,3,7},
97 {3,0,8},{0,4,9},{1,4,10},
98 {2,4,11}, {3,4,12} };
99
100 static double MidPoints[1][3] = { {0.5,0.5,0.0} };
101 //----------------------------------------------------------------------------
GetEdgeArray(int edgeId)102 int *vtkQuadraticPyramid::GetEdgeArray(int edgeId)
103 {
104 return PyramidEdges[edgeId];
105 }
106 //----------------------------------------------------------------------------
GetFaceArray(int faceId)107 int *vtkQuadraticPyramid::GetFaceArray(int faceId)
108 {
109 return PyramidFaces[faceId];
110 }
111
112 //----------------------------------------------------------------------------
GetEdge(int edgeId)113 vtkCell *vtkQuadraticPyramid::GetEdge(int edgeId)
114 {
115 edgeId = (edgeId < 0 ? 0 : (edgeId > 7 ? 7 : edgeId ));
116
117 for (int i=0; i<3; i++)
118 {
119 this->Edge->PointIds->SetId(i,this->PointIds->GetId(PyramidEdges[edgeId][i]));
120 this->Edge->Points->SetPoint(i,this->Points->GetPoint(PyramidEdges[edgeId][i]));
121 }
122
123 return this->Edge;
124 }
125
126 //----------------------------------------------------------------------------
GetFace(int faceId)127 vtkCell *vtkQuadraticPyramid::GetFace(int faceId)
128 {
129 faceId = (faceId < 0 ? 0 : (faceId > 4 ? 4 : faceId ));
130
131 // load point id's and coordinates
132 // be careful with the first one:
133 if(faceId > 0)
134 {
135 for (int i=0; i<6; i++)
136 {
137 this->TriangleFace->PointIds->SetId(i,this->PointIds->GetId(PyramidFaces[faceId][i]));
138 this->TriangleFace->Points->SetPoint(i,this->Points->GetPoint(PyramidFaces[faceId][i]));
139 }
140 return this->TriangleFace;
141 }
142 else
143 {
144 for (int i=0; i<8; i++)
145 {
146 this->Face->PointIds->SetId(i,this->PointIds->GetId(PyramidFaces[faceId][i]));
147 this->Face->Points->SetPoint(i,this->Points->GetPoint(PyramidFaces[faceId][i]));
148 }
149 return this->Face;
150 }
151 }
152
153 //----------------------------------------------------------------------------
154 static const double VTK_DIVERGED = 1.e6;
155 static const int VTK_PYRAMID_MAX_ITERATION=10;
156 static const double VTK_PYRAMID_CONVERGED=1.e-03;
157
EvaluatePosition(double * x,double * closestPoint,int & subId,double pcoords[3],double & dist2,double * weights)158 int vtkQuadraticPyramid::EvaluatePosition(double* x,
159 double* closestPoint,
160 int& subId, double pcoords[3],
161 double& dist2, double *weights)
162 {
163 int iteration, converged;
164 double params[3];
165 double fcol[3], rcol[3], scol[3], tcol[3];
166 int i, j;
167 double d, pt[3];
168 double derivs[3*13];
169
170 // set initial position for Newton's method
171 subId = 0;
172 pcoords[0] = pcoords[1] = pcoords[2] = params[0] = params[1] = params[2]=0.5;
173
174 // enter iteration loop
175 for (iteration=converged=0;
176 !converged && (iteration < VTK_PYRAMID_MAX_ITERATION); iteration++)
177 {
178 // calculate element interpolation functions and derivatives
179 this->InterpolationFunctions(pcoords, weights);
180 this->InterpolationDerivs(pcoords, derivs);
181
182 // calculate newton functions
183 for (i=0; i<3; i++)
184 {
185 fcol[i] = rcol[i] = scol[i] = tcol[i] = 0.0;
186 }
187 for (i=0; i<13; i++)
188 {
189 this->Points->GetPoint(i, pt);
190 for (j=0; j<3; j++)
191 {
192 fcol[j] += pt[j] * weights[i];
193 rcol[j] += pt[j] * derivs[i];
194 scol[j] += pt[j] * derivs[i+13];
195 tcol[j] += pt[j] * derivs[i+26];
196 }
197 }
198
199 for (i=0; i<3; i++)
200 {
201 fcol[i] -= x[i];
202 }
203
204 // compute determinants and generate improvements
205 d=vtkMath::Determinant3x3(rcol,scol,tcol);
206 if ( fabs(d) < 1.e-20)
207 {
208 return -1;
209 }
210
211 pcoords[0] = params[0] - 0.5*vtkMath::Determinant3x3 (fcol,scol,tcol) / d;
212 pcoords[1] = params[1] - 0.5*vtkMath::Determinant3x3 (rcol,fcol,tcol) / d;
213 pcoords[2] = params[2] - 0.5*vtkMath::Determinant3x3 (rcol,scol,fcol) / d;
214
215 // check for convergence
216 if ( ((fabs(pcoords[0]-params[0])) < VTK_PYRAMID_CONVERGED) &&
217 ((fabs(pcoords[1]-params[1])) < VTK_PYRAMID_CONVERGED) &&
218 ((fabs(pcoords[2]-params[2])) < VTK_PYRAMID_CONVERGED) )
219 {
220 converged = 1;
221 }
222
223 // Test for bad divergence (S.Hirschberg 11.12.2001)
224 else if ((fabs(pcoords[0]) > VTK_DIVERGED) ||
225 (fabs(pcoords[1]) > VTK_DIVERGED) ||
226 (fabs(pcoords[2]) > VTK_DIVERGED))
227 {
228 return -1;
229 }
230
231 // if not converged, repeat
232 else
233 {
234 params[0] = pcoords[0];
235 params[1] = pcoords[1];
236 params[2] = pcoords[2];
237 }
238 }
239
240 // if not converged, set the parametric coordinates to arbitrary values
241 // outside of element
242 if ( !converged )
243 {
244 return -1;
245 }
246
247 this->InterpolationFunctions(pcoords, weights);
248
249 if ( pcoords[0] >= -0.001 && pcoords[0] <= 1.001 &&
250 pcoords[1] >= -0.001 && pcoords[1] <= 1.001 &&
251 pcoords[2] >= -0.001 && pcoords[2] <= 1.001 )
252 {
253 if (closestPoint)
254 {
255 closestPoint[0] = x[0]; closestPoint[1] = x[1]; closestPoint[2] = x[2];
256 dist2 = 0.0; //inside pyramid
257 }
258 return 1;
259 }
260 else
261 {
262 double pc[3], w[13];
263 if (closestPoint)
264 {
265 for (i=0; i<3; i++) //only approximate, not really true for warped hexa
266 {
267 if (pcoords[i] < 0.0)
268 {
269 pc[i] = 0.0;
270 }
271 else if (pcoords[i] > 1.0)
272 {
273 pc[i] = 1.0;
274 }
275 else
276 {
277 pc[i] = pcoords[i];
278 }
279 }
280 this->EvaluateLocation(subId, pc, closestPoint,
281 static_cast<double *>(w));
282 dist2 = vtkMath::Distance2BetweenPoints(closestPoint,x);
283 }
284 return 0;
285 }
286 }
287
288 //----------------------------------------------------------------------------
EvaluateLocation(int & vtkNotUsed (subId),double pcoords[3],double x[3],double * weights)289 void vtkQuadraticPyramid::EvaluateLocation(int& vtkNotUsed(subId),
290 double pcoords[3],
291 double x[3], double *weights)
292 {
293 int i, j;
294 double pt[3];
295
296 this->InterpolationFunctions(pcoords, weights);
297
298 x[0] = x[1] = x[2] = 0.0;
299 for (i=0; i<13; i++)
300 {
301 this->Points->GetPoint(i, pt);
302 for (j=0; j<3; j++)
303 {
304 x[j] += pt[j] * weights[i];
305 }
306 }
307 }
308
309 //----------------------------------------------------------------------------
CellBoundary(int subId,double pcoords[3],vtkIdList * pts)310 int vtkQuadraticPyramid::CellBoundary(int subId, double pcoords[3],
311 vtkIdList *pts)
312 {
313 return this->Pyramid->CellBoundary(subId, pcoords, pts);
314 }
315
316 //----------------------------------------------------------------------------
Subdivide(vtkPointData * inPd,vtkCellData * inCd,vtkIdType cellId,vtkDataArray * cellScalars)317 void vtkQuadraticPyramid::Subdivide(vtkPointData *inPd, vtkCellData *inCd,
318 vtkIdType cellId, vtkDataArray *cellScalars)
319 {
320 int numMidPts, i, j;
321 double weights[13];
322 double x[3];
323 double s;
324
325 //Copy point and cell attribute data, first make sure it's empty:
326 this->PointData->Initialize();
327 this->CellData->Initialize();
328 // Make sure to copy ALL arrays. These field data have to be
329 // identical to the input field data. Otherwise, CopyData
330 // that occurs later may not work because the output field
331 // data was initialized (CopyAllocate) with the input field
332 // data.
333 this->PointData->CopyAllOn();
334 this->CellData->CopyAllOn();
335 this->PointData->CopyAllocate(inPd,14);
336 this->CellData->CopyAllocate(inCd,6);
337 for (i=0; i<13; i++)
338 {
339 this->PointData->CopyData(inPd,this->PointIds->GetId(i),i);
340 this->CellScalars->SetValue( i, cellScalars->GetTuple1(i));
341 }
342 for (i=0; i<6; i++)
343 {
344 this->CellData->CopyData(inCd,cellId,i);
345 }
346
347 //Interpolate new values
348 double p[3];
349 this->Points->Resize(14);
350 this->CellScalars->Resize(14);
351 for ( numMidPts=0; numMidPts < 1; numMidPts++ )
352 {
353 this->InterpolationFunctions(MidPoints[numMidPts], weights);
354
355 x[0] = x[1] = x[2] = 0.0;
356 s = 0.0;
357 for (i=0; i<13; i++)
358 {
359 this->Points->GetPoint(i, p);
360 for (j=0; j<3; j++)
361 {
362 x[j] += p[j] * weights[i];
363 }
364 s += cellScalars->GetTuple1(i) * weights[i];
365 }
366 this->Points->SetPoint(13+numMidPts,x);
367 this->CellScalars->SetValue(13+numMidPts,s);
368 this->PointData->InterpolatePoint(inPd, 13+numMidPts,
369 this->PointIds, weights);
370 }
371 }
372
373 //----------------------------------------------------------------------------
Contour(double value,vtkDataArray * cellScalars,vtkIncrementalPointLocator * locator,vtkCellArray * verts,vtkCellArray * lines,vtkCellArray * polys,vtkPointData * inPd,vtkPointData * outPd,vtkCellData * inCd,vtkIdType cellId,vtkCellData * outCd)374 void vtkQuadraticPyramid::Contour(double value,
375 vtkDataArray* cellScalars,
376 vtkIncrementalPointLocator* locator,
377 vtkCellArray *verts,
378 vtkCellArray* lines,
379 vtkCellArray* polys,
380 vtkPointData* inPd,
381 vtkPointData* outPd,
382 vtkCellData* inCd,
383 vtkIdType cellId,
384 vtkCellData* outCd)
385 {
386 int i;
387 //subdivide into 6 linear pyramids
388 this->Subdivide(inPd,inCd,cellId,cellScalars);
389
390 //contour each linear pyramid separately
391 this->Scalars->SetNumberOfTuples(5); //num of vertices
392 for (i=0; i<6; i++) //for each pyramid
393 {
394 for (int j=0; j<5; j++) //for each point of pyramid
395 {
396 this->Pyramid->Points->SetPoint(j,this->Points->GetPoint(LinearPyramids[i][j]));
397 this->Pyramid->PointIds->SetId(j,LinearPyramids[i][j]);
398 this->Scalars->SetValue(j,this->CellScalars->GetValue(LinearPyramids[i][j]));
399 }
400 this->Pyramid->Contour(value,this->Scalars,locator,verts,lines,polys,
401 this->PointData,outPd,this->CellData,cellId,outCd);
402 }
403
404 //contour each linear tetra separately
405 this->Scalars->SetNumberOfTuples(4); //num of vertices
406 for (i=6; i<10; i++) //for each tetra
407 {
408 for (int j=0; j<4; j++) //for each point of tetra
409 {
410 this->Tetra->Points->SetPoint(j,this->Points->GetPoint(LinearPyramids[i][j]));
411 this->Tetra->PointIds->SetId(j,LinearPyramids[i][j]);
412 this->Scalars->SetTuple(j,this->CellScalars->GetTuple(LinearPyramids[i][j]));
413 }
414 this->Tetra->Contour(value,this->Scalars,locator,verts,lines,polys,
415 this->PointData,outPd,this->CellData,cellId,outCd);
416 }
417 }
418
419 //----------------------------------------------------------------------------
420 // Line-hex intersection. Intersection has to occur within [0,1] parametric
421 // coordinates and with specified tolerance.
422 //
IntersectWithLine(double * p1,double * p2,double tol,double & t,double * x,double * pcoords,int & subId)423 int vtkQuadraticPyramid::IntersectWithLine(double* p1, double* p2,
424 double tol, double& t,
425 double* x, double* pcoords, int& subId)
426 {
427 int intersection=0;
428 double tTemp;
429 double pc[3], xTemp[3];
430 int faceNum;
431 int inter;
432
433 t = VTK_DOUBLE_MAX;
434 for (faceNum=0; faceNum<5; faceNum++)
435 {
436 // We have 8 nodes on rect face
437 // and 6 on triangle faces
438 if(faceNum > 0)
439 {
440 for (int i=0; i<6; i++)
441 {
442 this->TriangleFace->PointIds->SetId(i,
443 this->PointIds->GetId(PyramidFaces[faceNum][i]));
444 }
445 inter = this->TriangleFace->IntersectWithLine(p1, p2, tol, tTemp,
446 xTemp, pc, subId);
447 }
448 else
449 {
450 for (int i=0; i<8; i++)
451 {
452 this->Face->Points->SetPoint(i,
453 this->Points->GetPoint(PyramidFaces[faceNum][i]));
454 }
455 inter = this->Face->IntersectWithLine(p1, p2, tol, tTemp,
456 xTemp, pc, subId);
457 }
458 if ( inter )
459 {
460 intersection = 1;
461 if ( tTemp < t )
462 {
463 t = tTemp;
464 x[0] = xTemp[0]; x[1] = xTemp[1]; x[2] = xTemp[2];
465 switch (faceNum)
466 {
467 case 0:
468 pcoords[0] = 0.0; pcoords[1] = pc[1]; pcoords[2] = pc[0];
469 break;
470
471 case 1:
472 pcoords[0] = 1.0; pcoords[1] = pc[0]; pcoords[2] = pc[1];
473 break;
474
475 case 2:
476 pcoords[0] = pc[0]; pcoords[1] = 0.0; pcoords[2] = pc[1];
477 break;
478
479 case 3:
480 pcoords[0] = pc[1]; pcoords[1] = 1.0; pcoords[2] = pc[0];
481 break;
482
483 case 4:
484 pcoords[0] = pc[1]; pcoords[1] = pc[0]; pcoords[2] = 0.0;
485 break;
486
487 case 5:
488 pcoords[0] = pc[0]; pcoords[1] = pc[1]; pcoords[2] = 1.0;
489 break;
490 }
491 }
492 }
493 }
494 return intersection;
495 }
496
497 //----------------------------------------------------------------------------
Triangulate(int vtkNotUsed (index),vtkIdList * ptIds,vtkPoints * pts)498 int vtkQuadraticPyramid::Triangulate(int vtkNotUsed(index), vtkIdList *ptIds,
499 vtkPoints *pts)
500 {
501 int i;
502 int ii;
503 pts->Reset();
504 ptIds->Reset();
505
506 for (i=0; i < 6; i++)
507 {
508 for ( int j=0; j < 5; j++)
509 {
510 ptIds->InsertId(5*i+j,this->PointIds->GetId(LinearPyramids[i][j]));
511 pts->InsertPoint(5*i+j,this->Points->GetPoint(LinearPyramids[i][j]));
512 }
513 }
514
515 for (ii=0, i=6 ; i < 10; i++, ii++)
516 {
517 for ( int j=0; j < 4; j++)
518 {
519 ptIds->InsertId(4*ii+j+30,this->PointIds->GetId(LinearPyramids[i][j]));
520 pts->InsertPoint(4*ii+j+30,this->Points->GetPoint(LinearPyramids[i][j]));
521 }
522 }
523
524 return 1;
525 }
526
527 //----------------------------------------------------------------------------
528 // Given parametric coordinates compute inverse Jacobian transformation
529 // matrix. Returns 9 elements of 3x3 inverse Jacobian plus interpolation
530 // function derivatives.
531 //
JacobianInverse(double pcoords[3],double ** inverse,double derivs[39])532 void vtkQuadraticPyramid::JacobianInverse(double pcoords[3], double **inverse,
533 double derivs[39])
534 {
535 int i, j;
536 double *m[3], m0[3], m1[3], m2[3];
537 double x[3];
538
539 // compute interpolation function derivatives
540 this->InterpolationDerivs(pcoords, derivs);
541
542 // create Jacobian matrix
543 m[0] = m0; m[1] = m1; m[2] = m2;
544 for (i=0; i < 3; i++) //initialize matrix
545 {
546 m0[i] = m1[i] = m2[i] = 0.0;
547 }
548
549 for ( j=0; j < 13; j++ )
550 {
551 this->Points->GetPoint(j, x);
552 for ( i=0; i < 3; i++ )
553 {
554 m0[i] += x[i] * derivs[j];
555 m1[i] += x[i] * derivs[13 + j];
556 m2[i] += x[i] * derivs[26 + j];
557 }
558 }
559
560 // now find the inverse
561 if ( vtkMath::InvertMatrix(m,inverse,3) == 0 )
562 {
563 vtkErrorMacro(<<"Jacobian inverse not found");
564 return;
565 }
566 }
567
568 //----------------------------------------------------------------------------
Derivatives(int vtkNotUsed (subId),double pcoords[3],double * values,int dim,double * derivs)569 void vtkQuadraticPyramid::Derivatives(int vtkNotUsed(subId),
570 double pcoords[3], double *values,
571 int dim, double *derivs)
572 {
573 double *jI[3], j0[3], j1[3], j2[3];
574 double functionDerivs[3*13], sum[3];
575 int i, j, k;
576
577 // compute inverse Jacobian and interpolation function derivatives
578 jI[0] = j0; jI[1] = j1; jI[2] = j2;
579 this->JacobianInverse(pcoords, jI, functionDerivs);
580
581 // now compute derivates of values provided
582 for (k=0; k < dim; k++) //loop over values per vertex
583 {
584 sum[0] = sum[1] = sum[2] = 0.0;
585 for ( i=0; i < 13; i++) //loop over interp. function derivatives
586 {
587 sum[0] += functionDerivs[i] * values[dim*i + k];
588 sum[1] += functionDerivs[13 + i] * values[dim*i + k];
589 sum[2] += functionDerivs[26 + i] * values[dim*i + k];
590 }
591 for (j=0; j < 3; j++) //loop over derivative directions
592 {
593 derivs[3*k + j] = sum[0]*jI[j][0] + sum[1]*jI[j][1] + sum[2]*jI[j][2];
594 }
595 }
596 }
597
598 //----------------------------------------------------------------------------
599 // Clip this quadratic pyramid using scalar value provided. Like contouring,
600 // except that it cuts the pyramid to produce tetrahedra.
601 //
Clip(double value,vtkDataArray * cellScalars,vtkIncrementalPointLocator * locator,vtkCellArray * tets,vtkPointData * inPd,vtkPointData * outPd,vtkCellData * inCd,vtkIdType cellId,vtkCellData * outCd,int insideOut)602 void vtkQuadraticPyramid::Clip(double value, vtkDataArray* cellScalars,
603 vtkIncrementalPointLocator* locator, vtkCellArray* tets,
604 vtkPointData* inPd, vtkPointData* outPd,
605 vtkCellData* inCd, vtkIdType cellId,
606 vtkCellData* outCd, int insideOut)
607 {
608 int i;
609 // create six linear pyramid + 4 tetra
610 this->Subdivide(inPd,inCd,cellId,cellScalars);
611
612 //contour each linear pyramid separately
613 this->Scalars->SetNumberOfTuples(5); //num of vertices
614 for (i=0; i<6; i++) //for each subdivided pyramid
615 {
616 for (int j=0; j<5; j++) //for each of the five vertices of the pyramid
617 {
618 this->Pyramid->Points->SetPoint(j,this->Points->GetPoint(LinearPyramids[i][j]));
619 this->Pyramid->PointIds->SetId(j,LinearPyramids[i][j]);
620 this->Scalars->SetValue(j,this->CellScalars->GetValue(LinearPyramids[i][j]));
621 }
622 this->Pyramid->Clip(value,this->Scalars,locator,tets,this->PointData,outPd,
623 this->CellData,cellId,outCd,insideOut);
624 }
625
626 this->Scalars->SetNumberOfTuples(4); //num of vertices
627 for (i=6; i<10; i++) //for each subdivided tetra
628 {
629 for (int j=0; j<4; j++) //for each of the four vertices of the tetra
630 {
631 this->Tetra->Points->SetPoint(j,this->Points->GetPoint(LinearPyramids[i][j]));
632 this->Tetra->PointIds->SetId(j,LinearPyramids[i][j]);
633 this->Scalars->SetValue(j,this->CellScalars->GetValue(LinearPyramids[i][j]));
634 }
635 this->Tetra->Clip(value,this->Scalars,locator,tets,this->PointData,outPd,
636 this->CellData,cellId,outCd,insideOut);
637 }
638 }
639
640 //----------------------------------------------------------------------------
641 // Compute interpolation functions for the fifteen nodes.
642 //
InterpolationFunctions(double pcoords[3],double weights[13])643 void vtkQuadraticPyramid::InterpolationFunctions(double pcoords[3],
644 double weights[13])
645 {
646 // VTK needs parametric coordinates to be between (0,1). Isoparametric
647 // shape functions are formulated between (-1,1). Here we do a
648 // coordinate system conversion from (0,1) to (-1,1).
649 const double r = 2*pcoords[0] - 1;
650 const double s = 2*pcoords[1] - 1;
651 const double t = 2*pcoords[2] - 1;
652
653 const double rm = 1.0 - r;
654 const double rp = 1.0 + r;
655 const double sm = 1.0 - s;
656 const double sp = 1.0 + s;
657 const double tm = 1.0 - t;
658 const double tp = 1.0 + t;
659 const double r2 = 1.0 - r*r;
660 const double s2 = 1.0 - s*s;
661 const double t2 = 1.0 - t*t;
662
663 // corners
664 weights[0] = 0.125 * rm * sm * tm * (-r - s - t - 2.0);
665 weights[1] = 0.125 * rp * sm * tm * ( r - s - t - 2.0);
666 weights[2] = 0.125 * rp * sp * tm * ( r + s - t - 2.0);
667 weights[3] = 0.125 * rm * sp * tm * (-r + s - t - 2.0);
668 weights[4] = 0.5 * t * tp;
669
670 // midsides of rectangles
671 weights[5] = 0.25 * r2 * sm * tm;
672 weights[6] = 0.25 * s2 * rp * tm;
673 weights[7] = 0.25 * r2 * sp * tm;
674 weights[8] = 0.25 * s2 * rm * tm;
675
676
677 // midsides of triangles
678 weights[9] = 0.25 * ( 1 - r ) * (1 - s ) * t2;
679 weights[10] = 0.25 * ( 1 + r ) * (1 - s ) * t2;
680 weights[11] = 0.25 * ( 1 + r ) * (1 + s ) * t2;
681 weights[12] = 0.25 * ( 1 - r ) * (1 + s ) * t2;
682 }
683
684 //----------------------------------------------------------------------------
685 // Derivatives in parametric space.
686 //
InterpolationDerivs(double pcoords[3],double derivs[39])687 void vtkQuadraticPyramid::InterpolationDerivs(double pcoords[3],
688 double derivs[39])
689 {
690 //VTK needs parametric coordinates to be between (0,1). Isoparametric
691 //shape functions are formulated between (-1,1). Here we do a
692 //coordinate system conversion from (0,1) to (-1,1).
693 double r = 2*pcoords[0] - 1;
694 double s = 2*pcoords[1] - 1;
695 double t = 2*pcoords[2] - 1;
696
697 double rm = 1.0 - r;
698 double rp = 1.0 + r;
699 double sm = 1.0 - s;
700 double sp = 1.0 + s;
701 double tm = 1.0 - t;
702 //double tp = 1.0 + t;
703 double r2 = 1.0 - r*r;
704 double t2 = 1.0 - t*t;
705
706 //r-derivatives
707 // corners
708 derivs[0] = -0.125*(sm*tm - 2.0*r*sm*tm - s*sm*tm - t*sm*tm - 2.0*sm*tm);
709 derivs[1] = 0.125*(sm*tm + 2.0*r*sm*tm - s*sm*tm - t*sm*tm - 2.0*sm*tm);
710 derivs[2] = 0.125*(sp*tm + 2.0*r*sp*tm + s*sp*tm - t*sp*tm - 2.0*sp*tm);
711 derivs[3] = -0.125*(sp*tm - 2.0*r*sp*tm + s*sp*tm - t*sp*tm - 2.0*sp*tm);
712 derivs[4] = 0.0;
713
714 // midsides of rectangles
715 derivs[5] = -0.5*r*sm*tm;
716 derivs[6] = 0.25*(tm - s*s*tm);
717 derivs[7] = -0.5*r*sp*tm;
718 derivs[8] = -0.25*(tm - s*s*tm);
719
720 // midsides of triangles
721 derivs[9] = -0.25 * (1 - s ) * (1 - t*t );
722 derivs[10] = 0.25 * (1 - s ) * (1 - t*t );
723 derivs[11] = 0.25 * (1 + s ) * (1 - t*t );
724 derivs[12] = -0.25 * (1 + s ) * (1 - t*t );
725
726 //s-derivatives
727 // corners
728 derivs[13] = -0.125*(rm*tm - 2.0*s*rm*tm - r*rm*tm - t*rm*tm - 2.0*rm*tm);
729 derivs[14] = -0.125*(rp*tm - 2.0*s*rp*tm + r*rp*tm - t*rp*tm - 2.0*rp*tm);
730 derivs[15] = 0.125*(rp*tm + 2.0*s*rp*tm + r*rp*tm - t*rp*tm - 2.0*rp*tm);
731 derivs[16] = 0.125*(rm*tm + 2.0*s*rm*tm - r*rm*tm - t*rm*tm - 2.0*rm*tm);
732 derivs[17] = 0.0;
733
734 // midsides of rectangles
735 derivs[18] = -0.25 * tm * r2;
736 derivs[19] = -0.5 * tm * s * rp;
737 derivs[20] = 0.25 * tm * r2;
738 derivs[21] = -0.5 * tm * s * rm;
739
740 // midsides of triangles
741 derivs[22] = -0.25 * rm * t2;
742 derivs[23] = -0.25 * rp * t2;
743 derivs[24] = 0.25 * rp * t2;
744 derivs[25] = 0.25 * rm * t2;
745
746 //t-derivatives
747 // corners
748 derivs[26] = -0.125*(rm*sm - 2.0*t*rm*sm - r*rm*sm - s*rm*sm - 2.0*rm*sm);
749 derivs[27] = -0.125*(rp*sm - 2.0*t*rp*sm + r*rp*sm - s*rp*sm - 2.0*rp*sm);
750 derivs[28] = -0.125*(rp*sp - 2.0*t*rp*sp + r*rp*sp + s*rp*sp - 2.0*rp*sp);
751 derivs[29] = -0.125*(rm*sp - 2.0*t*rm*sp - r*rm*sp + s*rm*sp - 2.0*rm*sp);
752 derivs[30] = 0.5 + t;
753
754 // midsides of rectangles
755 derivs[31] = -0.25*(sm - r*r*sm);
756 derivs[32] = -0.25*(rp - s*s*rp);
757 derivs[33] = -0.25*(sp - r*r*sp);
758 derivs[34] = -0.25*(rm - s*s*rm);
759
760
761 // midsides of triangles
762 derivs[35] = -0.5 * ( 1 - r ) * (1 - s ) * t;
763 derivs[36] = -0.5 * ( 1 + r ) * (1 - s ) * t;
764 derivs[37] = -0.5 * ( 1 + r ) * (1 + s ) * t;
765 derivs[38] = -0.5 * ( 1 - r ) * (1 + s ) * t;
766
767 // we compute derivatives in in [-1; 1] but we need them in [ 0; 1]
768 for(int i = 0; i < 39; i++)
769 derivs[i] *= 2;
770 }
771
772 static double vtkQPyramidCellPCoords[39] = {0.0,0.0,0.0, 1.0,0.0,0.0, 1.0,1.0,0.0,
773 0.0,1.0,0.0, 0.0,0.0,1.0, 0.5,0.0,0.0,
774 1.0,0.5,0.0, 0.5,1.0,0.0, 0.0,0.5,0.0,
775 0.0,0.0,0.5, 1.0,0.0,0.5,
776 1.0,1.0,0.5, 0.0,1.0,0.5 };
777
778 //----------------------------------------------------------------------------
GetParametricCoords()779 double *vtkQuadraticPyramid::GetParametricCoords()
780 {
781 return vtkQPyramidCellPCoords;
782 }
783
784 //----------------------------------------------------------------------------
PrintSelf(ostream & os,vtkIndent indent)785 void vtkQuadraticPyramid::PrintSelf(ostream& os, vtkIndent indent)
786 {
787 this->Superclass::PrintSelf(os,indent);
788
789 os << indent << "Edge:\n";
790 this->Edge->PrintSelf(os,indent.GetNextIndent());
791 os << indent << "TriangleFace:\n";
792 this->TriangleFace->PrintSelf(os,indent.GetNextIndent());
793 os << indent << "Face:\n";
794 this->Face->PrintSelf(os,indent.GetNextIndent());
795 os << indent << "Tetra:\n";
796 this->Tetra->PrintSelf(os,indent.GetNextIndent());
797 os << indent << "Pyramid:\n";
798 this->Pyramid->PrintSelf(os,indent.GetNextIndent());
799 os << indent << "PointData:\n";
800 this->PointData->PrintSelf(os,indent.GetNextIndent());
801 os << indent << "CellData:\n";
802 this->CellData->PrintSelf(os,indent.GetNextIndent());
803 os << indent << "Scalars:\n";
804 this->Scalars->PrintSelf(os,indent.GetNextIndent());
805 }
806