VTK  9.3.1
vtkMoleculeMapper.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
15 #ifndef vtkMoleculeMapper_h
16 #define vtkMoleculeMapper_h
17 
18 #include "vtkDomainsChemistryModule.h" // For export macro
19 #include "vtkMapper.h"
20 #include "vtkNew.h" // For vtkNew
21 
22 VTK_ABI_NAMESPACE_BEGIN
23 class vtkActor;
24 class vtkGlyph3DMapper;
25 class vtkIdTypeArray;
26 class vtkMolecule;
27 class vtkPeriodicTable;
28 class vtkPolyData;
29 class vtkPolyDataMapper;
30 class vtkRenderer;
31 class vtkSelection;
32 class vtkSphereSource;
33 class vtkTrivialProducer;
34 
35 class VTKDOMAINSCHEMISTRY_EXPORT vtkMoleculeMapper : public vtkMapper
36 {
37 public:
38  static vtkMoleculeMapper* New();
39  vtkTypeMacro(vtkMoleculeMapper, vtkMapper);
40  void PrintSelf(ostream& os, vtkIndent indent) override;
41 
43 
46  void SetInputData(vtkMolecule* in);
49 
61  void UseBallAndStickSettings();
62 
74  void UseVDWSpheresSettings();
75 
87  void UseLiquoriceStickSettings();
88 
103  void UseFastSettings();
104 
106 
109  vtkGetMacro(RenderAtoms, bool);
110  vtkSetMacro(RenderAtoms, bool);
111  vtkBooleanMacro(RenderAtoms, bool);
113 
115 
118  vtkGetMacro(RenderBonds, bool);
119  vtkSetMacro(RenderBonds, bool);
120  vtkBooleanMacro(RenderBonds, bool);
122 
124 
128  vtkGetMacro(RenderLattice, bool);
129  vtkSetMacro(RenderLattice, bool);
130  vtkBooleanMacro(RenderLattice, bool);
132 
133  enum
134  {
135  CovalentRadius = 0,
138  CustomArrayRadius
139  };
140 
142 
147  vtkGetMacro(AtomicRadiusType, int);
148  vtkSetMacro(AtomicRadiusType, int);
149  const char* GetAtomicRadiusTypeAsString();
150  void SetAtomicRadiusTypeToCovalentRadius() { this->SetAtomicRadiusType(CovalentRadius); }
151  void SetAtomicRadiusTypeToVDWRadius() { this->SetAtomicRadiusType(VDWRadius); }
152  void SetAtomicRadiusTypeToUnitRadius() { this->SetAtomicRadiusType(UnitRadius); }
153  void SetAtomicRadiusTypeToCustomArrayRadius() { this->SetAtomicRadiusType(CustomArrayRadius); }
155 
157 
162  vtkGetMacro(AtomicRadiusScaleFactor, float);
163  vtkSetMacro(AtomicRadiusScaleFactor, float);
165 
167 
171  vtkGetMacro(UseMultiCylindersForBonds, bool);
172  vtkSetMacro(UseMultiCylindersForBonds, bool);
173  vtkBooleanMacro(UseMultiCylindersForBonds, bool);
175 
176  enum
177  {
178  SingleColor = 0,
179  DiscreteByAtom
180  };
181 
183 
193  vtkGetMacro(BondColorMode, int);
194  vtkSetClampMacro(BondColorMode, int, SingleColor, DiscreteByAtom);
195  void SetBondColorModeToSingleColor() { this->SetBondColorMode(SingleColor); }
196  void SetBondColorModeToDiscreteByAtom() { this->SetBondColorMode(DiscreteByAtom); }
197  const char* GetBondColorModeAsString();
199 
201 
210  vtkGetMacro(AtomColorMode, int);
211  vtkSetClampMacro(AtomColorMode, int, SingleColor, DiscreteByAtom);
213 
215 
219  vtkGetVector3Macro(AtomColor, unsigned char);
220  vtkSetVector3Macro(AtomColor, unsigned char);
222 
224 
228  vtkGetVector3Macro(BondColor, unsigned char);
229  vtkSetVector3Macro(BondColor, unsigned char);
231 
233 
236  vtkGetMacro(BondRadius, float);
237  vtkSetMacro(BondRadius, float);
239 
241 
245  vtkGetVector3Macro(LatticeColor, unsigned char);
246  vtkSetVector3Macro(LatticeColor, unsigned char);
248 
250 
254  virtual void GetSelectedAtomsAndBonds(
255  vtkSelection* selection, vtkIdTypeArray* atomIds, vtkIdTypeArray* bondIds);
256  virtual void GetSelectedAtoms(vtkSelection* selection, vtkIdTypeArray* atomIds)
257  {
258  this->GetSelectedAtomsAndBonds(selection, atomIds, nullptr);
259  }
260  virtual void GetSelectedBonds(vtkSelection* selection, vtkIdTypeArray* bondIds)
261  {
262  this->GetSelectedAtomsAndBonds(selection, nullptr, bondIds);
263  }
265 
267 
270  void Render(vtkRenderer*, vtkActor*) override;
271  void ReleaseGraphicsResources(vtkWindow*) override;
272  double* GetBounds() override;
273  void GetBounds(double bounds[6]) override { Superclass::GetBounds(bounds); }
274  int FillInputPortInformation(int port, vtkInformation* info) override;
275  bool GetSupportsSelection() override { return true; }
277 
279 
283  vtkGetStringMacro(AtomicRadiusArrayName);
284  vtkSetStringMacro(AtomicRadiusArrayName);
286 
291  virtual void SetMapScalars(bool map);
292 
296  vtkPeriodicTable* GetPeriodicTable() { return this->PeriodicTable; }
297 
298 protected:
300  ~vtkMoleculeMapper() override;
301 
303 
311  unsigned char AtomColor[3];
313 
315 
321  float BondRadius;
322  unsigned char BondColor[3];
324 
326 
330  void GlyphRender(vtkRenderer* ren, vtkActor* act);
331 
333 
341  virtual void UpdateGlyphPolyData();
342  virtual void UpdateAtomGlyphPolyData();
343  virtual void UpdateBondGlyphPolyData();
345 
347 
353 
354  unsigned char LatticeColor[3];
357  virtual void UpdateLatticePolyData();
358 
363 
364 private:
365  vtkMoleculeMapper(const vtkMoleculeMapper&) = delete;
366  void operator=(const vtkMoleculeMapper&) = delete;
367 };
368 
369 VTK_ABI_NAMESPACE_END
370 #endif
Access to information about the elements.
represents an object (geometry & properties) in a rendered scene
Definition: vtkActor.h:40
vtkNew< vtkGlyph3DMapper > BondGlyphMapper
Internal mappers.
bool RenderBonds
Customize bond rendering.
vtkNew< vtkPeriodicTable > PeriodicTable
Periodic table for lookups.
bool GetSupportsSelection() override
Reimplemented from base class.
Store vtkAlgorithm input/output information.
class describing a molecule
Definition: vtkMolecule.h:83
char * AtomicRadiusArrayName
Customize atom rendering.
float AtomicRadiusScaleFactor
Customize atom rendering.
bool GlyphDataInitialized
Cached variables and update methods.
void SetAtomicRadiusTypeToCovalentRadius()
Get/Set the type of radius used to generate the atoms.
bool UseMultiCylindersForBonds
Customize bond rendering.
virtual void GetSelectedBonds(vtkSelection *selection, vtkIdTypeArray *bondIds)
Extract the ids atoms and/or bonds rendered by this molecule from a vtkSelection object.
vtkNew< vtkPolyData > LatticePolyData
abstract specification for renderers
Definition: vtkRenderer.h:61
vtkNew< vtkPolyData > AtomGlyphPolyData
Cached variables and update methods.
data object that represents a "selection" in VTK.
Definition: vtkSelection.h:49
dynamic, self-adjusting array of vtkIdType
bool RenderAtoms
Customize atom rendering.
vtkGlyph3D on the GPU.
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:79
void SetAtomicRadiusTypeToVDWRadius()
Get/Set the type of radius used to generate the atoms.
void SetAtomicRadiusTypeToCustomArrayRadius()
Get/Set the type of radius used to generate the atoms.
int AtomColorMode
Customize atom rendering.
Mapper that draws vtkMolecule objects.
int AtomicRadiusType
Customize atom rendering.
window superclass for vtkRenderWindow
Definition: vtkWindow.h:27
create a polygonal sphere centered at the origin
vtkNew< vtkPolyData > BondGlyphPolyData
Cached variables and update methods.
Producer for stand-alone data objects.
vtkNew< vtkPolyDataMapper > LatticeMapper
a simple class to control print indentation
Definition: vtkIndent.h:28
vtkNew< vtkTrivialProducer > AtomGlyphPointOutput
Cached variables and update methods.
vtkNew< vtkGlyph3DMapper > AtomGlyphMapper
Internal mappers.
void SetAtomicRadiusTypeToUnitRadius()
Get/Set the type of radius used to generate the atoms.
virtual double * GetBounds()=0
Return bounding box (array of six doubles) of data expressed as (xmin,xmax, ymin,ymax, zmin,zmax).
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
virtual int FillInputPortInformation(int port, vtkInformation *info)
Fill the input port information objects for this algorithm.
vtkNew< vtkTrivialProducer > BondGlyphPointOutput
Cached variables and update methods.
abstract class specifies interface to map data to graphics primitives
Definition: vtkMapper.h:76
vtkPeriodicTable * GetPeriodicTable()
Accessor to internal structure.
map vtkPolyData to graphics primitives
void SetBondColorModeToDiscreteByAtom()
Get/Set the method by which bonds are colored.
float BondRadius
Customize bond rendering.
virtual void GetSelectedAtoms(vtkSelection *selection, vtkIdTypeArray *atomIds)
Extract the ids atoms and/or bonds rendered by this molecule from a vtkSelection object.
static vtkAlgorithm * New()
int BondColorMode
Customize bond rendering.
double * GetBounds() override
Return bounding box (array of six doubles) of data expressed as (xmin,xmax, ymin,ymax, zmin,zmax).
vtkDataSet * GetInput()
Get the input as a vtkDataSet.
virtual void Render(vtkRenderer *ren, vtkActor *a)=0
Method initiates the mapping process.
void SetBondColorModeToSingleColor()
Get/Set the method by which bonds are colored.
void ReleaseGraphicsResources(vtkWindow *) override
Release any graphics resources that are being consumed by this mapper.
Definition: vtkMapper.h:104
void GetBounds(double bounds[6]) override
Reimplemented from base class.