[vtkusers] bug or intentional?

Philipp.Batchelor Philipp.Batchelor at kcl.ac.uk
Fri Jun 8 10:34:34 EDT 2001


Hello,
I have been getting some strange results from vtkSubdivideTetra, namely,
starting with the simplest cube made of 5 vtkTetras (4 congruent simplices, one 
regular) see 'cube.ug.vtk', I get as output 43 points (60=5*12 cells), and 
clearly some points are duplicates, cf points 12 and 17 in attached datafile 
'cube.ug.r1.vtk', marked by >>...<<.

Looking at the source of vtkSubdivideTetra, the 'vtkMergePoints locator' is 
using the method :
'locator->InsertNextPoint(...)'
where I would rather expect:
'locator->InsertUniquePoint(...)'.
I modified vtkSubdivideTetra in that way, and afterwards get only 31 points, cf 
'cube.ug.r1_clean.vtk', which is much more satisfying for my application.

Is is intentional to use 'locator->InsertNextPoint'? or should 
vtkSubdivideTetra be modified? In case of interest, my modified version is 
attached too.
Sorry if the question is stupid.
Ph
-------------- next part --------------
# vtk DataFile Version 3.0
vtk output
ASCII
DATASET UNSTRUCTURED_GRID
POINTS 8 float
0 0 0 1 0 1 1 1 0
0 1 1 0 0 1 1 1 1
1 0 0 0 1 0
CELLS 5 25
4 0 1 2 3
4 1 0 3 4
4 1 3 5 2
4 0 6 2 1
4 0 2 7 3
CELL_TYPES 5
10
10
10
10
10
-------------- next part --------------
# vtk DataFile Version 3.0
vtk output
ASCII
DATASET UNSTRUCTURED_GRID
POINTS 43 float
0 0 0 1 0 1 1 1 0 
0 1 1 0 0 1 1 1 1 
1 0 0 0 1 0 0.5 0.5 0.5 
0.5 0 0.5 1 0.5 0.5 0.5 0.5 0 
>>0 0.5 0.5<< 0.5 0.5 1 0.5 1 0.5 
0.25 0.25 0.75 0.5 0 0.5 >>0 0.5 0.5<< 
0.5 0.5 1 0.5 0 1 0 0 0.5 
0 0.5 1 0.75 0.75 0.75 0.5 0.5 1 
0.5 1 1 1 0.5 1 1 0.5 0.5 
0.5 1 0.5 1 1 0.5 0.75 0.25 0.25 
0.5 0 0 1 0.5 0 0.5 0.5 0 
0.5 0 0.5 1 0 0.5 1 0.5 0.5 
0.25 0.75 0.25 0.5 0.5 0 0.5 1 0 
0 0.5 0 0 0.5 0.5 0.5 1 0.5 
0 1 0.5 
CELLS 60 300
4 0 9 11 12 
4 1 9 10 13 
4 2 11 10 14 
4 3 12 13 14 
4 8 9 11 12 
4 8 9 10 13 
4 8 11 10 14 
4 8 12 13 14 
4 8 9 10 11 
4 8 9 13 12 
4 8 10 14 13 
4 8 11 14 12 
4 1 16 18 19 
4 0 16 17 20 
4 3 18 17 21 
4 4 19 20 21 
4 15 16 18 19 
4 15 16 17 20 
4 15 18 17 21 
4 15 19 20 21 
4 15 16 17 18 
4 15 16 20 19 
4 15 17 21 20 
4 15 18 21 19 
4 1 23 25 26 
4 3 23 24 27 
4 5 25 24 28 
4 2 26 27 28 
4 22 23 25 26 
4 22 23 24 27 
4 22 25 24 28 
4 22 26 27 28 
4 22 23 24 25 
4 22 23 27 26 
4 22 24 28 27 
4 22 25 28 26 
4 0 30 32 33 
4 6 30 31 34 
4 2 32 31 35 
4 1 33 34 35 
4 29 30 32 33 
4 29 30 31 34 
4 29 32 31 35 
4 29 33 34 35 
4 29 30 31 32 
4 29 30 34 33 
4 29 31 35 34 
4 29 32 35 33 
4 0 37 39 40 
4 2 37 38 41 
4 7 39 38 42 
4 3 40 41 42 
4 36 37 39 40 
4 36 37 38 41 
4 36 39 38 42 
4 36 40 41 42 
4 36 37 38 39 
4 36 37 41 40 
4 36 38 42 41 
4 36 39 42 40 

CELL_TYPES 60
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10

-------------- next part --------------
# vtk DataFile Version 3.0
vtk output
ASCII
DATASET UNSTRUCTURED_GRID
POINTS 31 float
0 0 0 1 0 1 1 1 0 
0 1 1 0 0 1 1 1 1 
1 0 0 0 1 0 0.5 0.5 0.5 
0.5 0 0.5 1 0.5 0.5 0.5 0.5 0 
0 0.5 0.5 0.5 0.5 1 0.5 1 0.5 
0.25 0.25 0.75 0.5 0 1 0 0 0.5 
0 0.5 1 0.75 0.75 0.75 0.5 1 1 
1 0.5 1 1 1 0.5 0.75 0.25 0.25 
0.5 0 0 1 0.5 0 1 0 0.5 
0.25 0.75 0.25 0.5 1 0 0 0.5 0 
0 1 0.5 
CELLS 60 300
4 0 9 11 12 
4 1 9 10 13 
4 2 11 10 14 
4 3 12 13 14 
4 8 9 11 12 
4 8 9 10 13 
4 8 11 10 14 
4 8 12 13 14 
4 8 9 10 11 
4 8 9 13 12 
4 8 10 14 13 
4 8 11 14 12 
4 1 9 13 16 
4 0 9 12 17 
4 3 13 12 18 
4 4 16 17 18 
4 15 9 13 16 
4 15 9 12 17 
4 15 13 12 18 
4 15 16 17 18 
4 15 9 12 13 
4 15 9 17 16 
4 15 12 18 17 
4 15 13 18 16 
4 1 13 21 10 
4 3 13 20 14 
4 5 21 20 22 
4 2 10 14 22 
4 19 13 21 10 
4 19 13 20 14 
4 19 21 20 22 
4 19 10 14 22 
4 19 13 20 21 
4 19 13 14 10 
4 19 20 22 14 
4 19 21 22 10 
4 0 24 11 9 
4 6 24 25 26 
4 2 11 25 10 
4 1 9 26 10 
4 23 24 11 9 
4 23 24 25 26 
4 23 11 25 10 
4 23 9 26 10 
4 23 24 25 11 
4 23 24 26 9 
4 23 25 10 26 
4 23 11 10 9 
4 0 11 29 12 
4 2 11 28 14 
4 7 29 28 30 
4 3 12 14 30 
4 27 11 29 12 
4 27 11 28 14 
4 27 29 28 30 
4 27 12 14 30 
4 27 11 28 29 
4 27 11 14 12 
4 27 28 30 14 
4 27 29 30 12 

CELL_TYPES 60
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10
10

-------------- next part --------------
/*=========================================================================

  Program:   Visualization Toolkit
  Module:    $RCSfile: vtkSubdivideTetra.cxx,v $
  Language:  C++
  Date:      $Date: 2000/12/10 20:08:26 $
  Version:   $Revision: 1.14 $


Copyright (c) 1993-2001 Ken Martin, Will Schroeder, Bill Lorensen 
All rights reserved.

Redistribution and use in source and binary forms, with or without
modification, are permitted provided that the following conditions are met:

 * Redistributions of source code must retain the above copyright notice,
   this list of conditions and the following disclaimer.

 * Redistributions in binary form must reproduce the above copyright notice,
   this list of conditions and the following disclaimer in the documentation
   and/or other materials provided with the distribution.

 * Neither name of Ken Martin, Will Schroeder, or Bill Lorensen nor the names
   of any contributors may be used to endorse or promote products derived
   from this software without specific prior written permission.

 * Modified source versions must be plainly marked as such, and must not be
   misrepresented as being the original software.

THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS ``AS IS''
AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
ARE DISCLAIMED. IN NO EVENT SHALL THE AUTHORS OR CONTRIBUTORS BE LIABLE FOR
ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL
DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR
SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY,
OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.

=========================================================================*/
#include "vtkSubdivideTetraClean.h"
#include "vtkMergePoints.h"
#include "vtkObjectFactory.h"



//------------------------------------------------------------------------------
vtkSubdivideTetraClean* vtkSubdivideTetraClean::New()
{
  // First try to create the object from the vtkObjectFactory
  vtkObject* ret = vtkObjectFactory::CreateInstance("vtkSubdivideTetraClean");
  if(ret)
    {
    return (vtkSubdivideTetraClean*)ret;
    }
  // If the factory was unable to create the object, then create it here.
  return new vtkSubdivideTetraClean;
}




// Description:
// Construct with all types of clipping turned off.
vtkSubdivideTetraClean::vtkSubdivideTetraClean()
{
}

void vtkSubdivideTetraClean::Execute()
{
  vtkUnstructuredGrid *input=(vtkUnstructuredGrid *)this->GetInput();
  int numPts = input->GetNumberOfPoints();
  int numCells = input->GetNumberOfCells();
  vtkPoints *inPts=input->GetPoints();
  int cellId, i, pts[4];
  vtkCell *cell;
  vtkPointData *pd = input->GetPointData();
  vtkUnstructuredGrid *output = this->GetOutput();
  vtkPointData *outputPD = output->GetPointData();
  vtkPoints *newPts;
  int ptId;
  vtkCellTypes *types=vtkCellTypes::New();
  float weights[4], x0[3], x1[3], x2[3], x3[3], x[3];
  int p0, p1, p2, p3, center;
  int e01, e02, e03, e12, e13, e23;
  vtkMergePoints *locator;
  
  vtkDebugMacro(<<"Executing mesh subdivide");

  input->GetCellTypes(types);
  if ( types->GetNumberOfTypes() != 1 || 
  types->GetCellType(0) != VTK_TETRA )
    {
    vtkErrorMacro(<<"Must be tetrahedra");
    types->Delete();
    return;
    }
  input->GetCells();
  
  // Copy original points and point data
  newPts = vtkPoints::New();
  newPts->Allocate(5*numPts,numPts);
  outputPD->InterpolateAllocate(pd,5*numPts,numPts);
  
  output->Allocate(numCells);
  output->SetPoints(newPts);

  locator = vtkMergePoints::New();
  locator->InitPointInsertion (newPts, input->GetBounds());

  for (ptId=0; ptId < numPts; ptId++)
    {
    locator->InsertNextPoint(inPts->GetPoint(ptId));
    outputPD->CopyData(pd,ptId,ptId);
    }

  // loop over tetrahedra, generating sixteen new ones for each. This is
  // done by introducing mid-edge nodes and a single mid-tetra node.
  for(cellId=0; cellId < numCells; cellId++)
    {
    cell = input->GetCell(cellId);

    // get tetra points
    cell->Points->GetPoint(0,x0);
    cell->Points->GetPoint(1,x1);
    cell->Points->GetPoint(2,x2);
    cell->Points->GetPoint(3,x3);

    p0 = cell->PointIds->GetId(0);
    p1 = cell->PointIds->GetId(1);
    p2 = cell->PointIds->GetId(2);
    p3 = cell->PointIds->GetId(3);

        // compute center point
    weights[0] = weights[1] = weights[2] = weights[3] = 0.25;
    for (i=0; i<3; i++)
      {
      x[i] = 0.25*(x0[i] + x1[i] + x2[i] + x3[i]);
      }
    locator->InsertUniquePoint(x,center );
    outputPD->InterpolatePoint(pd, center, cell->PointIds, weights);
    
    // compute edge points
    // edge 0-1
    for (i=0; i<3; i++)
      {
      x[i] = 0.5 * (x1[i] + x0[i]);
      }
    locator->InsertUniquePoint(x,e01);
    outputPD->InterpolateEdge(pd, e01, p0, p1, 0.5);
    
    // edge 1-2
    for (i=0; i<3; i++)
      {
      x[i] = 0.5 * (x2[i] + x1[i]);
      }
    locator->InsertUniquePoint(x,e12);
    outputPD->InterpolateEdge(pd, e12, p1, p2, 0.5);
    
    // edge 2-0
    for (i=0; i<3; i++)
      {
      x[i] = 0.5 * (x2[i] + x0[i]);
      }
    locator->InsertUniquePoint(x,e02);
    outputPD->InterpolateEdge(pd, e02, p2, p0, 0.5);
    
    // edge 0-3
    for (i=0; i<3; i++)
      {
      x[i] = 0.5 * (x3[i] + x0[i]);
      }
    locator->InsertUniquePoint(x,e03);
    outputPD->InterpolateEdge(pd, e03, p0, p3, 0.5);
    
    // edge 1-3
    for (i=0; i<3; i++)
      {
      x[i] = 0.5 * (x3[i] + x1[i]);
      }
    locator->InsertUniquePoint(x,e13);
    outputPD->InterpolateEdge(pd, e13, p1, p3, 0.5);
    
    // edge 2-3
    for (i=0; i<3; i++)
      {
      x[i] = 0.5 * (x3[i] + x2[i]);
      }
    locator->InsertUniquePoint(x,e23);
    outputPD->InterpolateEdge(pd, e23, p2, p3, 0.5);

    // Now create tetrahedra
    // First, four tetra from each vertex
    pts[0] = p0;
    pts[1] = e01;
    pts[2] = e02;
    pts[3] = e03;
    output->InsertNextCell(VTK_TETRA, 4, pts);
    pts[0] = p1;
    pts[1] = e01;
    pts[2] = e12;
    pts[3] = e13;
    output->InsertNextCell(VTK_TETRA, 4, pts);
    pts[0] = p2;
    pts[1] = e02;
    pts[2] = e12;
    pts[3] = e23;
    output->InsertNextCell(VTK_TETRA, 4, pts);
    pts[0] = p3;
    pts[1] = e03;
    pts[2] = e13;
    pts[3] = e23;
    output->InsertNextCell(VTK_TETRA, 4, pts);

    // Now four tetra from cut-off tetra corners
    pts[0] = center;
    pts[1] = e01;
    pts[2] = e02;
    pts[3] = e03;
    output->InsertNextCell(VTK_TETRA, 4, pts);
    pts[1] = e01;
    pts[2] = e12;
    pts[3] = e13;
    output->InsertNextCell(VTK_TETRA, 4, pts);
    pts[1] = e02;
    pts[2] = e12;
    pts[3] = e23;
    output->InsertNextCell(VTK_TETRA, 4, pts);
    pts[1] = e03;
    pts[2] = e13;
    pts[3] = e23;
    output->InsertNextCell(VTK_TETRA, 4, pts);

    // Now four tetra from triangles on tetra faces
    pts[0] = center;
    pts[1] = e01;
    pts[2] = e12;
    pts[3] = e02;
    output->InsertNextCell(VTK_TETRA, 4, pts);
    pts[1] = e01;
    pts[2] = e13;
    pts[3] = e03;
    output->InsertNextCell(VTK_TETRA, 4, pts);
    pts[1] = e12;
    pts[2] = e23;
    pts[3] = e13;
    output->InsertNextCell(VTK_TETRA, 4, pts);
    pts[1] = e02;
    pts[2] = e23;
    pts[3] = e03;
    output->InsertNextCell(VTK_TETRA, 4, pts);

    } //for all cells

  vtkDebugMacro(<<"Subdivided " << numCells << " cells");

  locator->Delete();
  types->Delete();
  newPts->Delete();
  output->Squeeze();
}

void vtkSubdivideTetraClean::PrintSelf(ostream& os, vtkIndent indent)
{
  vtkUnstructuredGridToUnstructuredGridFilter::PrintSelf(os,indent);
}



More information about the vtkusers mailing list