1 /*=========================================================================
2 
3   Program:   Visualization Toolkit
4   Module:    vtkBiQuadraticQuadraticWedge.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 
16 //Thanks to Soeren Gebbert  who developed this class and
17 //integrated it into VTK 5.0.
18 
19 #include "vtkBiQuadraticQuadraticWedge.h"
20 
21 #include "vtkObjectFactory.h"
22 #include "vtkDoubleArray.h"
23 #include "vtkWedge.h"
24 #include "vtkMath.h"
25 #include "vtkQuadraticEdge.h"
26 #include "vtkBiQuadraticQuad.h"
27 #include "vtkQuadraticTriangle.h"
28 #include "vtkPoints.h"
29 
30 #include <cassert>
31 
32 vtkStandardNewMacro (vtkBiQuadraticQuadraticWedge);
33 
34 //----------------------------------------------------------------------------
35 // Construct the biquadratic quadratic wedge with 18 points
36 
vtkBiQuadraticQuadraticWedge()37 vtkBiQuadraticQuadraticWedge::vtkBiQuadraticQuadraticWedge ()
38 {
39   this->Points->SetNumberOfPoints (18);
40   this->PointIds->SetNumberOfIds (18);
41   for (int i = 0; i < 18; i++)
42     {
43     this->Points->SetPoint(i, 0.0, 0.0, 0.0);
44     this->PointIds->SetId(i,0);
45     }
46 
47   this->Edge = vtkQuadraticEdge::New ();
48   this->Face = vtkBiQuadraticQuad::New ();
49   this->TriangleFace = vtkQuadraticTriangle::New ();
50   this->Wedge = vtkWedge::New ();
51 
52   this->Scalars = vtkDoubleArray::New ();
53   this->Scalars->SetNumberOfTuples (6); //Number of vertices from a linear wedge
54 
55 }
56 
57 //----------------------------------------------------------------------------
~vtkBiQuadraticQuadraticWedge()58 vtkBiQuadraticQuadraticWedge::~vtkBiQuadraticQuadraticWedge ()
59 {
60   this->Edge->Delete ();
61   this->Face->Delete ();
62   this->TriangleFace->Delete ();
63   this->Wedge->Delete ();
64   this->Scalars->Delete ();
65 }
66 
67 //----------------------------------------------------------------------------
68 // We are using 8 linear wedge
69 static int LinearWedges[8][6] = { {0,6,8,12,15,17},
70                                   {6,7,8,15,16,17},
71                                   {6,1,7,15,13,16},
72                                   {8,7,2,17,16,14},
73                                   {12,15,17,3,9,11},
74                                   {15,16,17,9,10,11},
75                                   {15,13,16,9,4,10},
76                                   {17,16,14,11,10,5} };
77 
78 // We use 2 quadratic triangles and 3 quadratic-linear quads
79 static int WedgeFaces[5][9] = {
80     {0, 1, 2,  6,  7,  8,  0,  0,  0},  // first quad triangle
81     {3, 5, 4, 11, 10,  9,  0,  0,  0},  // second quad triangle
82     {0, 3, 4,  1, 12,  9, 13,  6, 15},  // 1. biquad quad
83     {1, 4, 5,  2, 13, 10, 14,  7, 16},  // 2. biquad quad
84     {2, 5, 3,  0, 14, 11, 12,  8, 17}   // 3. biquad quad
85 };
86 
87 // We have 9 quadratic edges
88 static int WedgeEdges[9][3] = {
89     {0, 1,  6}, {1, 2,  7}, {2, 0,  8},
90     {3, 4,  9}, {4, 5, 10}, {5, 3, 11},
91     {0, 3, 12}, {1, 4, 13}, {2, 5, 14}
92 };
93 //----------------------------------------------------------------------------
GetEdgeArray(int edgeId)94 int *vtkBiQuadraticQuadraticWedge::GetEdgeArray(int edgeId)
95 {
96   return WedgeEdges[edgeId];
97 }
98 //----------------------------------------------------------------------------
GetFaceArray(int faceId)99 int *vtkBiQuadraticQuadraticWedge::GetFaceArray(int faceId)
100 {
101   return WedgeFaces[faceId];
102 }
103 
104 //----------------------------------------------------------------------------
GetEdge(int edgeId)105 vtkCell *vtkBiQuadraticQuadraticWedge::GetEdge(int edgeId)
106 {
107   edgeId = (edgeId < 0 ? 0 : (edgeId > 8 ? 8 : edgeId));
108 
109   //We have 9 quadratic edges
110   for (int i = 0; i < 3; i++)
111     {
112     this->Edge->PointIds->SetId (i, this->PointIds->GetId (WedgeEdges[edgeId][i]));
113     this->Edge->Points->SetPoint (i, this->Points->GetPoint (WedgeEdges[edgeId][i]));
114     }
115 
116   return this->Edge;
117 }
118 
119 //----------------------------------------------------------------------------
GetFace(int faceId)120 vtkCell *vtkBiQuadraticQuadraticWedge::GetFace(int faceId)
121 {
122   faceId = (faceId < 0 ? 0 : (faceId > 4 ? 4 : faceId));
123 
124   // load point id's and coordinates
125   // be careful with the last two:
126   if (faceId < 2)
127     {
128     for (int i = 0; i < 6; i++)
129       {
130       this->TriangleFace->PointIds->SetId (i, this->PointIds->GetId (WedgeFaces[faceId][i]));
131       this->TriangleFace->Points->SetPoint (i, this->Points->GetPoint (WedgeFaces[faceId][i]));
132       }
133     return this->TriangleFace;
134     }
135   else
136     {
137     for (int i = 0; i < 9; i++)
138       {
139       this->Face->PointIds->SetId (i, this->PointIds->GetId (WedgeFaces[faceId][i]));
140       this->Face->Points->SetPoint (i, this->Points->GetPoint (WedgeFaces[faceId][i]));
141       }
142     return this->Face;
143     }
144 }
145 
146 //----------------------------------------------------------------------------
147 static const double VTK_DIVERGED = 1.e6;
148 static const int VTK_WEDGE_MAX_ITERATION = 20;
149 static const double VTK_WEDGE_CONVERGED = 1.e-03;
150 
EvaluatePosition(double * x,double * closestPoint,int & subId,double pcoords[3],double & dist2,double * weights)151 int vtkBiQuadraticQuadraticWedge::EvaluatePosition (double *x,
152   double *closestPoint,
153   int &subId, double pcoords[3], double &dist2, double *weights)
154 {
155   int iteration, converged;
156   double params[3];
157   double fcol[3], rcol[3], scol[3], tcol[3];
158   int i, j;
159   double d, pt[3];
160   double derivs[3 * 18];
161 
162   //  set initial position for Newton's method
163   subId = 0;
164   pcoords[0] = pcoords[1] = pcoords[2] = params[0] = params[1] = params[2] = 0.5;
165 
166   //  enter iteration loop
167   for (iteration = converged = 0; !converged && (iteration < VTK_WEDGE_MAX_ITERATION); iteration++)
168     {
169     //  calculate element interpolation functions and derivatives
170     this->InterpolationFunctions (pcoords, weights);
171     this->InterpolationDerivs (pcoords, derivs);
172 
173     //  calculate newton functions
174     for (i = 0; i < 3; i++)
175       {
176       fcol[i] = rcol[i] = scol[i] = tcol[i] = 0.0;
177       }
178     for (i = 0; i < 18; i++)
179       {
180       this->Points->GetPoint (i, pt);
181       for (j = 0; j < 3; j++)
182         {
183         fcol[j] += pt[j] * weights[i];
184         rcol[j] += pt[j] * derivs[i];
185         scol[j] += pt[j] * derivs[i + 18];
186         tcol[j] += pt[j] * derivs[i + 36];
187         }
188       }
189 
190     for (i = 0; i < 3; i++)
191       {
192       fcol[i] -= x[i];
193       }
194 
195     //  compute determinants and generate improvements
196     d = vtkMath::Determinant3x3 (rcol, scol, tcol);
197     if (fabs (d) < 1.e-20)
198       {
199       return -1;
200       }
201 
202     pcoords[0] = params[0] - 0.5 * vtkMath::Determinant3x3 (fcol, scol, tcol) / d;
203     pcoords[1] = params[1] - 0.5 * vtkMath::Determinant3x3 (rcol, fcol, tcol) / d;
204     pcoords[2] = params[2] - 0.5 * vtkMath::Determinant3x3 (rcol, scol, fcol) / d;
205 
206     //  check for convergence
207     if (((fabs (pcoords[0] - params[0])) < VTK_WEDGE_CONVERGED) &&
208       ((fabs (pcoords[1] - params[1])) < VTK_WEDGE_CONVERGED) &&
209       ((fabs (pcoords[2] - params[2])) < VTK_WEDGE_CONVERGED))
210       {
211       converged = 1;
212       }
213 
214     // Test for bad divergence (S.Hirschberg 11.12.2001)
215     else if ((fabs (pcoords[0]) > VTK_DIVERGED) ||
216       (fabs (pcoords[1]) > VTK_DIVERGED) || (fabs (pcoords[2]) > VTK_DIVERGED))
217       {
218       return -1;
219       }
220 
221     //  if not converged, repeat
222     else
223       {
224       params[0] = pcoords[0];
225       params[1] = pcoords[1];
226       params[2] = pcoords[2];
227       }
228     }
229 
230   //  if not converged, set the parametric coordinates to arbitrary values
231   //  outside of element
232   if (!converged)
233     {
234     return -1;
235     }
236 
237   this->InterpolationFunctions (pcoords, weights);
238 
239   if (pcoords[0] >= -0.001 && pcoords[0] <= 1.001 &&
240     pcoords[1] >= -0.001 && pcoords[1] <= 1.001 && pcoords[2] >= -0.001 && pcoords[2] <= 1.001)
241     {
242     if (closestPoint)
243       {
244       closestPoint[0] = x[0];
245       closestPoint[1] = x[1];
246       closestPoint[2] = x[2];
247       dist2 = 0.0;    //inside wedge
248       }
249     return 1;
250     }
251   else
252     {
253     double pc[3], w[18];
254     if (closestPoint)
255       {
256       for (i = 0; i < 3; i++)  //only approximate, not really true for warped hexa
257         {
258         if (pcoords[i] < 0.0)
259           {
260           pc[i] = 0.0;
261           }
262         else if (pcoords[i] > 1.0)
263           {
264           pc[i] = 1.0;
265           }
266         else
267           {
268           pc[i] = pcoords[i];
269           }
270         }
271       this->EvaluateLocation (subId, pc, closestPoint,
272                               static_cast<double *>(w));
273       dist2 = vtkMath::Distance2BetweenPoints (closestPoint, x);
274       }
275     return 0;
276     }
277 }
278 
279 //----------------------------------------------------------------------------
EvaluateLocation(int & vtkNotUsed (subId),double pcoords[3],double x[3],double * weights)280 void vtkBiQuadraticQuadraticWedge::EvaluateLocation (int &vtkNotUsed (subId),
281   double pcoords[3], double x[3], double *weights)
282 {
283   double pt[3];
284 
285   this->InterpolationFunctions (pcoords, weights);
286 
287   x[0] = x[1] = x[2] = 0.0;
288   for (int i = 0; i < 18; i++)
289     {
290     this->Points->GetPoint (i, pt);
291     for (int j = 0; j < 3; j++)
292       {
293       x[j] += pt[j] * weights[i];
294       }
295     }
296 }
297 
298 //----------------------------------------------------------------------------
CellBoundary(int subId,double pcoords[3],vtkIdList * pts)299 int vtkBiQuadraticQuadraticWedge::CellBoundary (int subId,
300   double pcoords[3], vtkIdList * pts)
301 {
302   return this->Wedge->CellBoundary (subId, pcoords, pts);
303 }
304 
305 //----------------------------------------------------------------------------
Contour(double value,vtkDataArray * cellScalars,vtkIncrementalPointLocator * locator,vtkCellArray * verts,vtkCellArray * lines,vtkCellArray * polys,vtkPointData * inPd,vtkPointData * outPd,vtkCellData * inCd,vtkIdType cellId,vtkCellData * outCd)306 void vtkBiQuadraticQuadraticWedge::Contour (double value,
307           vtkDataArray * cellScalars,
308           vtkIncrementalPointLocator * locator,
309           vtkCellArray * verts,
310           vtkCellArray * lines,
311           vtkCellArray * polys,
312           vtkPointData * inPd,
313           vtkPointData * outPd, vtkCellData * inCd,
314           vtkIdType cellId, vtkCellData * outCd)
315 {
316   //contour each linear wedge separately
317   for (int i=0; i<8; i++) //for each wedge
318     {
319     for (int j=0; j<6; j++) //for each point of wedge
320       {
321       this->Wedge->Points->SetPoint(j,this->Points->GetPoint(LinearWedges[i][j]));
322       this->Wedge->PointIds->SetId(j,this->PointIds->GetId(LinearWedges[i][j]));
323       this->Scalars->SetValue(j,cellScalars->GetTuple1(LinearWedges[i][j]));
324       }
325     this->Wedge->Contour(value,this->Scalars,locator,verts,lines,polys, inPd,outPd,inCd,cellId,outCd);
326     }
327 }
328 
329 //----------------------------------------------------------------------------
330 // Clip this biquadratic wedge using scalar value provided. Like contouring,
331 // except that it cuts the wedge to produce tetrahedra.
332 void
Clip(double value,vtkDataArray * cellScalars,vtkIncrementalPointLocator * locator,vtkCellArray * tets,vtkPointData * inPd,vtkPointData * outPd,vtkCellData * inCd,vtkIdType cellId,vtkCellData * outCd,int insideOut)333 vtkBiQuadraticQuadraticWedge::Clip (double value, vtkDataArray *cellScalars,
334              vtkIncrementalPointLocator * locator, vtkCellArray * tets,
335              vtkPointData * inPd, vtkPointData * outPd,
336              vtkCellData * inCd, vtkIdType cellId, vtkCellData * outCd, int insideOut)
337 {
338   //contour each linear wedge separately
339   for (int i=0; i<8; i++) //for each wedge
340     {
341     for (int j=0; j<6; j++) //for each of the six vertices of the wedge
342       {
343       this->Wedge->Points->SetPoint(j,this->Points->GetPoint(LinearWedges[i][j]));
344       this->Wedge->PointIds->SetId(j,this->PointIds->GetId(LinearWedges[i][j]));
345       this->Scalars->SetValue(j,cellScalars->GetTuple1(LinearWedges[i][j]));
346       }
347     this->Wedge->Clip(value,this->Scalars,locator,tets,inPd,outPd, inCd,cellId,outCd,insideOut);
348     }
349 }
350 
351 
352 //----------------------------------------------------------------------------
353 // Line-hex intersection. Intersection has to occur within [0,1] parametric
354 // coordinates and with specified tolerance.
IntersectWithLine(double * p1,double * p2,double tol,double & t,double * x,double * pcoords,int & subId)355 int vtkBiQuadraticQuadraticWedge::IntersectWithLine (double *p1, double *p2,
356                                                      double tol, double &t,
357                                                      double *x,
358                                                      double *pcoords,
359                                                      int &subId)
360 {
361   int intersection = 0;
362   double tTemp;
363   double pc[3], xTemp[3];
364   int faceNum;
365   int inter;
366 
367   t = VTK_DOUBLE_MAX;
368   for (faceNum = 0; faceNum < 5; faceNum++)
369     {
370     // We have 9 nodes on biquad face
371     // and 6 on triangle faces
372     if (faceNum < 2)
373       {
374       for (int i = 0; i < 6; i++)
375         {
376         this->TriangleFace->PointIds->SetId (i, this->PointIds->GetId (WedgeFaces[faceNum][i]));
377         this->TriangleFace->Points->SetPoint (i, this->Points->GetPoint (WedgeFaces[faceNum][i]));
378         }
379       inter = this->TriangleFace->IntersectWithLine (p1, p2, tol, tTemp, xTemp,
380                                                      pc, subId);
381       }
382     else
383       {
384       for (int i = 0; i < 9; i++)
385         {
386         this->Face->Points->SetPoint (i, this->Points->GetPoint (WedgeFaces[faceNum][i]));
387         }
388       inter = this->Face->IntersectWithLine (p1, p2, tol, tTemp, xTemp, pc,
389                                              subId);
390       }
391     if (inter)
392       {
393       intersection = 1;
394       if (tTemp < t)
395         {
396         t = tTemp;
397         x[0] = xTemp[0];
398         x[1] = xTemp[1];
399         x[2] = xTemp[2];
400         switch (faceNum)
401           {
402           case 0:
403             pcoords[0] = 0.0;
404             pcoords[1] = pc[1];
405             pcoords[2] = pc[0];
406             break;
407 
408           case 1:
409             pcoords[0] = 1.0;
410             pcoords[1] = pc[0];
411             pcoords[2] = pc[1];
412             break;
413 
414           case 2:
415             pcoords[0] = pc[0];
416             pcoords[1] = 0.0;
417             pcoords[2] = pc[1];
418             break;
419 
420           case 3:
421             pcoords[0] = pc[1];
422             pcoords[1] = 1.0;
423             pcoords[2] = pc[0];
424             break;
425 
426           case 4:
427             pcoords[0] = pc[1];
428             pcoords[1] = pc[0];
429             pcoords[2] = 0.0;
430             break;
431 
432           case 5:
433             pcoords[0] = pc[0];
434             pcoords[1] = pc[1];
435             pcoords[2] = 1.0;
436             break;
437           default:
438             assert("check:impossible case" && 0); // reaching this line is a bug
439             break;
440           }
441         }
442       }
443     }
444   return intersection;
445 }
446 
447 //----------------------------------------------------------------------------
Triangulate(int vtkNotUsed (index),vtkIdList * ptIds,vtkPoints * pts)448 int vtkBiQuadraticQuadraticWedge::Triangulate (int vtkNotUsed (index),
449   vtkIdList * ptIds, vtkPoints * pts)
450 {
451   pts->Reset ();
452   ptIds->Reset ();
453 
454   for (int i = 0; i < 8; i++)
455     {
456     for (int j = 0; j < 6; j++)
457       {
458       ptIds->InsertId (6 * i + j, this->PointIds->GetId (LinearWedges[i][j]));
459       pts->InsertPoint (6 * i + j, this->Points->GetPoint (LinearWedges[i][j]));
460       }
461     }
462 
463   return 1;
464 }
465 
466 //----------------------------------------------------------------------------
467 // Given parametric coordinates compute inverse Jacobian transformation
468 // matrix. Returns 9 elements of 3x3 inverse Jacobian plus interpolation
469 // function derivatives.
JacobianInverse(double pcoords[3],double ** inverse,double derivs[54])470 void vtkBiQuadraticQuadraticWedge::JacobianInverse(double pcoords[3],
471   double **inverse, double derivs[54])
472 {
473   int i, j;
474   double *m[3], m0[3], m1[3], m2[3];
475   double x[3];
476 
477   // compute interpolation function derivatives
478   this->InterpolationDerivs (pcoords, derivs);
479 
480   // create Jacobian matrix
481   m[0] = m0;
482   m[1] = m1;
483   m[2] = m2;
484 
485   for (i = 0; i < 3; i++)  //initialize matrix
486     {
487     m0[i] = m1[i] = m2[i] = 0.0;
488     }
489 
490   for (j = 0; j < 18; j++)
491     {
492     this->Points->GetPoint (j, x);
493     for (i = 0; i < 3; i++)
494       {
495       m0[i] += x[i] * derivs[j];
496       m1[i] += x[i] * derivs[18 + j];
497       m2[i] += x[i] * derivs[36 + j];
498       }
499     }
500 
501   // now find the inverse
502   if (vtkMath::InvertMatrix(m, inverse, 3) == 0)
503     {
504     vtkErrorMacro (<<"Jacobian inverse not found");
505     return;
506     }
507 }
508 
509 //----------------------------------------------------------------------------
Derivatives(int vtkNotUsed (subId),double pcoords[3],double * values,int dim,double * derivs)510 void vtkBiQuadraticQuadraticWedge::Derivatives (int vtkNotUsed (subId),
511   double pcoords[3], double *values, int dim, double *derivs)
512 {
513   double *jI[3], j0[3], j1[3], j2[3];
514   double functionDerivs[3 * 18], sum[3];
515   int i, j, k;
516 
517   // compute inverse Jacobian and interpolation function derivatives
518   jI[0] = j0;
519   jI[1] = j1;
520   jI[2] = j2;
521 
522   this->JacobianInverse (pcoords, jI, functionDerivs);
523 
524   // now compute derivates of values provided
525   for (k = 0; k < dim; k++)  //loop over values per vertex
526     {
527     sum[0] = sum[1] = sum[2] = 0.0;
528     for (i = 0; i < 18; i++)  //loop over interp. function derivatives
529       {
530       sum[0] += functionDerivs[i] * values[dim * i + k];
531       sum[1] += functionDerivs[18 + i] * values[dim * i + k];
532       sum[2] += functionDerivs[36 + i] * values[dim * i + k];
533       }
534     for (j = 0; j < 3; j++)  //loop over derivative directions
535       {
536       derivs[3 * k + j] = sum[0] * jI[j][0] + sum[1] * jI[j][1] + sum[2] * jI[j][2];
537       }
538     }
539 }
540 
541 //----------------------------------------------------------------------------
542 // Compute interpolation functions for the fifteen nodes.
InterpolationFunctions(double pcoords[3],double weights[18])543 void vtkBiQuadraticQuadraticWedge::InterpolationFunctions (double pcoords[3], double weights[18])
544 {
545   // VTK needs parametric coordinates to be between (0,1). Isoparametric
546   // shape functions are formulated between (-1,1). Here we do a
547   // coordinate system conversion from (0,1) to (-1,1).
548   double x = 2*(pcoords[0]- 0.5);
549   double y = 2*(pcoords[1]- 0.5);
550   double z = 2*(pcoords[2]- 0.5);
551 
552   // corners
553   weights[0] =-0.25 * (x + y) * (x + y + 1) * z * (1 - z);
554   weights[1] =-0.25 *  x      * (x + 1)     * z * (1 - z);
555   weights[2] =-0.25 *      y  * (1 + y)     * z * (1 - z);
556   weights[3] = 0.25 * (x + y) * (x + y + 1) * z * (1 + z);
557   weights[4] = 0.25 *  x      * (x + 1)     * z * (1 + z);
558   weights[5] = 0.25 *      y  * (1 + y)     * z * (1 + z);
559 
560   // midsides of quadratic triangles
561   weights[6] =  (x + 1)*(x + y) *  0.5 * z * (1 - z);
562   weights[7] = -(x + 1)*(y + 1) *  0.5 * z * (1 - z);
563   weights[8] =  (y + 1)*(x + y) *  0.5 * z * (1 - z);
564   weights[9] = -(x + 1)*(x + y) *  0.5 * z * (1 + z);
565   weights[10]=  (x + 1)*(y + 1) *  0.5 * z * (1 + z);
566   weights[11]= -(y + 1)*(x + y) *  0.5 * z * (1 + z);
567 
568   // midsides of edges between the two triangles
569   weights[12] = 0.5 * (x + y) * (x + y + 1) * (1 + z)*(1 - z);
570   weights[13] = 0.5 *  x      * (x + 1)     * (1 + z)*(1 - z);
571   weights[14] = 0.5 *      y  * (1 + y)     * (1 + z)*(1 - z);
572 
573   //Centerpoints of the biquadratic quads
574   weights[15] = -(x + 1)*(x + y) * (1 + z)*(1 - z);
575   weights[16] =  (x + 1)*(y + 1) * (1 + z)*(1 - z);
576   weights[17] = -(y + 1)*(x + y) * (1 + z)*(1 - z);
577 }
578 
579 //----------------------------------------------------------------------------
580 // Derivatives in parametric space.
InterpolationDerivs(double pcoords[3],double derivs[54])581 void vtkBiQuadraticQuadraticWedge::InterpolationDerivs (double pcoords[3], double derivs[54])
582 {
583   //VTK needs parametric coordinates to be between (0,1). Isoparametric
584   //shape functions are formulated between (-1,1). Here we do a
585   //coordinate system conversion from (0,1) to (-1,1).
586   double x = 2*(pcoords[0]- 0.5);
587   double y = 2*(pcoords[1]- 0.5);
588   double z = 2*(pcoords[2]- 0.5);
589 
590   //Derivatives in x-direction
591   // corners
592   derivs[0] = -0.25 * (2 * x + 2 * y + 1) * z * (1 - z);
593   derivs[1] = -0.25 * (2 * x + 1)         * z * (1 - z);
594   derivs[2] =  0;
595   derivs[3] =  0.25 * (2 * x + 2 * y + 1) * z * (1 + z);
596   derivs[4] =  0.25 * (2 * x + 1)         * z * (1 + z);
597   derivs[5] =  0;
598   // midsides of quadratic triangles
599   derivs[6] =  (2 * x + y + 1) *  0.5 * z * (1 - z);
600   derivs[7] = -(y + 1)         *  0.5 * z * (1 - z);
601   derivs[8] =  (y + 1)         *  0.5 * z * (1 - z);
602   derivs[9] = -(2 * x + y + 1) *  0.5 * z * (1 + z);
603   derivs[10] = (y + 1)         *  0.5 * z * (1 + z) ;
604   derivs[11] =-(y + 1)         *  0.5 * z * (1 + z) ;
605   // midsides of edges between the two triangles
606   derivs[12] = 0.5 * (2 * x + 2 * y + 1) * (1 + z)*(1 - z);
607   derivs[13] = 0.5 * (2 * x + 1)         * (1 + z)*(1 - z);
608   derivs[14] = 0;
609   //Centerpoints of the biquadratic quads
610   derivs[15] = -(2 * x + y + 1) * (1 + z)*(1 - z);
611   derivs[16] =  (y + 1)         * (1 + z)*(1 - z);
612   derivs[17] = -(y + 1)         * (1 + z)*(1 - z);
613 
614   //Derivatives in y-direction
615   // corners
616   derivs[18] = -0.25 * (2 * y + 2 * x + 1)   * z * (1 - z);
617   derivs[19] =  0;
618   derivs[20] = -0.25 * (2 * y + 1)           * z * (1 - z);
619   derivs[21] =  0.25 * (2 * y + 2 * x + 1)   * z * (1 + z);
620   derivs[22] =  0;
621   derivs[23] =  0.25 * (2 * y + 1)           * z * (1 + z);
622   // midsides of quadratic triangles
623   derivs[24] =  (x + 1)         *  0.5 * z * (1 - z);
624   derivs[25] = -(x + 1)         *  0.5 * z * (1 - z);
625   derivs[26] =  (2 * y + x + 1) *  0.5 * z * (1 - z);
626   derivs[27] = -(x + 1)         *  0.5 * z * (1 + z);
627   derivs[28] =  (x + 1)         *  0.5 * z * (1 + z);
628   derivs[29] = -(2 * y + x + 1) *  0.5 * z * (1 + z);
629   // midsides of edges between the two triangles
630   derivs[30] = 0.5 * (2 * y + 2 * x + 1) * (1 + z)*(1 - z);
631   derivs[31] = 0;
632   derivs[32] = 0.5 * (2 * y + 1)         * (1 + z)*(1 - z);
633   //Centerpoints of the biquadratic quads
634   derivs[33] = -(x + 1)         * (1 + z)*(1 - z);
635   derivs[34] =  (x + 1)         * (1 + z)*(1 - z);
636   derivs[35] = -(2 * y + x + 1) * (1 + z)*(1 - z);
637 
638   //Derivatives in z-direction
639   // corners
640   derivs[36] = -0.25 * (x + y) * (x + y + 1) * (1 - 2 * z);
641   derivs[37] = -0.25 *  x      * (x + 1)     * (1 - 2 * z);
642   derivs[38] = -0.25 *      y  * (1 + y)     * (1 - 2 * z);
643   derivs[39] =  0.25 * (x + y) * (x + y + 1) * (1 + 2 * z);
644   derivs[40] =  0.25 *  x      * (x + 1)     * (1 + 2 * z);
645   derivs[41] =  0.25 *      y  * (1 + y)     * (1 + 2 * z);
646   // midsides of quadratic triangles
647   derivs[42] =  (x + 1)*(x + y) *  0.5 * (1 - 2 * z);
648   derivs[43] = -(x + 1)*(y + 1) *  0.5 * (1 - 2 * z);
649   derivs[44] =  (y + 1)*(x + y) *  0.5 * (1 - 2 * z);
650   derivs[45] = -(x + 1)*(x + y) *  0.5 * (1 + 2 * z);
651   derivs[46] =  (x + 1)*(y + 1) *  0.5 * (1 + 2 * z);
652   derivs[47] = -(y + 1)*(x + y) *  0.5 * (1 + 2 * z);
653   // midsides of edges between the two triangles
654   derivs[48] = 0.5 * (x + y) * (x + y + 1) * (-2 * z);
655   derivs[49] = 0.5 *  x      * (x + 1)     * (-2 * z);
656   derivs[50] = 0.5 *      y  * (1 + y)     * (-2 * z);
657   //Centerpoints of the biquadratic quads
658   derivs[51] = -(x + 1)*(x + y) * (-2 * z);
659   derivs[52] =  (x + 1)*(y + 1) * (-2 * z);
660   derivs[53] = -(y + 1)*(x + y) * (-2 * z);
661 
662   // we compute derivatives in [-1; 1] but we need them in [ 0; 1]
663   for(int i = 0; i < 54; i++)
664     derivs[i] *= 2;
665 }
666 
667 //----------------------------------------------------------------------------
668 static double vtkQWedgeCellPCoords[54] = {
669   0.0,0.0,0.0, 1.0,0.0,0.0, 0.0,1.0,0.0,
670   0.0,0.0,1.0, 1.0,0.0,1.0, 0.0,1.0,1.0,
671   0.5,0.0,0.0, 0.5,0.5,0.0, 0.0,0.5,0.0,
672   0.5,0.0,1.0, 0.5,0.5,1.0, 0.0,0.5,1.0,
673   0.0,0.0,0.5, 1.0,0.0,0.5, 0.0,1.0,0.5,
674   0.5,0.0,0.5, 0.5,0.5,0.5, 0.0,0.5,0.5
675 };
676 
GetParametricCoords()677 double *vtkBiQuadraticQuadraticWedge::GetParametricCoords()
678 {
679   return vtkQWedgeCellPCoords;
680 }
681 
682 //----------------------------------------------------------------------------
PrintSelf(ostream & os,vtkIndent indent)683 void vtkBiQuadraticQuadraticWedge::PrintSelf(ostream & os, vtkIndent indent)
684 {
685   this->Superclass::PrintSelf(os, indent);
686 
687   os << indent << "Edge:\n";
688   this->Edge->PrintSelf (os, indent.GetNextIndent ());
689   os << indent << "TriangleFace:\n";
690   this->TriangleFace->PrintSelf (os, indent.GetNextIndent ());
691   os << indent << "Face:\n";
692   this->Face->PrintSelf (os, indent.GetNextIndent ());
693   os << indent << "Wedge:\n";
694   this->Wedge->PrintSelf (os, indent.GetNextIndent ());
695   os << indent << "Scalars:\n";
696   this->Scalars->PrintSelf (os, indent.GetNextIndent ());
697 }
698