1 /*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkQuadraticQuad.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 "vtkQuadraticQuad.h"
16
17 #include "vtkCellData.h"
18 #include "vtkDoubleArray.h"
19 #include "vtkMath.h"
20 #include "vtkObjectFactory.h"
21 #include "vtkPointData.h"
22 #include "vtkPoints.h"
23 #include "vtkQuad.h"
24 #include "vtkQuadraticEdge.h"
25
26 vtkStandardNewMacro(vtkQuadraticQuad);
27
28 //------------------------------------------------------------------------------
29 // Construct the quad with eight points.
vtkQuadraticQuad()30 vtkQuadraticQuad::vtkQuadraticQuad()
31 {
32 this->Edge = vtkQuadraticEdge::New();
33 this->Quad = vtkQuad::New();
34 this->PointData = vtkPointData::New();
35 this->CellData = vtkCellData::New();
36 this->CellScalars = vtkDoubleArray::New();
37 this->CellScalars->SetNumberOfTuples(9);
38 this->Scalars = vtkDoubleArray::New();
39 this->Scalars->SetNumberOfTuples(4);
40
41 // We add a fictitious ninth point in order to process the cell. The ninth
42 // point is in the center of the cell.
43 this->Points->SetNumberOfPoints(9);
44 this->PointIds->SetNumberOfIds(9);
45 for (int i = 0; i < 9; i++)
46 {
47 this->Points->SetPoint(i, 0.0, 0.0, 0.0);
48 this->PointIds->SetId(i, 0);
49 }
50 this->Points->SetNumberOfPoints(8);
51 this->PointIds->SetNumberOfIds(8);
52 }
53
54 //------------------------------------------------------------------------------
~vtkQuadraticQuad()55 vtkQuadraticQuad::~vtkQuadraticQuad()
56 {
57 this->Edge->Delete();
58 this->Quad->Delete();
59
60 this->Scalars->Delete();
61 this->PointData->Delete();
62 this->CellData->Delete();
63 this->CellScalars->Delete();
64 }
65 //------------------------------------------------------------------------------
GetEdge(int edgeId)66 vtkCell* vtkQuadraticQuad::GetEdge(int edgeId)
67 {
68 edgeId = (edgeId < 0 ? 0 : (edgeId > 3 ? 3 : edgeId));
69 int p = (edgeId + 1) % 4;
70
71 // load point id's
72 this->Edge->PointIds->SetId(0, this->PointIds->GetId(edgeId));
73 this->Edge->PointIds->SetId(1, this->PointIds->GetId(p));
74 this->Edge->PointIds->SetId(2, this->PointIds->GetId(edgeId + 4));
75
76 // load coordinates
77 this->Edge->Points->SetPoint(0, this->Points->GetPoint(edgeId));
78 this->Edge->Points->SetPoint(1, this->Points->GetPoint(p));
79 this->Edge->Points->SetPoint(2, this->Points->GetPoint(edgeId + 4));
80
81 return this->Edge;
82 }
83
84 //------------------------------------------------------------------------------
85
86 static int LinearQuads[4][4] = {
87 { 0, 4, 8, 7 },
88 { 4, 1, 5, 8 },
89 { 8, 5, 2, 6 },
90 { 7, 8, 6, 3 },
91 };
92
Subdivide(double * weights)93 void vtkQuadraticQuad::Subdivide(double* weights)
94 {
95 int i, j;
96 double pc[3], x[3];
97
98 pc[0] = pc[1] = 0.5;
99 vtkQuadraticQuad::InterpolationFunctions(pc, weights);
100
101 double p[3];
102 x[0] = x[1] = x[2] = 0.0;
103 for (i = 0; i < 8; i++)
104 {
105 this->Points->GetPoint(i, p);
106 for (j = 0; j < 3; j++)
107 {
108 x[j] += p[j] * weights[i];
109 }
110 }
111 this->Points->SetPoint(8, x);
112 }
113
114 //------------------------------------------------------------------------------
EvaluatePosition(const double * x,double closestPoint[3],int & subId,double pcoords[3],double & minDist2,double weights[])115 int vtkQuadraticQuad::EvaluatePosition(const double* x, double closestPoint[3], int& subId,
116 double pcoords[3], double& minDist2, double weights[])
117 {
118 double pc[3], dist2;
119 int ignoreId, i, returnStatus = 0, status;
120 double tempWeights[4];
121 double closest[3];
122
123 // compute the midquad node
124 this->Subdivide(weights);
125
126 // four linear quads are used
127 for (minDist2 = VTK_DOUBLE_MAX, i = 0; i < 4; i++)
128 {
129 this->Quad->Points->SetPoint(0, this->Points->GetPoint(LinearQuads[i][0]));
130 this->Quad->Points->SetPoint(1, this->Points->GetPoint(LinearQuads[i][1]));
131 this->Quad->Points->SetPoint(2, this->Points->GetPoint(LinearQuads[i][2]));
132 this->Quad->Points->SetPoint(3, this->Points->GetPoint(LinearQuads[i][3]));
133
134 status = this->Quad->EvaluatePosition(x, closest, ignoreId, pc, dist2, tempWeights);
135 if (status != -1 && dist2 < minDist2)
136 {
137 returnStatus = status;
138 minDist2 = dist2;
139 subId = i;
140 pcoords[0] = pc[0];
141 pcoords[1] = pc[1];
142 }
143 }
144
145 // adjust parametric coordinates
146 if (returnStatus != -1)
147 {
148 if (subId == 0)
149 {
150 pcoords[0] /= 2.0;
151 pcoords[1] /= 2.0;
152 }
153 else if (subId == 1)
154 {
155 pcoords[0] = 0.5 + (pcoords[0] / 2.0);
156 pcoords[1] /= 2.0;
157 }
158 else if (subId == 2)
159 {
160 pcoords[0] = 0.5 + (pcoords[0] / 2.0);
161 pcoords[1] = 0.5 + (pcoords[1] / 2.0);
162 }
163 else
164 {
165 pcoords[0] /= 2.0;
166 pcoords[1] = 0.5 + (pcoords[1] / 2.0);
167 }
168 pcoords[2] = 0.0;
169 if (closestPoint != nullptr)
170 {
171 // Compute both closestPoint and weights
172 this->EvaluateLocation(subId, pcoords, closestPoint, weights);
173 }
174 else
175 {
176 // Compute weights only
177 vtkQuadraticQuad::InterpolationFunctions(pcoords, weights);
178 }
179 }
180
181 return returnStatus;
182 }
183 //------------------------------------------------------------------------------
EvaluateLocation(int & vtkNotUsed (subId),const double pcoords[3],double x[3],double * weights)184 void vtkQuadraticQuad::EvaluateLocation(
185 int& vtkNotUsed(subId), const double pcoords[3], double x[3], double* weights)
186 {
187 int i, j;
188 double pt[3];
189
190 vtkQuadraticQuad::InterpolationFunctions(pcoords, weights);
191
192 x[0] = x[1] = x[2] = 0.0;
193 for (i = 0; i < 8; i++)
194 {
195 this->Points->GetPoint(i, pt);
196 for (j = 0; j < 3; j++)
197 {
198 x[j] += pt[j] * weights[i];
199 }
200 }
201 }
202
203 //------------------------------------------------------------------------------
CellBoundary(int subId,const double pcoords[3],vtkIdList * pts)204 int vtkQuadraticQuad::CellBoundary(int subId, const double pcoords[3], vtkIdList* pts)
205 {
206 return this->Quad->CellBoundary(subId, pcoords, pts);
207 }
208
209 static double MidPoints[1][3] = { { 0.5, 0.5, 0.0 } };
210
211 //------------------------------------------------------------------------------
InterpolateAttributes(vtkPointData * inPd,vtkCellData * inCd,vtkIdType cellId,vtkDataArray * cellScalars)212 void vtkQuadraticQuad::InterpolateAttributes(
213 vtkPointData* inPd, vtkCellData* inCd, vtkIdType cellId, vtkDataArray* cellScalars)
214 {
215 int numMidPts, i, j;
216 double weights[20];
217 double x[3];
218 double s;
219
220 // Copy point and cell attribute data, first make sure it's empty:
221 this->PointData->Initialize();
222 this->CellData->Initialize();
223 // Make sure to copy ALL arrays. These field data have to be
224 // identical to the input field data. Otherwise, CopyData
225 // that occurs later may not work because the output field
226 // data was initialized (CopyAllocate) with the input field
227 // data.
228 this->PointData->CopyAllOn();
229 this->CellData->CopyAllOn();
230 this->PointData->CopyAllocate(inPd, 9);
231 this->CellData->CopyAllocate(inCd, 4);
232
233 // copy the point data over into point ids 0->7
234 for (i = 0; i < 8; i++)
235 {
236 this->PointData->CopyData(inPd, this->PointIds->GetId(i), i);
237 this->CellScalars->SetValue(i, cellScalars->GetTuple1(i));
238 }
239 // copy the cell data over to the linear cell
240 this->CellData->CopyData(inCd, cellId, 0);
241
242 // Interpolate new values
243 double p[3];
244 this->Points->Resize(9);
245 this->CellScalars->Resize(9);
246 for (numMidPts = 0; numMidPts < 1; numMidPts++)
247 {
248 vtkQuadraticQuad::InterpolationFunctions(MidPoints[numMidPts], weights);
249
250 x[0] = x[1] = x[2] = 0.0;
251 s = 0.0;
252 for (i = 0; i < 8; i++)
253 {
254 this->Points->GetPoint(i, p);
255 for (j = 0; j < 3; j++)
256 {
257 x[j] += p[j] * weights[i];
258 }
259 s += cellScalars->GetTuple1(i) * weights[i];
260 }
261 this->Points->SetPoint(8 + numMidPts, x);
262 this->CellScalars->SetValue(8 + numMidPts, s);
263 this->PointData->InterpolatePoint(inPd, 8 + numMidPts, this->PointIds, weights);
264 }
265 }
266
267 //------------------------------------------------------------------------------
Contour(double value,vtkDataArray * cellScalars,vtkIncrementalPointLocator * locator,vtkCellArray * verts,vtkCellArray * lines,vtkCellArray * polys,vtkPointData * inPd,vtkPointData * outPd,vtkCellData * inCd,vtkIdType cellId,vtkCellData * outCd)268 void vtkQuadraticQuad::Contour(double value, vtkDataArray* cellScalars,
269 vtkIncrementalPointLocator* locator, vtkCellArray* verts, vtkCellArray* lines,
270 vtkCellArray* polys, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId,
271 vtkCellData* outCd)
272 {
273 // interpolate point and cell data
274 this->InterpolateAttributes(inPd, inCd, cellId, cellScalars);
275
276 // contour each linear quad separately
277 for (int i = 0; i < 4; i++)
278 {
279 for (int j = 0; j < 4; j++) // for each of the four vertices of the linear quad
280 {
281 this->Quad->Points->SetPoint(j, this->Points->GetPoint(LinearQuads[i][j]));
282 this->Quad->PointIds->SetId(j, LinearQuads[i][j]);
283 this->Scalars->SetValue(j, this->CellScalars->GetValue(LinearQuads[i][j]));
284 }
285
286 this->Quad->Contour(value, this->Scalars, locator, verts, lines, polys, this->PointData, outPd,
287 this->CellData, i, outCd);
288 }
289 }
290 //------------------------------------------------------------------------------
291 // Clip this quadratic quad using scalar value provided. Like contouring,
292 // except that it cuts the quad to produce other quads and triangles.
Clip(double value,vtkDataArray * cellScalars,vtkIncrementalPointLocator * locator,vtkCellArray * polys,vtkPointData * inPd,vtkPointData * outPd,vtkCellData * inCd,vtkIdType cellId,vtkCellData * outCd,int insideOut)293 void vtkQuadraticQuad::Clip(double value, vtkDataArray* cellScalars,
294 vtkIncrementalPointLocator* locator, vtkCellArray* polys, vtkPointData* inPd, vtkPointData* outPd,
295 vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd, int insideOut)
296 {
297 // interpolate point and cell data
298 this->InterpolateAttributes(inPd, inCd, cellId, cellScalars);
299
300 // contour each linear quad separately
301 for (int i = 0; i < 4; i++)
302 {
303 for (int j = 0; j < 4; j++) // for each of the four vertices of the linear quad
304 {
305 this->Quad->Points->SetPoint(j, this->Points->GetPoint(LinearQuads[i][j]));
306 this->Quad->PointIds->SetId(j, LinearQuads[i][j]);
307 this->Scalars->SetValue(j, this->CellScalars->GetValue(LinearQuads[i][j]));
308 }
309
310 this->Quad->Clip(value, this->Scalars, locator, polys, this->PointData, outPd, this->CellData,
311 i, outCd, insideOut);
312 }
313 }
314
315 //------------------------------------------------------------------------------
316 // Line-line intersection. Intersection has to occur within [0,1] parametric
317 // coordinates and with specified tolerance.
IntersectWithLine(const double * p1,const double * p2,double tol,double & t,double * x,double * pcoords,int & subId)318 int vtkQuadraticQuad::IntersectWithLine(
319 const double* p1, const double* p2, double tol, double& t, double* x, double* pcoords, int& subId)
320 {
321 int subTest, i;
322 subId = 0;
323 double weights[8];
324
325 // first define the midquad point
326 this->Subdivide(weights);
327
328 // intersect the four linear quads
329 for (i = 0; i < 4; i++)
330 {
331 this->Quad->Points->SetPoint(0, this->Points->GetPoint(LinearQuads[i][0]));
332 this->Quad->Points->SetPoint(1, this->Points->GetPoint(LinearQuads[i][1]));
333 this->Quad->Points->SetPoint(2, this->Points->GetPoint(LinearQuads[i][2]));
334 this->Quad->Points->SetPoint(3, this->Points->GetPoint(LinearQuads[i][3]));
335
336 if (this->Quad->IntersectWithLine(p1, p2, tol, t, x, pcoords, subTest))
337 {
338 return 1;
339 }
340 }
341
342 return 0;
343 }
344
345 //------------------------------------------------------------------------------
Triangulate(int vtkNotUsed (index),vtkIdList * ptIds,vtkPoints * pts)346 int vtkQuadraticQuad::Triangulate(int vtkNotUsed(index), vtkIdList* ptIds, vtkPoints* pts)
347 {
348 pts->Reset();
349 ptIds->Reset();
350
351 // Create six linear triangles: one at each corner and two
352 // to cover the remaining quadrilateral.
353
354 // First the corner vertices
355 ptIds->InsertId(0, this->PointIds->GetId(0));
356 ptIds->InsertId(1, this->PointIds->GetId(4));
357 ptIds->InsertId(2, this->PointIds->GetId(7));
358 pts->InsertPoint(0, this->Points->GetPoint(0));
359 pts->InsertPoint(1, this->Points->GetPoint(4));
360 pts->InsertPoint(2, this->Points->GetPoint(7));
361
362 ptIds->InsertId(3, this->PointIds->GetId(4));
363 ptIds->InsertId(4, this->PointIds->GetId(1));
364 ptIds->InsertId(5, this->PointIds->GetId(5));
365 pts->InsertPoint(3, this->Points->GetPoint(4));
366 pts->InsertPoint(4, this->Points->GetPoint(1));
367 pts->InsertPoint(5, this->Points->GetPoint(5));
368
369 ptIds->InsertId(6, this->PointIds->GetId(5));
370 ptIds->InsertId(7, this->PointIds->GetId(2));
371 ptIds->InsertId(8, this->PointIds->GetId(6));
372 pts->InsertPoint(6, this->Points->GetPoint(5));
373 pts->InsertPoint(7, this->Points->GetPoint(2));
374 pts->InsertPoint(8, this->Points->GetPoint(6));
375
376 ptIds->InsertId(9, this->PointIds->GetId(6));
377 ptIds->InsertId(10, this->PointIds->GetId(3));
378 ptIds->InsertId(11, this->PointIds->GetId(7));
379 pts->InsertPoint(9, this->Points->GetPoint(6));
380 pts->InsertPoint(10, this->Points->GetPoint(3));
381 pts->InsertPoint(11, this->Points->GetPoint(7));
382
383 // Now the two remaining triangles
384 // Choose the triangulation that minimizes the edge length
385 // across the cell.
386 double x4[3], x5[3], x6[3], x7[3];
387 this->Points->GetPoint(4, x4);
388 this->Points->GetPoint(5, x5);
389 this->Points->GetPoint(6, x6);
390 this->Points->GetPoint(7, x7);
391
392 if (vtkMath::Distance2BetweenPoints(x4, x6) <= vtkMath::Distance2BetweenPoints(x5, x7))
393 {
394 ptIds->InsertId(12, this->PointIds->GetId(4));
395 ptIds->InsertId(13, this->PointIds->GetId(6));
396 ptIds->InsertId(14, this->PointIds->GetId(7));
397 pts->InsertPoint(12, this->Points->GetPoint(4));
398 pts->InsertPoint(13, this->Points->GetPoint(6));
399 pts->InsertPoint(14, this->Points->GetPoint(7));
400
401 ptIds->InsertId(15, this->PointIds->GetId(4));
402 ptIds->InsertId(16, this->PointIds->GetId(5));
403 ptIds->InsertId(17, this->PointIds->GetId(6));
404 pts->InsertPoint(15, this->Points->GetPoint(4));
405 pts->InsertPoint(16, this->Points->GetPoint(5));
406 pts->InsertPoint(17, this->Points->GetPoint(6));
407 }
408 else
409 {
410 ptIds->InsertId(12, this->PointIds->GetId(5));
411 ptIds->InsertId(13, this->PointIds->GetId(6));
412 ptIds->InsertId(14, this->PointIds->GetId(7));
413 pts->InsertPoint(12, this->Points->GetPoint(5));
414 pts->InsertPoint(13, this->Points->GetPoint(6));
415 pts->InsertPoint(14, this->Points->GetPoint(7));
416
417 ptIds->InsertId(15, this->PointIds->GetId(5));
418 ptIds->InsertId(16, this->PointIds->GetId(7));
419 ptIds->InsertId(17, this->PointIds->GetId(4));
420 pts->InsertPoint(15, this->Points->GetPoint(5));
421 pts->InsertPoint(16, this->Points->GetPoint(7));
422 pts->InsertPoint(17, this->Points->GetPoint(4));
423 }
424
425 return 1;
426 }
427
428 //------------------------------------------------------------------------------
Derivatives(int vtkNotUsed (subId),const double pcoords[3],const double * values,int dim,double * derivs)429 void vtkQuadraticQuad::Derivatives(
430 int vtkNotUsed(subId), const double pcoords[3], const double* values, int dim, double* derivs)
431 {
432 double sum[2], p[3];
433 double functionDerivs[16];
434 double *J[3], J0[3], J1[3], J2[3];
435 double *JI[3], JI0[3], JI1[3], JI2[3];
436
437 vtkQuadraticQuad::InterpolationDerivs(pcoords, functionDerivs);
438
439 // Compute transposed Jacobian and inverse Jacobian
440 J[0] = J0;
441 J[1] = J1;
442 J[2] = J2;
443 JI[0] = JI0;
444 JI[1] = JI1;
445 JI[2] = JI2;
446 for (int k = 0; k < 3; k++)
447 {
448 J0[k] = J1[k] = 0.0;
449 }
450
451 for (int i = 0; i < 8; i++)
452 {
453 this->Points->GetPoint(i, p);
454 for (int j = 0; j < 2; j++)
455 {
456 for (int k = 0; k < 3; k++)
457 {
458 J[j][k] += p[k] * functionDerivs[j * 8 + i];
459 }
460 }
461 }
462
463 // Compute third row vector in transposed Jacobian and normalize it, so that Jacobian determinant
464 // stays the same.
465 vtkMath::Cross(J0, J1, J2);
466 if (vtkMath::Normalize(J2) == 0.0 || !vtkMath::InvertMatrix(J, JI, 3)) // degenerate
467 {
468 for (int j = 0; j < dim; j++)
469 {
470 for (int i = 0; i < 3; i++)
471 {
472 derivs[j * dim + i] = 0.0;
473 }
474 }
475 return;
476 }
477
478 // Loop over "dim" derivative values. For each set of values,
479 // compute derivatives
480 // in local system and then transform into modelling system.
481 // First compute derivatives in local x'-y' coordinate system
482 for (int j = 0; j < dim; j++)
483 {
484 sum[0] = sum[1] = 0.0;
485 for (int i = 0; i < 8; i++) // loop over interp. function derivatives
486 {
487 sum[0] += functionDerivs[i] * values[dim * i + j];
488 sum[1] += functionDerivs[8 + i] * values[dim * i + j];
489 }
490 // dBydx = sum[0]*JI[0][0] + sum[1]*JI[0][1];
491 // dBydy = sum[0]*JI[1][0] + sum[1]*JI[1][1];
492
493 // Transform into global system (dot product with global axes)
494 derivs[3 * j] = sum[0] * JI[0][0] + sum[1] * JI[0][1];
495 derivs[3 * j + 1] = sum[0] * JI[1][0] + sum[1] * JI[1][1];
496 derivs[3 * j + 2] = sum[0] * JI[2][0] + sum[1] * JI[2][1];
497 }
498 }
499
500 //------------------------------------------------------------------------------
501 // Compute interpolation functions. The first four nodes are the corner
502 // vertices; the others are mid-edge nodes.
InterpolationFunctions(const double pcoords[3],double weights[8])503 void vtkQuadraticQuad::InterpolationFunctions(const double pcoords[3], double weights[8])
504 {
505 double r = pcoords[0];
506 double s = pcoords[1];
507
508 // midedge weights
509 weights[4] = 4 * r * (1.0 - r) * (1.0 - s);
510 weights[5] = 4 * r * (1.0 - s) * s;
511 weights[6] = 4 * r * (1.0 - r) * s;
512 weights[7] = 4 * (1.0 - r) * (1.0 - s) * s;
513
514 // corner
515 weights[0] = (1.0 - r) * (1.0 - s) - 0.5 * (weights[4] + weights[7]);
516 weights[1] = r * (1.0 - s) - 0.5 * (weights[4] + weights[5]);
517 weights[2] = r * s - 0.5 * (weights[5] + weights[6]);
518 weights[3] = (1.0 - r) * s - 0.5 * (weights[6] + weights[7]);
519 }
520
521 //------------------------------------------------------------------------------
522 // Derivatives in parametric space.
InterpolationDerivs(const double pcoords[3],double derivs[16])523 void vtkQuadraticQuad::InterpolationDerivs(const double pcoords[3], double derivs[16])
524 {
525 // Coordinate conversion
526 double r = pcoords[0];
527 double s = pcoords[1];
528
529 // Derivatives in the r-direction
530 // midedge
531 derivs[4] = 4 * (1.0 - s) * (1.0 - 2 * r);
532 derivs[5] = 4 * (1.0 - s) * s;
533 derivs[6] = 4 * s * (1.0 - 2 * r);
534 derivs[7] = -4 * (1.0 - s) * s;
535 derivs[0] = -(1.0 - s) - 0.5 * (derivs[4] + derivs[7]);
536 derivs[1] = (1.0 - s) - 0.5 * (derivs[4] + derivs[5]);
537 derivs[2] = s - 0.5 * (derivs[5] + derivs[6]);
538 derivs[3] = -s - 0.5 * (derivs[6] + derivs[7]);
539
540 // Derivatives in the s-direction
541 // midedge
542 derivs[12] = -4 * r * (1.0 - r);
543 derivs[13] = 4 * r * (1.0 - 2 * s);
544 derivs[14] = 4 * r * (1.0 - r);
545 derivs[15] = 4 * (1.0 - r) * (1.0 - 2 * s);
546 derivs[8] = -(1.0 - r) - 0.5 * (derivs[12] + derivs[15]);
547 derivs[9] = -r - 0.5 * (derivs[12] + derivs[13]);
548 derivs[10] = r - 0.5 * (derivs[13] + derivs[14]);
549 derivs[11] = (1.0 - r) - 0.5 * (derivs[14] + derivs[15]);
550 }
551
552 //------------------------------------------------------------------------------
553 static double vtkQQuadCellPCoords[24] = {
554 0.0, 0.0, 0.0, //
555 1.0, 0.0, 0.0, //
556 1.0, 1.0, 0.0, //
557 0.0, 1.0, 0.0, //
558 0.5, 0.0, 0.0, //
559 1.0, 0.5, 0.0, //
560 0.5, 1.0, 0.0, //
561 0.0, 0.5, 0.0 //
562 };
GetParametricCoords()563 double* vtkQuadraticQuad::GetParametricCoords()
564 {
565 return vtkQQuadCellPCoords;
566 }
567
568 //------------------------------------------------------------------------------
PrintSelf(ostream & os,vtkIndent indent)569 void vtkQuadraticQuad::PrintSelf(ostream& os, vtkIndent indent)
570 {
571 this->Superclass::PrintSelf(os, indent);
572
573 os << indent << "Edge:\n";
574 this->Edge->PrintSelf(os, indent.GetNextIndent());
575 os << indent << "Quad:\n";
576 this->Quad->PrintSelf(os, indent.GetNextIndent());
577 os << indent << "Scalars:\n";
578 this->Scalars->PrintSelf(os, indent.GetNextIndent());
579 }
580