Skip to content
Snippets Groups Projects
Commit 521898a2 authored by Nicolas Vuaille's avatar Nicolas Vuaille
Browse files

Molecule to polydata conversion

 * introduce vtkMoleculeToLinesFilter
 ** Create polydata with atom == point, bond == line cell. Copy data.

 * Fix AtomBall and BondStick conversion:
 ** set array name for bond orders and atomic numbers
 ** fix allocated size
parent 3ce232d9
No related merge requests found
......@@ -7,6 +7,7 @@ set(Module_SRCS
vtkMoleculeMapper.cxx
vtkMoleculeToAtomBallFilter.cxx
vtkMoleculeToBondStickFilter.cxx
vtkMoleculeToLinesFilter.cxx
vtkMoleculeToPolyDataFilter.cxx
vtkPeriodicTable.cxx
vtkPointSetToMoleculeFilter.cxx
......
......@@ -28,6 +28,7 @@ vtk_add_test_cxx(vtkDomainsChemistryCxxTests tests
TestMoleculeIOLegacy.cxx
TestMoleculeSelection.cxx,NO_VALID
TestMoleculeMapperPropertyUpdate.cxx
TestMoleculeToLines.cxx,NO_VALID
TestMultiCylinderOn.cxx
TestMultiCylinderOff.cxx
TestPeriodicTable.cxx,NO_VALID
......
/*=========================================================================
Program: Visualization Toolkit
Module: TestMoleculeToLines.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 "vtkTestUtilities.h"
#include "vtkCMLMoleculeReader.h"
#include "vtkCellData.h"
#include "vtkDataSetAttributes.h"
#include "vtkMolecule.h"
#include "vtkMoleculeToLinesFilter.h"
#include "vtkNew.h"
#include "vtkPointData.h"
#include "vtkPolyData.h"
#define CheckNumbers(name, first, second) \
if (first != second) \
{ \
cerr << "Error : wrong number of " << #name << ". Got " << first << " but expects " << second \
<< endl; \
return EXIT_FAILURE; \
}
int TestMoleculeToLines(int argc, char* argv[])
{
char* fileName = vtkTestUtilities::ExpandDataFileName(argc, argv, "Data/porphyrin.cml");
// read molecule from cml file
vtkNew<vtkCMLMoleculeReader> reader;
reader->SetFileName(fileName);
reader->Update();
delete[] fileName;
vtkMolecule* molecule = reader->GetOutput();
// convert
vtkNew<vtkMoleculeToLinesFilter> converter;
converter->SetInputConnection(reader->GetOutputPort());
converter->Update();
vtkPolyData* poly = converter->GetOutput();
// check number of points, lines and associated data
CheckNumbers("points", poly->GetNumberOfPoints(), molecule->GetNumberOfAtoms());
CheckNumbers("lines", poly->GetNumberOfLines(), molecule->GetNumberOfBonds());
CheckNumbers("pointData",
poly->GetPointData()->GetNumberOfArrays(),
molecule->GetAtomData()->GetNumberOfArrays());
CheckNumbers("cellData",
poly->GetCellData()->GetNumberOfArrays(),
molecule->GetBondData()->GetNumberOfArrays());
return EXIT_SUCCESS;
}
......@@ -59,6 +59,7 @@ int vtkMoleculeToAtomBallFilter::RequestData(
vtkCellArray *polys = vtkCellArray::New();
vtkPoints *points = vtkPoints::New();
vtkUnsignedShortArray *atomicNums = vtkUnsignedShortArray::New();
atomicNums->SetName(input->GetAtomicNumberArrayName());
// Initialize a SphereSource
vtkSphereSource *sphereSource = vtkSphereSource::New();
......@@ -71,7 +72,8 @@ int vtkMoleculeToAtomBallFilter::RequestData(
GetNumberOfPoints());
polys->Allocate(numAtoms * sphereSource->GetOutput()->GetPolys()->
GetNumberOfCells());
atomicNums->Allocate(points->GetNumberOfPoints());
atomicNums->Allocate(numAtoms * sphereSource->GetOutput()->GetPoints()->
GetNumberOfPoints());
// Initialize some variables for later
vtkIdType numCellPoints, *cellPoints;
......
......@@ -52,6 +52,7 @@ int vtkMoleculeToBondStickFilter::RequestData(
vtkCellArray *polys = vtkCellArray::New();
vtkPoints *points = vtkPoints::New();
vtkUnsignedShortArray *bondOrders = vtkUnsignedShortArray::New();
bondOrders->SetName(input->GetBondOrdersArrayName());
// Initialize a CylinderSource
vtkCylinderSource *cylSource = vtkCylinderSource::New();
......@@ -64,7 +65,8 @@ int vtkMoleculeToBondStickFilter::RequestData(
GetNumberOfPoints());
polys->Allocate(3 * numBonds * cylSource->GetOutput()->GetPolys()->
GetNumberOfCells());
bondOrders->Allocate(points->GetNumberOfPoints());
bondOrders->Allocate(3 * numBonds * cylSource->GetOutput()->GetPoints()->
GetNumberOfPoints());
// Create a transform object to map the cylinder source to the bond
vtkTransform *xform = vtkTransform::New();
......
/*=========================================================================
Program: Visualization Toolkit
Module: vtkMoleculeToLinesFilter.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 "vtkMoleculeToLinesFilter.h"
#include "vtkCellArray.h"
#include "vtkCellData.h"
#include "vtkInformation.h"
#include "vtkMolecule.h"
#include "vtkPointData.h"
vtkStandardNewMacro(vtkMoleculeToLinesFilter);
//----------------------------------------------------------------------------
int vtkMoleculeToLinesFilter::RequestData(vtkInformation*,
vtkInformationVector** inputVector,
vtkInformationVector* outputVector)
{
vtkMolecule* input = vtkMolecule::SafeDownCast(vtkDataObject::GetData(inputVector[0]));
vtkPolyData* output = vtkPolyData::SafeDownCast(vtkDataObject::GetData(outputVector));
vtkNew<vtkCellArray> bonds;
// 2 point ids + 1 VTKCellType = 3 values per bonds
bonds->Allocate(3 * input->GetNumberOfBonds());
for (vtkIdType bondInd = 0; bondInd < input->GetNumberOfBonds(); ++bondInd)
{
vtkBond bond = input->GetBond(bondInd);
vtkIdType ids[2] = { bond.GetBeginAtomId(), bond.GetEndAtomId() };
bonds->InsertNextCell(2, ids);
}
output->SetPoints(input->GetAtomicPositionArray());
output->SetLines(bonds);
output->GetPointData()->DeepCopy(input->GetAtomData());
output->GetCellData()->DeepCopy(input->GetBondData());
return 1;
}
/*=========================================================================
Program: Visualization Toolkit
Module: vtkMoleculeToLinesFilter.h
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.
=========================================================================*/
/**
* @class vtkMoleculeToLinesFilter
* @brief Convert a molecule into a simple polydata with lines.
*
* vtkMoleculeToLinesFilter is a filter class that takes vtkMolecule as input and
* generates polydata on output.
* Conversion is done following this rules:
* - 1 atom == 1 point
* - 1 bond == 1 line (cell of type VTK_LINE)
* - atom data is copied as point data
* - bond data is copied as cell data
*/
#ifndef vtkMoleculeToLinesFilter_h
#define vtkMoleculeToLinesFilter_h
#include "vtkDomainsChemistryModule.h" // For export macro
#include "vtkMoleculeToPolyDataFilter.h"
class VTKDOMAINSCHEMISTRY_EXPORT vtkMoleculeToLinesFilter : public vtkMoleculeToPolyDataFilter
{
public:
static vtkMoleculeToLinesFilter* New();
vtkTypeMacro(vtkMoleculeToLinesFilter, vtkMoleculeToPolyDataFilter);
protected:
vtkMoleculeToLinesFilter() = default;
~vtkMoleculeToLinesFilter() = default;
int RequestData(vtkInformation*, vtkInformationVector**, vtkInformationVector*) override;
private:
vtkMoleculeToLinesFilter(const vtkMoleculeToLinesFilter&) = delete;
void operator=(const vtkMoleculeToLinesFilter&) = delete;
};
#endif
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment