source: vtkvis/trunk/Molecule.cpp

Last change on this file was 6409, checked in by ldelgass, 8 years ago

Set default const color on atoms, set default colormode to constant or scalar
if no element array.

  • Property svn:eol-style set to native
File size: 35.1 KB
Line 
1/* -*- mode: c++; c-basic-offset: 4; indent-tabs-mode: nil -*- */
2/*
3 * Copyright (C) 2004-2012  HUBzero Foundation, LLC
4 *
5 * Author: Leif Delgass <ldelgass@purdue.edu>
6 */
7
8#include <cstdio>
9#include <cfloat>
10#include <cassert>
11
12#include <vtkDataSet.h>
13#include <vtkCellArray.h>
14#include <vtkPointData.h>
15#include <vtkFloatArray.h>
16#include <vtkDoubleArray.h>
17#include <vtkIntArray.h>
18#include <vtkStringArray.h>
19#include <vtkPolyData.h>
20#include <vtkPolyDataMapper.h>
21#include <vtkActor.h>
22#include <vtkProperty.h>
23#include <vtkSphereSource.h>
24#include <vtkCylinderSource.h>
25#include <vtkTransform.h>
26#include <vtkTransformPolyDataFilter.h>
27#include <vtkGlyph3DMapper.h>
28#include <vtkPointSetToLabelHierarchy.h>
29#include <vtkLabelPlacementMapper.h>
30#include <vtkTextProperty.h>
31#include <vtkTransformPolyDataFilter.h>
32
33#include "Molecule.h"
34#include "MoleculeData.h"
35#include "Renderer.h"
36#include "Trace.h"
37
38using namespace VtkVis;
39
40Molecule::Molecule() :
41    GraphicsObject(),
42    _radiusScale(0.3),
43    _atomScaling(COVALENT_RADIUS),
44    _labelsOn(false),
45    _colorMap(NULL),
46    _colorMode(COLOR_BY_ELEMENTS),
47    _colorFieldType(DataSet::POINT_DATA)
48{
49    _bondColor[0] = _bondColor[1] = _bondColor[2] = 1.0f;
50    _colorFieldRange[0] = DBL_MAX;
51    _colorFieldRange[1] = -DBL_MAX;
52    _faceCulling = true;
53}
54
55Molecule::~Molecule()
56{
57#ifdef WANT_TRACE
58    if (_dataSet != NULL)
59        TRACE("Deleting Molecule for %s", _dataSet->getName().c_str());
60    else
61        TRACE("Deleting Molecule with NULL DataSet");
62#endif
63}
64
65void Molecule::setDataSet(DataSet *dataSet,
66                          Renderer *renderer)
67{
68    if (_dataSet != dataSet) {
69        _dataSet = dataSet;
70
71        _renderer = renderer;
72
73        if (renderer->getUseCumulativeRange()) {
74            renderer->getCumulativeDataRange(_dataRange,
75                                             _dataSet->getActiveScalarsName(),
76                                             1);
77            renderer->getCumulativeDataRange(_vectorMagnitudeRange,
78                                             _dataSet->getActiveVectorsName(),
79                                             3);
80            for (int i = 0; i < 3; i++) {
81                renderer->getCumulativeDataRange(_vectorComponentRange[i],
82                                                 _dataSet->getActiveVectorsName(),
83                                                 3, i);
84            }
85        } else {
86            _dataSet->getScalarRange(_dataRange);
87            _dataSet->getVectorRange(_vectorMagnitudeRange);
88            for (int i = 0; i < 3; i++) {
89                _dataSet->getVectorRange(_vectorComponentRange[i], i);
90            }
91        }
92
93        update();
94    }
95}
96
97/**
98 * \brief Create and initialize VTK Props to render a Molecule
99 */
100void Molecule::initProp()
101{
102    if (_atomProp == NULL) {
103        _atomProp = vtkSmartPointer<vtkActor>::New();
104        if (_faceCulling && _opacity == 1.0)
105            setCulling(_atomProp->GetProperty(), true);
106        _atomProp->GetProperty()->EdgeVisibilityOff();
107        _atomProp->GetProperty()->SetColor(_color[0], _color[1], _color[2]);
108        _atomProp->GetProperty()->SetEdgeColor(_edgeColor[0], _edgeColor[1], _edgeColor[2]);
109        _atomProp->GetProperty()->SetLineWidth(_edgeWidth);
110        _atomProp->GetProperty()->SetOpacity(_opacity);
111        _atomProp->GetProperty()->SetAmbient(.2);
112        _atomProp->GetProperty()->SetSpecular(.2);
113        _atomProp->GetProperty()->SetSpecularPower(80.0);
114        if (!_lighting)
115            _atomProp->GetProperty()->LightingOff();
116    }
117    if (_bondProp == NULL) {
118        _bondProp = vtkSmartPointer<vtkActor>::New();
119        if (_faceCulling && _opacity == 1.0)
120            setCulling(_bondProp->GetProperty(), true);
121        _bondProp->GetProperty()->EdgeVisibilityOff();
122        _bondProp->GetProperty()->SetColor(_bondColor[0], _bondColor[1], _bondColor[2]);
123        _bondProp->GetProperty()->SetEdgeColor(_edgeColor[0], _edgeColor[1], _edgeColor[2]);
124        _bondProp->GetProperty()->SetLineWidth(_edgeWidth);
125        _bondProp->GetProperty()->SetOpacity(_opacity);
126        _bondProp->GetProperty()->SetAmbient(.2);
127        _bondProp->GetProperty()->SetSpecular(.2);
128        _bondProp->GetProperty()->SetSpecularPower(80.0);
129        if (!_lighting)
130            _bondProp->GetProperty()->LightingOff();
131    }
132    if (_labelProp == NULL) {
133        _labelProp = vtkSmartPointer<vtkActor2D>::New();
134    }
135    if (_prop == NULL) {
136        _prop = vtkSmartPointer<vtkAssembly>::New();
137    }
138}
139
140/**
141 * \brief Internal method to set up pipeline after a state change
142 */
143void Molecule::update()
144{
145    if (_dataSet == NULL) {
146        return;
147    }
148    vtkDataSet *ds = _dataSet->getVtkDataSet();
149
150    if (_atomMapper == NULL) {
151        _atomMapper = vtkSmartPointer<vtkGlyph3DMapper>::New();
152        _atomMapper->SetResolveCoincidentTopologyToPolygonOffset();
153        _atomMapper->ScalarVisibilityOn();
154    }
155    if (_bondMapper == NULL) {
156        _bondMapper = vtkSmartPointer<vtkGlyph3DMapper>::New();
157        _bondMapper->SetResolveCoincidentTopologyToPolygonOffset();
158        _bondMapper->ScalarVisibilityOn();
159    }
160    if (_labelMapper == NULL) {
161        _labelMapper = vtkSmartPointer<vtkLabelPlacementMapper>::New();
162        _labelMapper->SetShapeToRoundedRect();
163        _labelMapper->SetBackgroundColor(1.0, 1.0, 0.7);
164        _labelMapper->SetBackgroundOpacity(0.8);
165        _labelMapper->SetMargin(3);
166    }
167
168    if (_lut == NULL) {
169        if (ds->GetPointData() == NULL ||
170            ds->GetPointData()->GetScalars() == NULL ||
171            strcmp(ds->GetPointData()->GetScalars()->GetName(), "element") != 0) {
172            WARN("No element array in dataset %s", _dataSet->getName().c_str());
173            if (ds->GetPointData()->GetScalars() == NULL) {
174                _colorMode = COLOR_CONSTANT;
175            } else {
176                _colorMode = COLOR_BY_SCALAR;
177            }
178            _atomScaling = NO_ATOM_SCALING;
179            setColorMap(ColorMap::getDefault());
180        } else {
181            TRACE("Using element default colormap");
182            setColorMap(ColorMap::getElementDefault());
183        }
184    }
185
186    initProp();
187
188    addLabelArray(ds);
189    addRadiusArray(ds, _atomScaling, _radiusScale);
190    computeBonds(ds);
191
192    vtkPolyData *pd = vtkPolyData::SafeDownCast(ds);
193    if (pd) {
194        TRACE("Points: %d Verts: %d Lines: %d Polys: %d Strips: %d",
195              pd->GetNumberOfPoints(),
196              pd->GetNumberOfVerts(),
197              pd->GetNumberOfLines(),
198              pd->GetNumberOfPolys(),
199              pd->GetNumberOfStrips());
200        // DataSet is a vtkPolyData
201        if (pd->GetNumberOfLines() > 0) {
202            // Bonds
203            setupBondPolyData();
204
205            if (_cylinderSource == NULL) {
206                _cylinderSource = vtkSmartPointer<vtkCylinderSource>::New();
207                _cylinderSource->SetRadius(0.075);
208                _cylinderSource->SetHeight(1.0);
209                _cylinderSource->SetResolution(12);
210                _cylinderSource->CappingOff();
211                _cylinderTrans = vtkSmartPointer<vtkTransformPolyDataFilter>::New();
212                _cylinderTrans->SetInputConnection(_cylinderSource->GetOutputPort());
213                vtkSmartPointer<vtkTransform> trans = vtkSmartPointer<vtkTransform>::New();
214                trans->RotateZ(-90.0);
215                _cylinderTrans->SetTransform(trans);
216            }
217            if (_lineSource == NULL) {
218                _lineSource = vtkSmartPointer<vtkLineSource>::New();
219                _lineSource->SetPoint1(-0.5, 0, 0);
220                _lineSource->SetPoint2(0.5, 0, 0);
221            }
222
223            _bondMapper->SetSourceConnection(_cylinderTrans->GetOutputPort());
224            _bondMapper->SetInputData(_bondPD);
225            _bondMapper->SetOrientationArray("bond_orientations");
226            _bondMapper->SetOrientationModeToDirection();
227            _bondMapper->OrientOn();
228            _bondMapper->SetScaleArray("bond_scales");
229            _bondMapper->SetScaleModeToScaleByVectorComponents();
230            _bondMapper->ScalingOn();
231            _bondMapper->ClampingOff();
232
233            _bondProp->SetMapper(_bondMapper);
234            getAssembly()->AddPart(_bondProp);
235        }
236        if (pd->GetNumberOfPoints() > 0) {
237            if (_labelHierarchy == NULL) {
238                _labelHierarchy = vtkSmartPointer<vtkPointSetToLabelHierarchy>::New();
239            }
240            if (_labelTransform == NULL) {
241                _labelTransform = vtkSmartPointer<vtkTransform>::New();
242            }
243            vtkSmartPointer<vtkTransformPolyDataFilter> transformFilter = vtkSmartPointer<vtkTransformPolyDataFilter>::New();
244            transformFilter->SetInputData(pd);
245            transformFilter->SetTransform(_labelTransform);
246            _labelHierarchy->SetInputConnection(transformFilter->GetOutputPort());
247            _labelHierarchy->SetLabelArrayName("_atom_labels");
248            _labelHierarchy->GetTextProperty()->SetColor(0, 0, 0);
249            _labelMapper->SetInputConnection(_labelHierarchy->GetOutputPort());
250            _labelProp->SetMapper(_labelMapper);
251
252            // Atoms
253            if (_sphereSource == NULL) {
254                _sphereSource = vtkSmartPointer<vtkSphereSource>::New();
255                _sphereSource->SetRadius(1.0);
256                _sphereSource->SetThetaResolution(14);
257                _sphereSource->SetPhiResolution(14);
258            }
259
260            _atomMapper->SetSourceConnection(_sphereSource->GetOutputPort());
261            _atomMapper->SetInputData(pd);
262            _atomMapper->SetScaleArray("_radii");
263            _atomMapper->SetScaleModeToScaleByMagnitude();
264            _atomMapper->ScalingOn();
265            _atomMapper->OrientOff();
266
267            _atomProp->SetMapper(_atomMapper);
268            getAssembly()->AddPart(_atomProp);
269        }
270    } else {
271        // DataSet is NOT a vtkPolyData
272        ERROR("DataSet is not a PolyData");
273        return;
274    }
275
276    setAtomLabelVisibility(_labelsOn);
277
278    if (pd->GetNumberOfPoints() > 0) {
279        _atomMapper->Update();
280        _labelMapper->Update();
281    }
282    if (pd->GetNumberOfLines() > 0) {
283        _bondMapper->Update();
284    }
285}
286
287void Molecule::updateLabelTransform()
288{
289    if (_labelTransform != NULL && getAssembly() != NULL) {
290        _labelTransform->SetMatrix(getAssembly()->GetMatrix());
291    }
292}
293
294void Molecule::updateRanges(Renderer *renderer)
295{
296    if (_dataSet == NULL) {
297        ERROR("called before setDataSet");
298        return;
299    }
300
301    if (renderer->getUseCumulativeRange()) {
302        renderer->getCumulativeDataRange(_dataRange,
303                                         _dataSet->getActiveScalarsName(),
304                                         1);
305        renderer->getCumulativeDataRange(_vectorMagnitudeRange,
306                                         _dataSet->getActiveVectorsName(),
307                                         3);
308        for (int i = 0; i < 3; i++) {
309            renderer->getCumulativeDataRange(_vectorComponentRange[i],
310                                             _dataSet->getActiveVectorsName(),
311                                             3, i);
312        }
313    } else {
314        _dataSet->getScalarRange(_dataRange);
315        _dataSet->getVectorRange(_vectorMagnitudeRange);
316        for (int i = 0; i < 3; i++) {
317            _dataSet->getVectorRange(_vectorComponentRange[i], i);
318        }
319    }
320
321    // Need to update color map ranges
322    double *rangePtr = _colorFieldRange;
323    if (_colorFieldRange[0] > _colorFieldRange[1]) {
324        rangePtr = NULL;
325    }
326    setColorMode(_colorMode, _colorFieldType, _colorFieldName.c_str(), rangePtr);
327}
328
329void Molecule::setColorMode(ColorMode mode)
330{
331    _colorMode = mode;
332    if (_dataSet == NULL)
333        return;
334
335    switch (mode) {
336    case COLOR_BY_ELEMENTS: // Assume default scalar is "element" array
337        setColorMode(mode,
338                     DataSet::POINT_DATA,
339                     "element");
340        break;
341    case COLOR_BY_SCALAR:
342        setColorMode(mode,
343                     _dataSet->getActiveScalarsType(),
344                     _dataSet->getActiveScalarsName());
345        break;
346    case COLOR_BY_VECTOR_MAGNITUDE:
347        setColorMode(mode,
348                     _dataSet->getActiveVectorsType(),
349                     _dataSet->getActiveVectorsName());
350        break;
351    case COLOR_BY_VECTOR_X:
352        setColorMode(mode,
353                     _dataSet->getActiveVectorsType(),
354                     _dataSet->getActiveVectorsName());
355        break;
356    case COLOR_BY_VECTOR_Y:
357        setColorMode(mode,
358                     _dataSet->getActiveVectorsType(),
359                     _dataSet->getActiveVectorsName());
360        break;
361    case COLOR_BY_VECTOR_Z:
362        setColorMode(mode,
363                     _dataSet->getActiveVectorsType(),
364                     _dataSet->getActiveVectorsName());
365        break;
366    case COLOR_CONSTANT:
367    default:
368        setColorMode(mode, DataSet::POINT_DATA, NULL, NULL);
369        break;
370    }
371}
372
373void Molecule::setColorMode(ColorMode mode,
374                            const char *name, double range[2])
375{
376    if (_dataSet == NULL)
377        return;
378    DataSet::DataAttributeType type = DataSet::POINT_DATA;
379    int numComponents = 1;
380    if (name != NULL && strlen(name) > 0 &&
381        !_dataSet->getFieldInfo(name, &type, &numComponents)) {
382        ERROR("Field not found: %s", name);
383        return;
384    }
385    setColorMode(mode, type, name, range);
386}
387
388void Molecule::setColorMode(ColorMode mode, DataSet::DataAttributeType type,
389                            const char *name, double range[2])
390{
391    _colorMode = mode;
392    _colorFieldType = type;
393    if (name == NULL)
394        _colorFieldName.clear();
395    else
396        _colorFieldName = name;
397    if (range == NULL) {
398        _colorFieldRange[0] = DBL_MAX;
399        _colorFieldRange[1] = -DBL_MAX;
400    } else {
401        memcpy(_colorFieldRange, range, sizeof(double)*2);
402    }
403
404    if (_dataSet == NULL || (_atomMapper == NULL && _bondMapper == NULL))
405        return;
406
407    switch (type) {
408    case DataSet::POINT_DATA:
409        if (_atomMapper != NULL)
410            _atomMapper->SetScalarModeToUsePointFieldData();
411        if (_bondMapper != NULL)
412            _bondMapper->SetScalarModeToUsePointFieldData();
413        break;
414    case DataSet::CELL_DATA:
415        if (_atomMapper != NULL)
416            _atomMapper->SetScalarModeToUseCellFieldData();
417        if (_bondMapper != NULL)
418            _bondMapper->SetScalarModeToUseCellFieldData();
419        break;
420    default:
421        ERROR("Unsupported DataAttributeType: %d", type);
422        return;
423    }
424
425    if (name != NULL && strlen(name) > 0) {
426        if (_atomMapper != NULL)
427            _atomMapper->SelectColorArray(name);
428        if (_bondMapper != NULL)
429            _bondMapper->SelectColorArray(name);
430    } else {
431        if (_atomMapper != NULL)
432            _atomMapper->SetScalarModeToDefault();
433        if (_bondMapper != NULL)
434            _bondMapper->SetScalarModeToDefault();
435    }
436
437    if (_lut != NULL) {
438        if (range != NULL) {
439            _lut->SetRange(range);
440        } else if (mode == COLOR_BY_ELEMENTS) {
441            double range[2];
442            range[0] = 0;
443            range[1] = NUM_ELEMENTS;
444            _lut->SetRange(range);
445        } else if (name != NULL && strlen(name) > 0) {
446            double r[2];
447            int comp = -1;
448            if (mode == COLOR_BY_VECTOR_X)
449                comp = 0;
450            else if (mode == COLOR_BY_VECTOR_Y)
451                comp = 1;
452            else if (mode == COLOR_BY_VECTOR_Z)
453                comp = 2;
454
455            if (_renderer->getUseCumulativeRange()) {
456                int numComponents;
457                if  (!_dataSet->getFieldInfo(name, type, &numComponents)) {
458                    ERROR("Field not found: %s, type: %d", name, type);
459                    return;
460                } else if (numComponents < comp+1) {
461                    ERROR("Request for component %d in field with %d components",
462                          comp, numComponents);
463                    return;
464                }
465                _renderer->getCumulativeDataRange(r, name, type, numComponents, comp);
466            } else {
467                _dataSet->getDataRange(r, name, type, comp);
468            }
469            _lut->SetRange(r);
470        } else {
471            switch (mode) {
472            case COLOR_CONSTANT:
473            case COLOR_BY_SCALAR:
474                _lut->SetRange(_dataRange);
475                break;
476            case COLOR_BY_VECTOR_MAGNITUDE:
477                _lut->SetRange(_vectorMagnitudeRange);
478                break;
479            case COLOR_BY_VECTOR_X:
480                _lut->SetRange(_vectorComponentRange[0]);
481                break;
482            case COLOR_BY_VECTOR_Y:
483                _lut->SetRange(_vectorComponentRange[1]);
484                break;
485            case COLOR_BY_VECTOR_Z:
486                _lut->SetRange(_vectorComponentRange[2]);
487                break;
488            default:
489                break;
490            }
491        }
492    }
493
494    switch (mode) {
495    case COLOR_BY_ELEMENTS:
496        if (_atomMapper != NULL)
497            _atomMapper->ScalarVisibilityOn();
498        if (_bondMapper != NULL)
499            _bondMapper->ScalarVisibilityOn();
500        break;
501    case COLOR_BY_SCALAR:
502        if (_atomMapper != NULL)
503            _atomMapper->ScalarVisibilityOn();
504        if (_bondMapper != NULL)
505            _bondMapper->ScalarVisibilityOn();
506        break;
507    case COLOR_BY_VECTOR_MAGNITUDE:
508        if (_atomMapper != NULL)
509            _atomMapper->ScalarVisibilityOn();
510        if (_bondMapper != NULL)
511            _bondMapper->ScalarVisibilityOn();
512        if (_lut != NULL) {
513            _lut->SetVectorModeToMagnitude();
514        }
515        break;
516    case COLOR_BY_VECTOR_X:
517        if (_atomMapper != NULL)
518            _atomMapper->ScalarVisibilityOn();
519        if (_bondMapper != NULL)
520            _bondMapper->ScalarVisibilityOn();
521        if (_lut != NULL) {
522            _lut->SetVectorModeToComponent();
523            _lut->SetVectorComponent(0);
524        }
525        break;
526    case COLOR_BY_VECTOR_Y:
527        if (_atomMapper != NULL)
528            _atomMapper->ScalarVisibilityOn();
529        if (_bondMapper != NULL)
530            _bondMapper->ScalarVisibilityOn();
531        if (_lut != NULL) {
532            _lut->SetVectorModeToComponent();
533            _lut->SetVectorComponent(1);
534        }
535        break;
536    case COLOR_BY_VECTOR_Z:
537        if (_atomMapper != NULL)
538            _atomMapper->ScalarVisibilityOn();
539        if (_bondMapper != NULL)
540            _bondMapper->ScalarVisibilityOn();
541        if (_lut != NULL) {
542            _lut->SetVectorModeToComponent();
543            _lut->SetVectorComponent(2);
544        }
545        break;
546    case COLOR_CONSTANT:
547    default:
548        if (_atomMapper != NULL)
549            _atomMapper->ScalarVisibilityOff();
550        if (_bondMapper != NULL)
551            _bondMapper->ScalarVisibilityOff();
552        break;
553    }
554}
555
556/**
557 * \brief Called when the color map has been edited
558 */
559void Molecule::updateColorMap()
560{
561    setColorMap(_colorMap);
562}
563
564/**
565 * \brief Associate a colormap lookup table with the DataSet
566 */
567void Molecule::setColorMap(ColorMap *cmap)
568{
569    if (cmap == NULL)
570        return;
571
572    _colorMap = cmap;
573 
574    if (_lut == NULL) {
575        _lut = vtkSmartPointer<vtkLookupTable>::New();
576        if (_atomMapper != NULL) {
577            _atomMapper->UseLookupTableScalarRangeOn();
578            _atomMapper->SetLookupTable(_lut);
579        }
580        if (_bondMapper != NULL) {
581            _bondMapper->UseLookupTableScalarRangeOn();
582            _bondMapper->SetLookupTable(_lut);
583        }
584        _lut->DeepCopy(cmap->getLookupTable());
585        switch (_colorMode) {
586        case COLOR_BY_ELEMENTS: {
587            double range[2];
588            range[0] = 0;
589            range[1] = NUM_ELEMENTS;
590            _lut->SetRange(range);
591        }
592            break;
593        case COLOR_CONSTANT:
594        case COLOR_BY_SCALAR:
595            _lut->SetRange(_dataRange);
596            break;
597        case COLOR_BY_VECTOR_MAGNITUDE:
598            _lut->SetRange(_vectorMagnitudeRange);
599            break;
600        case COLOR_BY_VECTOR_X:
601            _lut->SetRange(_vectorComponentRange[0]);
602            break;
603        case COLOR_BY_VECTOR_Y:
604            _lut->SetRange(_vectorComponentRange[1]);
605            break;
606        case COLOR_BY_VECTOR_Z:
607            _lut->SetRange(_vectorComponentRange[2]);
608            break;
609        default:
610            break;
611        }
612    } else {
613        double range[2];
614        _lut->GetTableRange(range);
615        _lut->DeepCopy(cmap->getLookupTable());
616        _lut->SetRange(range);
617        _lut->Modified();
618    }
619}
620
621void Molecule::setAtomLabelField(const char *fieldName)
622{
623    if (_labelHierarchy != NULL) {
624        if (strcmp(fieldName, "default") == 0) {
625            _labelHierarchy->SetLabelArrayName("_atom_labels");
626        } else {
627            _labelHierarchy->SetLabelArrayName(fieldName);
628        }
629    }
630}
631
632/**
633 * \brief Turn on/off rendering of atom labels
634 */
635void Molecule::setAtomLabelVisibility(bool state)
636{
637    _labelsOn = state;
638    if (_labelProp != NULL) {
639        _labelProp->SetVisibility((state ? 1 : 0));
640    }
641}
642
643/**
644 * \brief Turn on/off rendering of the atoms
645 */
646void Molecule::setAtomVisibility(bool state)
647{
648    if (_atomProp != NULL) {
649        _atomProp->SetVisibility((state ? 1 : 0));
650    }
651}
652
653/**
654 * \brief Turn on/off rendering of the bonds
655 */
656void Molecule::setBondVisibility(bool state)
657{
658    if (_bondProp != NULL) {
659        _bondProp->SetVisibility((state ? 1 : 0));
660    }
661}
662
663/**
664 * \brief Toggle visibility of the prop
665 */
666void Molecule::setVisibility(bool state)
667{
668    GraphicsObject::setVisibility(state);
669    if (_labelProp != NULL) {
670        if (!state)
671            _labelProp->SetVisibility(0);
672        else
673            setAtomLabelVisibility(_labelsOn);
674    }
675}
676
677/**
678 * \brief Set opacity of molecule
679 */
680void Molecule::setOpacity(double opacity)
681{
682    GraphicsObject::setOpacity(opacity);
683    if (_labelMapper != NULL) {
684        _labelMapper->SetBackgroundOpacity(opacity);
685    }
686}
687
688void Molecule::setAtomQuality(double quality)
689{
690    if (_sphereSource == NULL)
691        return;
692
693    if (quality > 10.0)
694        quality = 10.0;
695
696    int thetaRes = (int)(quality * 14.0);
697    int phiRes = (int)(quality * 14.0);
698    if (thetaRes < 4) thetaRes = 4;
699    if (phiRes < 3) phiRes = 3;
700
701    _sphereSource->SetThetaResolution(thetaRes);
702    _sphereSource->SetPhiResolution(phiRes);
703
704    if (_atomMapper != NULL) {
705        _atomMapper->Modified();
706        _atomMapper->Update();
707    }
708}
709
710void Molecule::setBondQuality(double quality)
711{
712    if (_cylinderSource == NULL)
713        return;
714
715    if (quality > 10.0)
716        quality = 10.0;
717
718    int res = (int)(quality * 12.0);
719    if (res < 3) res = 3;
720
721    _cylinderSource->SetResolution(res);
722
723    if (_bondMapper != NULL) {
724        _bondMapper->Modified();
725        _bondMapper->Update();
726    }
727}
728
729void Molecule::setBondStyle(BondStyle style)
730{
731    switch (style) {
732    case BOND_STYLE_CYLINDER:
733        if (_bondProp != NULL) {
734            _bondProp->GetProperty()->SetLineWidth(_edgeWidth);
735            _bondProp->GetProperty()->SetLighting(_lighting ? 1 : 0);
736        }
737        if (_bondMapper != NULL && _cylinderTrans != NULL &&
738            _bondMapper->GetInputConnection(1, 0) != _cylinderTrans->GetOutputPort()) {
739            _cylinderTrans->Modified();
740            _bondMapper->SetSourceConnection(_cylinderTrans->GetOutputPort());
741            _bondMapper->Modified();
742        }
743        break;
744    case BOND_STYLE_LINE:
745        if (_bondProp != NULL) {
746            _bondProp->GetProperty()->LightingOff();
747        }
748        if (_bondMapper != NULL && _lineSource != NULL &&
749            _bondMapper->GetInputConnection(1, 0) != _lineSource->GetOutputPort()) {
750            _lineSource->Modified();
751            _bondMapper->SetSourceConnection(_lineSource->GetOutputPort());
752            _bondMapper->Modified();
753        }
754         break;
755    default:
756        WARN("Unknown bond style");   
757    }
758}
759
760/**
761 * \brief Set constant bond color
762 */
763void Molecule::setBondColor(float color[3])
764{
765    _bondColor[0] = color[0];
766    _bondColor[1] = color[1];
767    _bondColor[2] = color[2];
768    if (_bondProp != NULL) {
769        _bondProp->GetProperty()->SetColor(_bondColor[0], _bondColor[1], _bondColor[2]);
770    }
771}
772
773/**
774 * \brief Set mode to determine how bonds are colored
775 */
776void Molecule::setBondColorMode(BondColorMode mode)
777{
778    if (_bondMapper == NULL) return;
779
780    switch (mode) {
781    case BOND_COLOR_BY_ELEMENTS:
782        _bondMapper->ScalarVisibilityOn();
783        break;
784    case BOND_COLOR_CONSTANT:
785        _bondMapper->ScalarVisibilityOff();
786        break;
787    default:
788        WARN("Unknown bond color mode");
789    }
790}
791
792/**
793 * \brief Set a group of world coordinate planes to clip rendering
794 *
795 * Passing NULL for planes will remove all cliping planes
796 */
797void Molecule::setClippingPlanes(vtkPlaneCollection *planes)
798{
799    if (_atomMapper != NULL) {
800        _atomMapper->SetClippingPlanes(planes);
801    }
802    if (_bondMapper != NULL) {
803        _bondMapper->SetClippingPlanes(planes);
804    }
805}
806
807/**
808 * \brief Set the radius type used for scaling atoms
809 */
810void Molecule::setAtomScaling(AtomScaling state)
811{
812    if (state == _atomScaling)
813        return;
814
815    _atomScaling = state;
816    if (_dataSet != NULL) {
817        vtkDataSet *ds = _dataSet->getVtkDataSet();
818        addRadiusArray(ds, _atomScaling, _radiusScale);
819        if (_atomMapper != NULL) {
820             assert(ds->GetPointData() != NULL &&
821                    ds->GetPointData()->GetArray("_radii") != NULL);
822            _atomMapper->SetScaleModeToScaleByMagnitude();
823            _atomMapper->SetScaleArray("_radii");
824            _atomMapper->ScalingOn();
825        }
826    }
827}
828
829/**
830 * \brief Set the constant radius scaling factor for atoms.  This
831 * can be used to convert from Angstroms to atom coordinates units.
832 */
833void Molecule::setAtomRadiusScale(double scale)
834{
835    if (scale == _radiusScale)
836        return;
837
838    _radiusScale = scale;
839    if (_dataSet != NULL) {
840        vtkDataSet *ds = _dataSet->getVtkDataSet();
841        addRadiusArray(ds, _atomScaling, _radiusScale);
842    }
843}
844
845/**
846 * \brief Set the constant radius scaling factor for bonds.
847 */
848void Molecule::setBondRadiusScale(double scale)
849{
850    if (_cylinderSource != NULL &&
851        _cylinderSource->GetRadius() != scale) {
852        _cylinderSource->SetRadius(scale);
853        // Workaround bug with source modification not causing
854        // mapper to be updated
855        if (_bondMapper != NULL) {
856            _bondMapper->Modified();
857        }
858    }
859}
860
861void Molecule::setupBondPolyData()
862{
863    if (_dataSet == NULL)
864        return;
865    if (_bondPD == NULL) {
866        _bondPD = vtkSmartPointer<vtkPolyData>::New();
867    } else {
868        _bondPD->Initialize();
869    }
870    vtkPolyData *pd = vtkPolyData::SafeDownCast(_dataSet->getVtkDataSet());
871    if (pd == NULL)
872        return;
873    vtkCellArray *lines = pd->GetLines();
874    lines->InitTraversal();
875    vtkIdType npts, *pts;
876    vtkSmartPointer<vtkPoints> bondPoints = vtkSmartPointer<vtkPoints>::New();
877    vtkSmartPointer<vtkIntArray> bondElements = vtkSmartPointer<vtkIntArray>::New();
878    vtkSmartPointer<vtkDoubleArray> bondVectors = vtkSmartPointer<vtkDoubleArray>::New();
879    vtkSmartPointer<vtkDoubleArray> bondScales = vtkSmartPointer<vtkDoubleArray>::New();
880    bondElements->SetName("element");
881    bondVectors->SetName("bond_orientations");
882    bondVectors->SetNumberOfComponents(3);
883    bondScales->SetName("bond_scales");
884    bondScales->SetNumberOfComponents(3);
885    vtkDataArray *elements = NULL;
886    if (pd->GetPointData() != NULL &&
887        pd->GetPointData()->GetScalars() != NULL &&
888        strcmp(pd->GetPointData()->GetScalars()->GetName(), "element") == 0) {
889        elements = pd->GetPointData()->GetScalars();
890    }
891    for (int i = 0; lines->GetNextCell(npts, pts); i++) {
892        assert(npts == 2);
893        double pt0[3], pt1[3];
894        double newPt0[3], newPt1[3];
895        //double *pt0 = pd->GetPoint(pts[0]);
896        //double *pt1 = pd->GetPoint(pts[1]);
897        pd->GetPoint(pts[0], pt0);
898        pd->GetPoint(pts[1], pt1);
899        double center[3];
900
901        for (int j = 0; j < 3; j++)
902            center[j] = pt0[j] + (pt1[j] - pt0[j]) * 0.5;
903        for (int j = 0; j < 3; j++)
904            newPt0[j] = pt0[j] + (pt1[j] - pt0[j]) * 0.25;
905        for (int j = 0; j < 3; j++)
906            newPt1[j] = pt0[j] + (pt1[j] - pt0[j]) * 0.75;
907
908        bondPoints->InsertNextPoint(newPt0);
909        bondPoints->InsertNextPoint(newPt1);
910
911#ifdef DEBUG
912        TRACE("Bond %d: (%g,%g,%g)-(%g,%g,%g)-(%g,%g,%g)", i,
913              pt0[0], pt0[1], pt0[2],
914              center[0], center[1], center[2],
915              pt1[0], pt1[1], pt1[2]);
916#endif
917
918        double vec[3];
919        for (int j = 0; j < 3; j++)
920            vec[j] = center[j] - pt0[j];
921
922        bondVectors->InsertNextTupleValue(vec);
923        bondVectors->InsertNextTupleValue(vec);
924#ifdef DEBUG
925        TRACE("Bond %d, vec: %g,%g,%g", i, vec[0], vec[1], vec[2]);
926#endif
927
928        double scale[3];
929        scale[0] = sqrt(vec[0]*vec[0] + vec[1]*vec[1] + vec[2]*vec[2]);
930        scale[1] = 1.0;
931        scale[2] = 1.0;
932
933        bondScales->InsertNextTupleValue(scale);
934        bondScales->InsertNextTupleValue(scale);
935
936        if (elements != NULL) {
937            int element = (int)elements->GetComponent(pts[0], 0);
938#ifdef DEBUG
939            TRACE("Bond %d, elt 0: %d", i, element);
940#endif
941            bondElements->InsertNextValue(element);
942            element = (int)elements->GetComponent(pts[1], 0);
943#ifdef DEBUG
944            TRACE("Bond %d, elt 1: %d", i, element);
945#endif
946            bondElements->InsertNextValue(element);
947        }
948    }
949    _bondPD->SetPoints(bondPoints);
950    if (elements != NULL)
951        _bondPD->GetPointData()->SetScalars(bondElements);
952    _bondPD->GetPointData()->AddArray(bondVectors);
953    _bondPD->GetPointData()->AddArray(bondScales);
954}
955
956void Molecule::addLabelArray(vtkDataSet *dataSet)
957{
958    vtkDataArray *elements = NULL;
959    if (dataSet->GetPointData() != NULL &&
960        dataSet->GetPointData()->GetScalars() != NULL &&
961        strcmp(dataSet->GetPointData()->GetScalars()->GetName(), "element") == 0) {
962        elements = dataSet->GetPointData()->GetScalars();
963    } else {
964        WARN("Can't label atoms without an element array");
965    }
966
967    vtkSmartPointer<vtkStringArray> labelArray = vtkSmartPointer<vtkStringArray>::New();
968    labelArray->SetName("_atom_labels");
969    vtkPolyData *pd = vtkPolyData::SafeDownCast(dataSet);
970    if (pd == NULL) {
971        ERROR("DataSet not a PolyData");
972        return;
973    }
974    for (int i = 0; i < pd->GetNumberOfPoints(); i++) {
975        char buf[32];
976        if (elements != NULL) {
977            int elt = (int)elements->GetComponent(i, 0);
978            sprintf(buf, "%s%d", g_elementNames[elt], i);
979        } else {
980            sprintf(buf, "%d", i);
981        }
982        labelArray->InsertNextValue(buf);
983    }
984    dataSet->GetPointData()->AddArray(labelArray);
985}
986
987/**
988 * \brief Add a scalar array to dataSet with sizes for the elements
989 * specified in the "element" scalar array
990 */
991void Molecule::addRadiusArray(vtkDataSet *dataSet, AtomScaling scaling, double scaleFactor)
992{
993    vtkDataArray *elements = NULL;
994    if (dataSet->GetPointData() != NULL &&
995        dataSet->GetPointData()->GetScalars() != NULL &&
996        strcmp(dataSet->GetPointData()->GetScalars()->GetName(), "element") == 0) {
997        elements = dataSet->GetPointData()->GetScalars();
998    } else if (scaling != NO_ATOM_SCALING) {
999        WARN("Can't use non-constant scaling without an element array");
1000    }
1001    const float *radiusSource = NULL;
1002    switch (scaling) {
1003    case VAN_DER_WAALS_RADIUS:
1004        radiusSource = g_vdwRadii;
1005        break;
1006    case COVALENT_RADIUS:
1007        radiusSource = g_covalentRadii;
1008        break;
1009    case ATOMIC_RADIUS:
1010        radiusSource = g_atomicRadii;
1011        break;
1012    case NO_ATOM_SCALING:
1013    default:
1014        ;
1015    }
1016    vtkPolyData *pd = vtkPolyData::SafeDownCast(dataSet);
1017    if (pd == NULL) {
1018        ERROR("DataSet not a PolyData");
1019        return;
1020    }
1021    vtkDataArray *array = dataSet->GetPointData()->GetArray("_radii");
1022    vtkSmartPointer<vtkFloatArray> radii;
1023    if (array == NULL) {
1024        radii = vtkSmartPointer<vtkFloatArray>::New();
1025        radii->SetName("_radii");
1026        radii->SetNumberOfComponents(1);
1027        radii->SetNumberOfValues(pd->GetNumberOfPoints());
1028    } else {
1029        radii = vtkFloatArray::SafeDownCast(array);
1030        assert(radii != NULL);
1031    }
1032    for (int i = 0; i < pd->GetNumberOfPoints(); i++) {
1033        float value;
1034        if (elements != NULL && radiusSource != NULL) {
1035            int elt = (int)elements->GetComponent(i, 0);
1036            value = radiusSource[elt] * scaleFactor;
1037        } else {
1038            value = scaleFactor;
1039        }
1040        radii->SetValue(i, value);
1041    }
1042    radii->Modified();
1043    if (array == NULL) {
1044        dataSet->GetPointData()->AddArray(radii);
1045    }
1046}
1047
1048/**
1049 * \brief Create a color map to map atomic numbers to element colors
1050 */
1051ColorMap *Molecule::createElementColorMap()
1052{
1053    ColorMap *elementCmap = new ColorMap("elementDefault");
1054    ColorMap::ControlPoint cp[NUM_ELEMENTS+1];
1055
1056    elementCmap->setNumberOfTableEntries(NUM_ELEMENTS+1);
1057    for (int i = 0; i <= NUM_ELEMENTS; i++) {
1058        cp[i].value = i/((double)NUM_ELEMENTS);
1059        for (int c = 0; c < 3; c++) {
1060            cp[i].color[c] = ((double)g_elementColors[i][c])/255.;
1061        }
1062        elementCmap->addControlPoint(cp[i]);
1063    }
1064    ColorMap::OpacityControlPoint ocp[2];
1065    ocp[0].value = 0;
1066    ocp[0].alpha = 1.0;
1067    ocp[1].value = 1.0;
1068    ocp[1].alpha = 1.0;
1069    elementCmap->addOpacityControlPoint(ocp[0]);
1070    elementCmap->addOpacityControlPoint(ocp[1]);
1071    elementCmap->build();
1072    double range[2];
1073    range[0] = 0;
1074    range[1] = NUM_ELEMENTS;
1075    elementCmap->getLookupTable()->SetRange(range);
1076
1077    return elementCmap;
1078}
1079
1080int Molecule::computeBonds(vtkDataSet *dataSet, float tolerance)
1081{
1082    vtkPolyData *pd = vtkPolyData::SafeDownCast(dataSet);
1083    if (pd == NULL) {
1084        ERROR("DataSet not a PolyData");
1085        return -1;
1086    }
1087
1088    if (pd->GetNumberOfPoints() < 2 ||
1089        pd->GetNumberOfLines() > 0) {
1090        return 0;
1091    }
1092
1093    vtkDataArray *elements = NULL;
1094    if (pd->GetPointData() != NULL &&
1095        pd->GetPointData()->GetScalars() != NULL &&
1096        strcmp(pd->GetPointData()->GetScalars()->GetName(), "element") == 0) {
1097        elements = pd->GetPointData()->GetScalars();
1098    } else {
1099        TRACE("Can't compute bonds without an element array");
1100        return -1;
1101    }
1102    int numAtoms = pd->GetNumberOfPoints();
1103    vtkSmartPointer<vtkCellArray> cells = vtkSmartPointer<vtkCellArray>::New();
1104    for (int i = 0; i < numAtoms; i++) {
1105        double pt1[3];
1106        int elem1 = (int)elements->GetComponent(i, 0);
1107        float r1 = g_covalentRadii[elem1];
1108        pd->GetPoint(i, pt1);
1109        for (int j = i+1; j < numAtoms; j++) {
1110            float r2, dist;
1111            double pt2[3];
1112            int elem2 = (int)elements->GetComponent(j, 0);
1113
1114            if (elem1 == 1 && elem2 == 1)
1115                continue;
1116
1117            r2 = g_covalentRadii[elem2];
1118            pd->GetPoint(j, pt2);
1119
1120            double x = pt1[0] - pt2[0];
1121            double y = pt1[1] - pt2[1];
1122            double z = pt1[2] - pt2[2];
1123            dist = (float)sqrt(x*x + y*y + z*z);
1124            //TRACE("i=%d,j=%d elem1=%d,elem2=%d r1=%g,r2=%g, dist=%g", i, j, elem1, elem2, r1, r2, dist);
1125            if (dist < (r1 + r2 + tolerance)) {
1126                vtkIdType pts[2];
1127                pts[0] = i;
1128                pts[1] = j;
1129                cells->InsertNextCell(2, pts);
1130            }
1131        }
1132    }
1133
1134    if (cells->GetNumberOfCells() > 0) {
1135        pd->SetLines(cells);
1136    }
1137
1138    TRACE("Generated %d bonds", cells->GetNumberOfCells());
1139    return cells->GetNumberOfCells();
1140}
Note: See TracBrowser for help on using the repository browser.