File: vtkQuadraturePointsGenerator.cxx

package info (click to toggle)
vtk7 7.1.1%2Bdfsg1-12
  • links: PTS, VCS
  • area: main
  • in suites: buster
  • size: 125,776 kB
  • sloc: cpp: 1,539,582; ansic: 106,521; python: 78,038; tcl: 47,013; xml: 8,142; yacc: 5,040; java: 4,439; perl: 3,132; lex: 1,926; sh: 1,500; makefile: 122; objc: 83
file content (336 lines) | stat: -rw-r--r-- 10,042 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
/*=========================================================================

  Program:   Visualization Toolkit
  Module:    vtkQuadraturePointsGenerator.cxx

  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
  All rights reserved.
  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.

     This software is distributed WITHOUT ANY WARRANTY; without even
     the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
     PURPOSE.  See the above copyright notice for more information.

=========================================================================*/

#include "vtkQuadraturePointsGenerator.h"

#include "vtkUnstructuredGrid.h"
#include "vtkPolyData.h"
#include "vtkPoints.h"
#include "vtkType.h"
#include "vtkDoubleArray.h"
#include "vtkDataArray.h"
#include "vtkPointData.h"
#include "vtkCellData.h"
#include "vtkCellType.h"
#include "vtkCellArray.h"
#include "vtkIdTypeArray.h"
#include "vtkIntArray.h"
#include "vtkUnstructuredGridAlgorithm.h"
#include "vtkInformation.h"
#include "vtkInformationVector.h"
#include "vtkInformationQuadratureSchemeDefinitionVectorKey.h"
#include "vtkObjectFactory.h"
#include "vtkDataSetAttributes.h"

#include "vtkQuadratureSchemeDefinition.h"
#include "vtkQuadraturePointsUtilities.hxx"

#include <sstream>
using std::ostringstream;

//-----------------------------------------------------------------------------
vtkStandardNewMacro(vtkQuadraturePointsGenerator);

//-----------------------------------------------------------------------------
vtkQuadraturePointsGenerator::vtkQuadraturePointsGenerator()
{
  this->SetNumberOfInputPorts(1);
  this->SetNumberOfOutputPorts(1);
}

//-----------------------------------------------------------------------------
vtkQuadraturePointsGenerator::~vtkQuadraturePointsGenerator()
{}

//-----------------------------------------------------------------------------
int vtkQuadraturePointsGenerator::FillInputPortInformation(
      int vtkNotUsed(port),
      vtkInformation *info)
{
  info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkUnstructuredGrid");
  return 1;
}

//-----------------------------------------------------------------------------
int vtkQuadraturePointsGenerator::RequestData(
      vtkInformation *,
      vtkInformationVector **input,
      vtkInformationVector *output)
{
  vtkDataObject *tmpDataObj;
  // Get the input.
  tmpDataObj
     = input[0]->GetInformationObject(0)->Get(vtkDataObject::DATA_OBJECT());
  vtkUnstructuredGrid *usgIn = vtkUnstructuredGrid::SafeDownCast(tmpDataObj);
  // Get the output.
  tmpDataObj
     = output->GetInformationObject(0)->Get(vtkDataObject::DATA_OBJECT());
  vtkPolyData *pdOut = vtkPolyData::SafeDownCast(tmpDataObj);

  // Quick sanity check.
  if (usgIn == NULL || pdOut == NULL || usgIn->GetNumberOfCells() == 0
      || usgIn->GetNumberOfPoints() == 0 || usgIn->GetCellData() == NULL
      || usgIn->GetCellData()->GetNumberOfArrays() == 0)
  {
    vtkErrorMacro("Filter data has not been configured correctly. Aborting.");
    return 1;
  }

  // Generate points for the selected data array.
  // user specified the offsets array.
  this->Generate(usgIn, this->GetInputArrayToProcess(0, input), pdOut);

  return 1;
}

// ----------------------------------------------------------------------------
int vtkQuadraturePointsGenerator::GenerateField(
        vtkUnstructuredGrid *usgIn,
        vtkDataArray* data,
        vtkDataArray* offsets,
        vtkPolyData* pdOut)
{
  vtkInformation *info = offsets->GetInformation();
  vtkInformationQuadratureSchemeDefinitionVectorKey *key
    = vtkQuadratureSchemeDefinition::DICTIONARY();
  if (!key->Has(info))
  {
    vtkErrorMacro(
      << "Dictionary is not present in array "
      << offsets->GetName() << " " << offsets
     << " Aborting.");
    return 0;
  }

  int dictSize = key->Size(info);
  vtkQuadratureSchemeDefinition **dict
    = new vtkQuadratureSchemeDefinition *[dictSize];
  key->GetRange(info, dict, 0, 0, dictSize);

  vtkIdType nVerts = pdOut->GetNumberOfPoints();

  vtkIdType cellId;
  vtkIdType ncell = usgIn->GetNumberOfCells();
  int cellType;
  // first loop through all cells to check if a shallow copy is possible
  bool shallowok = true;
  vtkIdType previous = -1;

  int offsetType = offsets->GetDataType();
  void *pOffsets = offsets->GetVoidPointer(0);
  switch(offsetType)
  {
    vtkTemplateMacro(
      for (cellId = 0; cellId < ncell; cellId++)
      {
        vtkIdType offset = static_cast<vtkIdType>(((VTK_TT*)pOffsets)[cellId]);
        if (offset != previous + 1)
        {
          shallowok = false;
          break;
        }
        cellType = usgIn->GetCellType(cellId);

        if (dict[cellType] == NULL)
        {
          previous = offset;
        }
        else
        {
          previous = offset + dict[cellType]->GetNumberOfQuadraturePoints();
        }
      }

      if ( (previous+1) != nVerts )
      {
        shallowok = false;
      }

      if (shallowok)
      {
        // ok, all the original cells are here, we can shallow copy the array
        // from input to output
        pdOut->GetPointData()->AddArray(data);
      }
      else
      {
        // in this case, we have to duplicate the valid tuples
        vtkDataArray *V_out = data->NewInstance();
        V_out->SetName(data->GetName());
        V_out->SetNumberOfComponents(data->GetNumberOfComponents());
        V_out->CopyComponentNames( data );
        for (cellId = 0; cellId < ncell; cellId++)
        {
          vtkIdType offset = static_cast<vtkIdType>(((VTK_TT*)pOffsets)[cellId]);
          cellType = usgIn->GetCellType(cellId);

          // a simple check to see if a scheme really exists for this cell type.
          // should not happen if the cell type has not been modified.
          if (dict[cellType] == NULL)
          {
            continue;
          }

          int np = dict[cellType]->GetNumberOfQuadraturePoints();
          for (int id = 0; id < np; id++)
          {
            V_out->InsertNextTuple(offset + id, data);
          }
        }
        V_out->Squeeze();
        pdOut->GetPointData()->AddArray(V_out);
        V_out->Delete();
      }
    );
  }
  delete[] dict;
  return 1;
}

//-----------------------------------------------------------------------------
int vtkQuadraturePointsGenerator::Generate(
        vtkUnstructuredGrid *usgIn,
        vtkDataArray* offsets,
        vtkPolyData *pdOut)
{
  if (usgIn == NULL || offsets == NULL || pdOut == NULL)
  {
    vtkErrorMacro("configuration error");
    return 0;
  }

  // Strategy:
  // create the points then move the FieldData to PointData

  const char* offsetName = offsets->GetName();
  if (offsetName == NULL)
  {
    vtkErrorMacro("offset array has no name, Skipping");
    return 1;
  }

  // get the dictionnary
  vtkInformation *info = offsets->GetInformation();
  vtkInformationQuadratureSchemeDefinitionVectorKey *key =
      vtkQuadratureSchemeDefinition::DICTIONARY();
  if (!key->Has(info))
  {
    vtkErrorMacro(
      << "Dictionary is not present in array "
      << offsets->GetName()
      << " Aborting.");
    return 0;
  }
  int dictSize = key->Size(info);
  vtkQuadratureSchemeDefinition **dict =
      new vtkQuadratureSchemeDefinition *[dictSize];
  key->GetRange(info, dict, 0, 0, dictSize);

  // Grab the point set.
  vtkDataArray *X = usgIn->GetPoints()->GetData();
  int X_type = X->GetDataType();
  void *pX = X->GetVoidPointer(0);

  // Create the result array.
  vtkDoubleArray *qPts = vtkDoubleArray::New();
  vtkIdType nCells = usgIn->GetNumberOfCells();
  qPts->Allocate(3* nCells ); // Expect at least one point per cell
  qPts->SetNumberOfComponents(3);

  // For all cells interpolate.
  switch (X_type)
  {
    vtkTemplateMacro(
      if (!Interpolate(
              usgIn,
              nCells,
              (VTK_TT*)pX,
              3,
              dict,
              qPts,
              (int*)NULL))
      {
        vtkWarningMacro("Failed to interpolate cell vertices "
            "to quadrature points. Aborting.");
      }
      );
  }
  delete[] dict;

  // Add the interpolated quadrature points to the output
  vtkIdType nVerts = qPts->GetNumberOfTuples();
  vtkPoints *p = vtkPoints::New();
  p->SetDataTypeToDouble();
  p->SetData(qPts);
  qPts->Delete();
  pdOut->SetPoints(p);
  p->Delete();
  // Generate vertices at the quadrature points
  vtkIdTypeArray *va = vtkIdTypeArray::New();
  va->SetNumberOfTuples(2* nVerts );
  vtkIdType *verts = va->GetPointer(0);
  for (int i = 0; i < nVerts; ++i)
  {
    verts[0] = 1;
    verts[1] = i;
    verts += 2;
  }
  vtkCellArray *cells = vtkCellArray::New();
  cells->SetCells(static_cast<vtkIdType> (nVerts), va);
  pdOut->SetVerts(cells);
  cells->Delete();
  va->Delete();

  // then loop over all fields to map the field array to the points
  int nArrays = usgIn->GetFieldData()->GetNumberOfArrays();
  for (int i = 0; i<nArrays; ++i)
  {
    vtkDataArray* array = usgIn->GetFieldData()->GetArray(i);
    if (array == NULL) continue;

    const char* arrayOffsetName = array->GetInformation()->Get(
        vtkQuadratureSchemeDefinition::QUADRATURE_OFFSET_ARRAY_NAME());
    if (arrayOffsetName == NULL)
    {
      // not an error, since non-quadrature point field data may
      //  be present.
      vtkDebugMacro(
        << "array " << array->GetName()
        << " has no offset array name, Skipping");
      continue;
    }

    if (strcmp(offsetName, arrayOffsetName) != 0)
    {
      // not an error, this array does not belong with the current
      // quadrature scheme definition.
      vtkDebugMacro(
        << "array " << array->GetName()
        << " has another offset array : "
        << arrayOffsetName << ", Skipping");
      continue;
    }

    this->GenerateField(usgIn, array, offsets, pdOut);
  }

  return 1;
}

//-----------------------------------------------------------------------------
void vtkQuadraturePointsGenerator::PrintSelf(ostream& os, vtkIndent indent)
{
  this->Superclass::PrintSelf(os, indent);
}