VTK  9.4.20241218
vtkUpdateCellsV8toV9.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
8#ifndef vtkUpdateCellsV8toV9_h
9#define vtkUpdateCellsV8toV9_h
10
11#include "vtkCellArray.h"
12#include "vtkCellData.h"
13#include "vtkCellType.h"
14#include "vtkCellTypes.h"
16#include "vtkIdList.h"
17#include "vtkIdTypeArray.h"
18#include "vtkNew.h"
19#include "vtkUnstructuredGrid.h"
20
21VTK_ABI_NAMESPACE_BEGIN
23{
24 vtkNew<vtkIdList> oldpts, newpts;
25
26 for (vtkIdType i = 0; i < output->GetNumberOfCells(); ++i)
27 {
28 vtkIdType type = output->GetCellTypesArray()->GetTypedComponent(i, 0);
31 {
32 output->GetCells()->GetCellAtId(i, oldpts);
33 newpts->DeepCopy(oldpts);
34
35 int degs[3];
36 if (output->GetCellData()->SetActiveAttribute(
38 {
40 double degs_double[3];
41 v->GetTuple(i, degs_double);
42 for (int ii = 0; ii < 3; ii++)
43 degs[ii] = static_cast<int>(degs_double[ii]);
44 }
45 else
46 {
47 int order =
48 static_cast<int>(round(std::cbrt(static_cast<int>(oldpts->GetNumberOfIds())))) - 1;
49 degs[0] = degs[1] = degs[2] = order;
50 }
51 for (int j = 0; j < oldpts->GetNumberOfIds(); j++)
52 {
54 if (j != newid)
55 {
56 newpts->SetId(j, oldpts->GetId(newid));
57 }
58 }
59 output->GetCells()->ReplaceCellAtId(i, newpts);
60 }
61 }
62}
63
65{
66 vtkIdType nCellTypes = distinctCellTypes->GetNumberOfValues();
67 for (vtkIdType i = 0; i < nCellTypes; ++i)
68 {
69 unsigned char type = distinctCellTypes->GetValue(i);
72 {
73 return true;
74 }
75 }
76 return false;
77}
78
79VTK_ABI_NAMESPACE_END
80#endif // vtkUpdateCellsV8toV9_h
81// VTK-HeaderTest-Exclude: vtkUpdateCellsV8toV9.h
ValueType GetTypedComponent(vtkIdType tupleIdx, int comp) const
Get component comp of the tuple at tupleIdx.
ValueType GetValue(vtkIdType valueIdx) const
Get the value at valueIdx.
vtkIdType GetNumberOfValues() const
Get the total number of values in the array.
void GetCellAtId(vtkIdType cellId, vtkIdType &cellSize, vtkIdType const *&cellPoints, vtkIdList *ptIds) override
Return the point ids for the cell at cellId.
void ReplaceCellAtId(vtkIdType cellId, vtkIdList *list)
Replaces the point ids for the specified cell with the supplied list.
abstract superclass for arrays of numeric data
virtual double * GetTuple(vtkIdType tupleIdx)=0
Get the data tuple at tupleIdx.
int SetActiveAttribute(const char *name, int attributeType)
Make the array with the given name the active attribute.
vtkDataArray * GetHigherOrderDegrees()
Set/Get the rational degrees data.
vtkCellData * GetCellData()
Return a pointer to this dataset's cell data.
Definition vtkDataSet.h:409
static vtkIdType NodeNumberingMappingFromVTK8To9(const int order[3], vtkIdType node_id_vtk8)
Allocate and hold a VTK object.
Definition vtkNew.h:167
dynamic, self-adjusting array of unsigned char
dataset represents arbitrary combinations of all possible cell types
vtkCellArray * GetCells()
Return the unstructured grid connectivity array.
vtkIdType GetNumberOfCells() override
Standard vtkDataSet methods; see vtkDataSet.h for documentation.
vtkUnsignedCharArray * GetCellTypesArray()
Get the array of all cell types in the grid.
Update cells from v8 node layout to v9 node layout.
@ VTK_BEZIER_HEXAHEDRON
@ VTK_LAGRANGE_HEXAHEDRON
@ VTK_HIGHER_ORDER_HEXAHEDRON
Definition vtkCellType.h:98
int vtkIdType
Definition vtkType.h:315
bool vtkNeedsNewFileVersionV8toV9(vtkUnsignedCharArray *distinctCellTypes)