vtkMoleculeMapper.h 9.95 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
/*=========================================================================

  Program:   Visualization Toolkit
  Module:    vtkMoleculeMapper.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.

=========================================================================*/
15
16
17
18
19
20
21
/**
 * @class   vtkMoleculeMapper
 * @brief   Mapper that draws vtkMolecule objects
 *
 *
 * vtkMoleculeMapper uses glyphs (display lists) to quickly render a
 * molecule.
Nicolas Vuaille's avatar
Nicolas Vuaille committed
22
 */
23

24
25
#ifndef vtkMoleculeMapper_h
#define vtkMoleculeMapper_h
26

27
#include "vtkDomainsChemistryModule.h" // For export macro
28
29
30
31
32
33
34
35
36
#include "vtkMapper.h"
#include "vtkNew.h" // For vtkNew

class vtkActor;
class vtkGlyph3DMapper;
class vtkIdTypeArray;
class vtkMolecule;
class vtkPeriodicTable;
class vtkPolyData;
37
class vtkPolyDataMapper;
38
39
40
41
42
class vtkRenderer;
class vtkSelection;
class vtkSphereSource;
class vtkTrivialProducer;

43
class VTKDOMAINSCHEMISTRY_EXPORT vtkMoleculeMapper : public vtkMapper
44
45
{
public:
Nicolas Vuaille's avatar
Nicolas Vuaille committed
46
47
  static vtkMoleculeMapper* New();
  vtkTypeMacro(vtkMoleculeMapper, vtkMapper);
48
  void PrintSelf(ostream& os, vtkIndent indent) override;
49

50
51
52
53
  //@{
  /**
   * Get/Set the input vtkMolecule.
   */
Nicolas Vuaille's avatar
Nicolas Vuaille committed
54
55
  void SetInputData(vtkMolecule* in);
  vtkMolecule* GetInput();
56
57
58
59
60
61
62
63
64
65
66
67
68
  //@}

  /**
   * Set ivars to default ball-and-stick settings. This is equivalent
   * to the following:
   * - SetRenderAtoms( true )
   * - SetRenderBonds( true )
   * - SetAtomicRadiusType( VDWRadius )
   * - SetAtomicRadiusScaleFactor( 0.3 )
   * - SetBondColorMode( DiscreteByAtom )
   * - SetUseMultiCylindersForBonds( true )
   * - SetBondRadius( 0.075 )
   */
69
70
  void UseBallAndStickSettings();

71
72
73
74
75
76
77
78
79
80
81
  /**
   * Set ivars to default van der Waals spheres settings. This is
   * equivalent to the following:
   * - SetRenderAtoms( true )
   * - SetRenderBonds( true )
   * - SetAtomicRadiusType( VDWRadius )
   * - SetAtomicRadiusScaleFactor( 1.0 )
   * - SetBondColorMode( DiscreteByAtom )
   * - SetUseMultiCylindersForBonds( true )
   * - SetBondRadius( 0.075 )
   */
82
83
  void UseVDWSpheresSettings();

84
85
86
87
88
89
90
91
92
93
94
  /**
   * Set ivars to default liquorice stick settings. This is
   * equivalent to the following:
   * - SetRenderAtoms( true )
   * - SetRenderBonds( true )
   * - SetAtomicRadiusType( UnitRadius )
   * - SetAtomicRadiusScaleFactor( 0.1 )
   * - SetBondColorMode( DiscreteByAtom )
   * - SetUseMultiCylindersForBonds( false )
   * - SetBondRadius( 0.1 )
   */
95
  void UseLiquoriceStickSettings();
96

97
98
99
100
101
102
103
104
105
106
107
108
109
110
  /**
   * Set ivars to use fast settings that may be useful for rendering
   * extremely large molecules where the overall shape is more
   * important than the details of the atoms/bond. This is equivalent
   * to the following:
   * - SetRenderAtoms( true )
   * - SetRenderBonds( true )
   * - SetAtomicRadiusType( UnitRadius )
   * - SetAtomicRadiusScaleFactor( 0.60 )
   * - SetBondColorMode( SingleColor )
   * - SetBondColor( 50, 50, 50 )
   * - SetUseMultiCylindersForBonds( false )
   * - SetBondRadius( 0.075 )
   */
111
112
  void UseFastSettings();

113
114
115
116
  //@{
  /**
   * Get/Set whether or not to render atoms. Default: On.
   */
117
118
119
  vtkGetMacro(RenderAtoms, bool);
  vtkSetMacro(RenderAtoms, bool);
  vtkBooleanMacro(RenderAtoms, bool);
120
  //@}
121

122
123
124
125
  //@{
  /**
   * Get/Set whether or not to render bonds. Default: On.
   */
126
127
128
  vtkGetMacro(RenderBonds, bool);
  vtkSetMacro(RenderBonds, bool);
  vtkBooleanMacro(RenderBonds, bool);
129
  //@}
130

131
132
133
134
135
  //@{
  /**
   * Get/Set whether or not to render the unit cell lattice, if present.
   * Default: On.
   */
Nicolas Vuaille's avatar
Nicolas Vuaille committed
136
137
138
  vtkGetMacro(RenderLattice, bool);
  vtkSetMacro(RenderLattice, bool);
  vtkBooleanMacro(RenderLattice, bool);
139
  //@}
140

Nicolas Vuaille's avatar
Nicolas Vuaille committed
141
142
  enum
  {
143
144
    CovalentRadius = 0,
    VDWRadius,
145
146
    UnitRadius,
    CustomArrayRadius
147
148
  };

149
150
151
152
153
154
  //@{
  /**
   * Get/Set the type of radius used to generate the atoms. Default:
   * VDWRadius. If CustomArrayRadius is used, the VertexData array named
   * 'radii' is used for per-atom radii.
   */
155
156
  vtkGetMacro(AtomicRadiusType, int);
  vtkSetMacro(AtomicRadiusType, int);
Nicolas Vuaille's avatar
Nicolas Vuaille committed
157
  const char* GetAtomicRadiusTypeAsString();
158
159
160
161
  void SetAtomicRadiusTypeToCovalentRadius() { this->SetAtomicRadiusType(CovalentRadius); }
  void SetAtomicRadiusTypeToVDWRadius() { this->SetAtomicRadiusType(VDWRadius); }
  void SetAtomicRadiusTypeToUnitRadius() { this->SetAtomicRadiusType(UnitRadius); }
  void SetAtomicRadiusTypeToCustomArrayRadius() { this->SetAtomicRadiusType(CustomArrayRadius); }
162
163
164
165
166
167
168
169
  //@}

  //@{
  /**
   * Get/Set the uniform scaling factor applied to the atoms.
   * This is ignored when AtomicRadiusType == CustomArrayRadius.
   * Default: 0.3.
   */
170
171
  vtkGetMacro(AtomicRadiusScaleFactor, float);
  vtkSetMacro(AtomicRadiusScaleFactor, float);
172
  //@}
173

174
175
176
177
178
  //@{
  /**
   * Get/Set whether multicylinders will be used to represent multiple
   * bonds. Default: On.
   */
179
180
181
  vtkGetMacro(UseMultiCylindersForBonds, bool);
  vtkSetMacro(UseMultiCylindersForBonds, bool);
  vtkBooleanMacro(UseMultiCylindersForBonds, bool);
182
  //@}
183

Nicolas Vuaille's avatar
Nicolas Vuaille committed
184
185
  enum
  {
186
    SingleColor = 0,
187
    DiscreteByAtom
188
189
  };

190
191
192
193
194
195
196
197
198
199
200
  //@{
  /**
   * Get/Set the method by which bonds are colored.

   * If 'SingleColor' is used, all bonds will be the same color. Use
   * SetBondColor to set the rgb values used.

   * If 'DiscreteByAtom' is selected, each bond is colored using the
   * same lookup table as the atoms at each end, with a sharp color
   * boundary at the bond center.
   */
201
  vtkGetMacro(BondColorMode, int);
Nicolas Vuaille's avatar
Nicolas Vuaille committed
202
  vtkSetClampMacro(BondColorMode, int, SingleColor, DiscreteByAtom);
203
204
  void SetBondColorModeToSingleColor() { this->SetBondColorMode(SingleColor); }
  void SetBondColorModeToDiscreteByAtom() { this->SetBondColorMode(DiscreteByAtom); }
Nicolas Vuaille's avatar
Nicolas Vuaille committed
205
  const char* GetBondColorModeAsString();
206
  //@}
207

Nicolas Vuaille's avatar
Nicolas Vuaille committed
208
209
210
211
212
213
214
215
216
217
218
  //@{
  /**
   * Get/Set the method by which atoms are colored.
   *
   * If 'SingleColor' is used, all atoms will have the same color. Use
   * SetAtomColor to set the rgb values to be used.
   *
   * If 'DiscreteByAtom' is selected, each atom is colored using the
   * internal lookup table.
   */
  vtkGetMacro(AtomColorMode, int);
219
  vtkSetClampMacro(AtomColorMode, int, SingleColor, DiscreteByAtom);
Nicolas Vuaille's avatar
Nicolas Vuaille committed
220
221
222
223
224
225
226
227
228
229
230
  //@}

  //@{
  /**
   * Get/Set the color of the atoms as an rgb tuple.
   * Default: {150, 150, 150} (grey)
   */
  vtkGetVector3Macro(AtomColor, unsigned char);
  vtkSetVector3Macro(AtomColor, unsigned char);
  //@}

231
232
233
234
235
  //@{
  /**
   * Get/Set the color of the bonds as an rgb tuple.
   * Default: {50, 50, 50} (dark grey)
   */
236
237
  vtkGetVector3Macro(BondColor, unsigned char);
  vtkSetVector3Macro(BondColor, unsigned char);
238
  //@}
239

240
241
242
243
  //@{
  /**
   * Get/Set the radius of the bond cylinders. Default: 0.075
   */
244
245
  vtkGetMacro(BondRadius, float);
  vtkSetMacro(BondRadius, float);
246
  //@}
247

248
249
250
251
252
  //@{
  /**
   * Get/Set the color of the bonds as an rgb tuple.
   * Default: {255, 255, 255} (white)
   */
Nicolas Vuaille's avatar
Nicolas Vuaille committed
253
254
  vtkGetVector3Macro(LatticeColor, unsigned char);
  vtkSetVector3Macro(LatticeColor, unsigned char);
255
  //@}
256

257
258
259
260
261
  //@{
  /**
   * Extract the ids atoms and/or bonds rendered by this molecule from a
   * vtkSelection object. The vtkIdTypeArray
   */
Nicolas Vuaille's avatar
Nicolas Vuaille committed
262
263
264
  virtual void GetSelectedAtomsAndBonds(
    vtkSelection* selection, vtkIdTypeArray* atomIds, vtkIdTypeArray* bondIds);
  virtual void GetSelectedAtoms(vtkSelection* selection, vtkIdTypeArray* atomIds)
265
  {
266
    this->GetSelectedAtomsAndBonds(selection, atomIds, nullptr);
267
  }
Nicolas Vuaille's avatar
Nicolas Vuaille committed
268
  virtual void GetSelectedBonds(vtkSelection* selection, vtkIdTypeArray* bondIds)
269
  {
270
    this->GetSelectedAtomsAndBonds(selection, nullptr, bondIds);
271
  }
272
  //@}
273

274
275
276
277
  //@{
  /**
   * Reimplemented from base class
   */
Nicolas Vuaille's avatar
Nicolas Vuaille committed
278
279
280
  void Render(vtkRenderer*, vtkActor*) override;
  void ReleaseGraphicsResources(vtkWindow*) override;
  double* GetBounds() override;
281
282
  void GetBounds(double bounds[6]) override { vtkAbstractMapper3D::GetBounds(bounds); }
  int FillInputPortInformation(int port, vtkInformation* info) override;
Nicolas Vuaille's avatar
Nicolas Vuaille committed
283
  bool GetSupportsSelection() override { return true; }
284
  //@}
285

Nicolas Vuaille's avatar
Nicolas Vuaille committed
286
  //@{
287
288
289
290
291
292
  /**
   * Get/Set the atomic radius array name. Default: "radii"
   * It is only used when AtomicRadiusType is set to CustomArrayRadius.
   */
  vtkGetStringMacro(AtomicRadiusArrayName);
  vtkSetStringMacro(AtomicRadiusArrayName);
Nicolas Vuaille's avatar
Nicolas Vuaille committed
293
  //@}
294

Nicolas Vuaille's avatar
Nicolas Vuaille committed
295
296
297
298
299
300
  /**
   * Helper method to set ScalarMode on both AtomGlyphMapper and BondGlyphMapper.
   * true means VTK_COLOR_MODE_MAP_SCALARS, false VTK_COLOR_MODE_DIRECT_SCALARS.
   */
  virtual void SetMapScalars(bool map);

David E. DeMarle's avatar
David E. DeMarle committed
301
302
  vtkPeriodicTable* GetPeriodicTable() { return this->PeriodicTable; }

303
304
protected:
  vtkMoleculeMapper();
305
  ~vtkMoleculeMapper() override;
306

307
308
309
310
  //@{
  /**
   * Customize atom rendering
   */
311
312
313
  bool RenderAtoms;
  int AtomicRadiusType;
  float AtomicRadiusScaleFactor;
314
  char* AtomicRadiusArrayName;
Nicolas Vuaille's avatar
Nicolas Vuaille committed
315
316
  int AtomColorMode;
  unsigned char AtomColor[3];
317
  //@}
318

319
320
321
322
  //@{
  /**
   * Customize bond rendering
   */
323
324
325
326
327
  bool RenderBonds;
  int BondColorMode;
  bool UseMultiCylindersForBonds;
  float BondRadius;
  unsigned char BondColor[3];
328
  //@}
329

330
331
  bool RenderLattice;

332
333
334
  /**
   * Internal render methods
   */
Nicolas Vuaille's avatar
Nicolas Vuaille committed
335
  void GlyphRender(vtkRenderer* ren, vtkActor* act);
336

337
338
339
340
  //@{
  /**
   * Cached variables and update methods
   */
341
342
343
344
345
  vtkNew<vtkPolyData> AtomGlyphPolyData;
  vtkNew<vtkTrivialProducer> AtomGlyphPointOutput;
  vtkNew<vtkPolyData> BondGlyphPolyData;
  vtkNew<vtkTrivialProducer> BondGlyphPointOutput;
  bool GlyphDataInitialized;
346
347
348
  virtual void UpdateGlyphPolyData();
  virtual void UpdateAtomGlyphPolyData();
  virtual void UpdateBondGlyphPolyData();
349
  //@}
350

351
352
353
354
  //@{
  /**
   * Internal mappers
   */
355
356
  vtkNew<vtkGlyph3DMapper> AtomGlyphMapper;
  vtkNew<vtkGlyph3DMapper> BondGlyphMapper;
357
  //@}
358

359
360
361
362
363
  unsigned char LatticeColor[3];
  vtkNew<vtkPolyData> LatticePolyData;
  vtkNew<vtkPolyDataMapper> LatticeMapper;
  virtual void UpdateLatticePolyData();

364
365
366
  /**
   * Periodic table for lookups
   */
367
368
369
  vtkNew<vtkPeriodicTable> PeriodicTable;

private:
370
371
  vtkMoleculeMapper(const vtkMoleculeMapper&) = delete;
  void operator=(const vtkMoleculeMapper&) = delete;
372
373
374
};

#endif