Changeset 3141 for trunk


Ignore:
Timestamp:
Aug 16, 2012 6:30:59 PM (12 years ago)
Author:
ldelgass
Message:

First pass at converting VTK Molecule to use display lists (Glyph3DMapper).
Still needs fixes for bond scaling.

Location:
trunk/packages/vizservers/vtkvis
Files:
2 edited

Legend:

Unmodified
Added
Removed
  • trunk/packages/vizservers/vtkvis/RpMolecule.cpp

    r3132 r3141  
    1010
    1111#include <vtkDataSet.h>
     12#include <vtkCellArray.h>
    1213#include <vtkPointData.h>
    1314#include <vtkFloatArray.h>
     15#include <vtkDoubleArray.h>
     16#include <vtkIntArray.h>
    1417#include <vtkStringArray.h>
    1518#include <vtkPolyData.h>
     
    1922#include <vtkTubeFilter.h>
    2023#include <vtkSphereSource.h>
     24#include <vtkCylinderSource.h>
     25#include <vtkTransform.h>
     26#include <vtkTransformPolyDataFilter.h>
    2127#include <vtkGlyph3D.h>
    2228#include <vtkPointSetToLabelHierarchy.h>
     
    2632#include "RpMolecule.h"
    2733#include "RpMoleculeData.h"
     34#include "RpVtkRenderer.h"
    2835#include "Trace.h"
    2936
     
    3441    _radiusScale(0.3),
    3542    _atomScaling(VAN_DER_WAALS_RADIUS),
    36     _colorMap(NULL),
    37     _labelsOn(false)
     43    _labelsOn(false),
     44    _colorMap(NULL)
    3845{
    3946    _bondColor[0] = _bondColor[1] = _bondColor[2] = 1.0f;
     
    4956        TRACE("Deleting Molecule with NULL DataSet");
    5057#endif
     58}
     59
     60void Molecule::setDataSet(DataSet *dataSet,
     61                          Renderer *renderer)
     62{
     63    if (_dataSet != dataSet) {
     64        _dataSet = dataSet;
     65
     66        _renderer = renderer;
     67
     68        if (renderer->getUseCumulativeRange()) {
     69            renderer->getCumulativeDataRange(_dataRange,
     70                                             _dataSet->getActiveScalarsName(),
     71                                             1);
     72            renderer->getCumulativeDataRange(_vectorMagnitudeRange,
     73                                             _dataSet->getActiveVectorsName(),
     74                                             3);
     75            for (int i = 0; i < 3; i++) {
     76                renderer->getCumulativeDataRange(_vectorComponentRange[i],
     77                                                 _dataSet->getActiveVectorsName(),
     78                                                 3, i);
     79            }
     80        } else {
     81            _dataSet->getScalarRange(_dataRange);
     82            _dataSet->getVectorRange(_vectorMagnitudeRange);
     83            for (int i = 0; i < 3; i++) {
     84                _dataSet->getVectorRange(_vectorComponentRange[i], i);
     85            }
     86        }
     87
     88        update();
     89    }
    5190}
    5291
     
    99138    vtkDataSet *ds = _dataSet->getVtkDataSet();
    100139
     140#ifndef MOLECULE_USE_GLYPH3D_MAPPER
    101141    if (_atomMapper == NULL) {
    102142        _atomMapper = vtkSmartPointer<vtkPolyDataMapper>::New();
     
    109149        _bondMapper->ScalarVisibilityOn();
    110150    }
     151#else
     152    if (_atomMapper == NULL) {
     153        _atomMapper = vtkSmartPointer<vtkGlyph3DMapper>::New();
     154        _atomMapper->SetResolveCoincidentTopologyToPolygonOffset();
     155        _atomMapper->ScalarVisibilityOn();
     156    }
     157    if (_bondMapper == NULL) {
     158        _bondMapper = vtkSmartPointer<vtkGlyph3DMapper>::New();
     159        _bondMapper->SetResolveCoincidentTopologyToPolygonOffset();
     160        _bondMapper->ScalarVisibilityOn();
     161    }
     162#endif
    111163    if (_labelMapper == NULL) {
    112164        _labelMapper = vtkSmartPointer<vtkLabelPlacementMapper>::New();
     
    146198        if (pd->GetNumberOfLines() > 0) {
    147199            // Bonds
     200#ifndef MOLECULE_USE_GLYPH3D_MAPPER
    148201            if (_tuber == NULL)
    149202                _tuber = vtkSmartPointer<vtkTubeFilter>::New();
     
    154207            _tuber->SetVaryRadiusToVaryRadiusOff();
    155208            _bondMapper->SetInputConnection(_tuber->GetOutputPort());
     209#else
     210            setupBondPolyData();
     211
     212            if (_cylinderSource == NULL) {
     213                _cylinderSource = vtkSmartPointer<vtkCylinderSource>::New();
     214                _cylinderSource->SetRadius(0.075);
     215                _cylinderSource->SetHeight(1.0);
     216                _cylinderSource->SetResolution(12);
     217                _cylinderSource->CappingOff();
     218            }
     219            vtkSmartPointer<vtkTransformPolyDataFilter> transPD = vtkSmartPointer<vtkTransformPolyDataFilter>::New();
     220            transPD->SetInputConnection(_cylinderSource->GetOutputPort());
     221            vtkSmartPointer<vtkTransform> trans = vtkSmartPointer<vtkTransform>::New();
     222            trans->RotateZ(-90.0);
     223            transPD->SetTransform(trans);
     224
     225            _bondMapper->SetSourceConnection(transPD->GetOutputPort());
     226            _bondMapper->SetInputConnection(_bondPD->GetProducerPort());
     227            _bondMapper->SetOrientationArray(vtkDataSetAttributes::VECTORS);
     228            _bondMapper->SetOrientationModeToDirection();
     229            _bondMapper->OrientOn();
     230            _bondMapper->SetScaleArray(vtkDataSetAttributes::VECTORS);
     231            _bondMapper->SetScaleModeToScaleByMagnitude();
     232            _bondMapper->ScalingOn();
     233            _bondMapper->ClampingOff();
     234#endif
    156235            _bondProp->SetMapper(_bondMapper);
    157236            getAssembly()->AddPart(_bondProp);
     
    170249            sphereSource->SetThetaResolution(14);
    171250            sphereSource->SetPhiResolution(14);
     251#ifndef MOLECULE_USE_GLYPH3D_MAPPER
    172252            if (_glypher == NULL)
    173253                _glypher = vtkSmartPointer<vtkGlyph3D>::New();
    174254            _glypher->SetSourceConnection(sphereSource->GetOutputPort());
    175255            _glypher->SetInput(pd);
    176             if (_atomScaling != NO_ATOM_SCALING &&
    177                 ds->GetPointData() != NULL &&
     256            if (ds->GetPointData() != NULL &&
    178257                ds->GetPointData()->GetVectors() != NULL) {
    179258                _glypher->SetScaleModeToScaleByVector();
     
    186265            _glypher->OrientOff();
    187266            _atomMapper->SetInputConnection(_glypher->GetOutputPort());
     267#else
     268            _atomMapper->SetSourceConnection(sphereSource->GetOutputPort());
     269            _atomMapper->SetInputConnection(pd->GetProducerPort());
     270            if (ds->GetPointData() != NULL &&
     271                ds->GetPointData()->GetVectors() != NULL) {
     272                _atomMapper->SetScaleArray(vtkDataSetAttributes::VECTORS);
     273                _atomMapper->SetScaleModeToScaleByMagnitude();
     274                _atomMapper->ScalingOn();
     275            } else {
     276                _atomMapper->SetScaleModeToNoDataScaling();
     277                _atomMapper->ScalingOff();
     278            }
     279            _atomMapper->OrientOff();
     280#endif
    188281            _atomProp->SetMapper(_atomMapper);
    189282            getAssembly()->AddPart(_atomProp);
     
    218311}
    219312
     313void Molecule::setColorMode(ColorMode mode)
     314{
     315    _colorMode = mode;
     316    if (_dataSet == NULL)
     317        return;
     318
     319    switch (mode) {
     320    case COLOR_BY_ELEMENTS: // Assume default scalar is "element" array
     321    case COLOR_BY_SCALAR:
     322        setColorMode(mode,
     323                     _dataSet->getActiveScalarsType(),
     324                     _dataSet->getActiveScalarsName(),
     325                     _dataRange);
     326        break;
     327    case COLOR_BY_VECTOR_MAGNITUDE:
     328        setColorMode(mode,
     329                     _dataSet->getActiveVectorsType(),
     330                     _dataSet->getActiveVectorsName(),
     331                     _vectorMagnitudeRange);
     332        break;
     333    case COLOR_BY_VECTOR_X:
     334        setColorMode(mode,
     335                     _dataSet->getActiveVectorsType(),
     336                     _dataSet->getActiveVectorsName(),
     337                     _vectorComponentRange[0]);
     338        break;
     339    case COLOR_BY_VECTOR_Y:
     340        setColorMode(mode,
     341                     _dataSet->getActiveVectorsType(),
     342                     _dataSet->getActiveVectorsName(),
     343                     _vectorComponentRange[1]);
     344        break;
     345    case COLOR_BY_VECTOR_Z:
     346        setColorMode(mode,
     347                     _dataSet->getActiveVectorsType(),
     348                     _dataSet->getActiveVectorsName(),
     349                     _vectorComponentRange[2]);
     350        break;
     351    case COLOR_CONSTANT:
     352    default:
     353        setColorMode(mode, DataSet::POINT_DATA, NULL, NULL);
     354        break;
     355    }
     356}
     357
     358void Molecule::setColorMode(ColorMode mode,
     359                            const char *name, double range[2])
     360{
     361    if (_dataSet == NULL)
     362        return;
     363    DataSet::DataAttributeType type = DataSet::POINT_DATA;
     364    int numComponents = 1;
     365    if (name != NULL && strlen(name) > 0 &&
     366        !_dataSet->getFieldInfo(name, &type, &numComponents)) {
     367        ERROR("Field not found: %s", name);
     368        return;
     369    }
     370    setColorMode(mode, type, name, range);
     371}
     372
     373void Molecule::setColorMode(ColorMode mode, DataSet::DataAttributeType type,
     374                            const char *name, double range[2])
     375{
     376    _colorMode = mode;
     377    _colorFieldType = type;
     378    if (name == NULL)
     379        _colorFieldName.clear();
     380    else
     381        _colorFieldName = name;
     382    if (range == NULL) {
     383        _colorFieldRange[0] = DBL_MAX;
     384        _colorFieldRange[1] = -DBL_MAX;
     385    } else {
     386        memcpy(_colorFieldRange, range, sizeof(double)*2);
     387    }
     388
     389    if (_dataSet == NULL || (_atomMapper == NULL && _bondMapper == NULL))
     390        return;
     391
     392    switch (type) {
     393    case DataSet::POINT_DATA:
     394        if (_atomMapper != NULL)
     395            _atomMapper->SetScalarModeToUsePointFieldData();
     396        if (_bondMapper != NULL)
     397            _bondMapper->SetScalarModeToUsePointFieldData();
     398        break;
     399    case DataSet::CELL_DATA:
     400        if (_atomMapper != NULL)
     401            _atomMapper->SetScalarModeToUseCellFieldData();
     402        if (_bondMapper != NULL)
     403            _bondMapper->SetScalarModeToUseCellFieldData();
     404        break;
     405    default:
     406        ERROR("Unsupported DataAttributeType: %d", type);
     407        return;
     408    }
     409
     410    if (name != NULL && strlen(name) > 0) {
     411        if (_atomMapper != NULL)
     412            _atomMapper->SelectColorArray(name);
     413        if (_bondMapper != NULL)
     414            _bondMapper->SelectColorArray(name);
     415    } else {
     416        if (_atomMapper != NULL)
     417            _atomMapper->SetScalarModeToDefault();
     418        if (_bondMapper != NULL)
     419            _bondMapper->SetScalarModeToDefault();
     420    }
     421
     422    if (_lut != NULL) {
     423        if (range != NULL) {
     424            _lut->SetRange(range);
     425        } else if (name != NULL && strlen(name) > 0) {
     426            double r[2];
     427            int comp = -1;
     428            if (mode == COLOR_BY_VECTOR_X)
     429                comp = 0;
     430            else if (mode == COLOR_BY_VECTOR_Y)
     431                comp = 1;
     432            else if (mode == COLOR_BY_VECTOR_Z)
     433                comp = 2;
     434
     435            if (_renderer->getUseCumulativeRange()) {
     436                int numComponents;
     437                if  (!_dataSet->getFieldInfo(name, type, &numComponents)) {
     438                    ERROR("Field not found: %s, type: %d", name, type);
     439                    return;
     440                } else if (numComponents < comp+1) {
     441                    ERROR("Request for component %d in field with %d components",
     442                          comp, numComponents);
     443                    return;
     444                }
     445                _renderer->getCumulativeDataRange(r, name, type, numComponents, comp);
     446            } else {
     447                _dataSet->getDataRange(r, name, type, comp);
     448            }
     449            _lut->SetRange(r);
     450        }
     451    }
     452
     453    switch (mode) {
     454    case COLOR_BY_ELEMENTS:
     455        if (_atomMapper != NULL)
     456            _atomMapper->ScalarVisibilityOn();
     457        if (_bondMapper != NULL)
     458            _bondMapper->ScalarVisibilityOn();
     459        break;
     460    case COLOR_BY_SCALAR:
     461        if (_atomMapper != NULL)
     462            _atomMapper->ScalarVisibilityOn();
     463        if (_bondMapper != NULL)
     464            _bondMapper->ScalarVisibilityOn();
     465        break;
     466    case COLOR_BY_VECTOR_MAGNITUDE:
     467        if (_atomMapper != NULL)
     468            _atomMapper->ScalarVisibilityOn();
     469        if (_bondMapper != NULL)
     470            _bondMapper->ScalarVisibilityOn();
     471        if (_lut != NULL) {
     472            _lut->SetVectorModeToMagnitude();
     473        }
     474        break;
     475    case COLOR_BY_VECTOR_X:
     476        if (_atomMapper != NULL)
     477            _atomMapper->ScalarVisibilityOn();
     478        if (_bondMapper != NULL)
     479            _bondMapper->ScalarVisibilityOn();
     480        if (_lut != NULL) {
     481            _lut->SetVectorModeToComponent();
     482            _lut->SetVectorComponent(0);
     483        }
     484        break;
     485    case COLOR_BY_VECTOR_Y:
     486        if (_atomMapper != NULL)
     487            _atomMapper->ScalarVisibilityOn();
     488        if (_bondMapper != NULL)
     489            _bondMapper->ScalarVisibilityOn();
     490        if (_lut != NULL) {
     491            _lut->SetVectorModeToComponent();
     492            _lut->SetVectorComponent(1);
     493        }
     494        break;
     495    case COLOR_BY_VECTOR_Z:
     496        if (_atomMapper != NULL)
     497            _atomMapper->ScalarVisibilityOn();
     498        if (_bondMapper != NULL)
     499            _bondMapper->ScalarVisibilityOn();
     500        if (_lut != NULL) {
     501            _lut->SetVectorModeToComponent();
     502            _lut->SetVectorComponent(2);
     503        }
     504        break;
     505    case COLOR_CONSTANT:
     506    default:
     507        if (_atomMapper != NULL)
     508            _atomMapper->ScalarVisibilityOff();
     509        if (_bondMapper != NULL)
     510            _bondMapper->ScalarVisibilityOff();
     511        break;
     512    }
     513}
     514
    220515/**
    221516 * \brief Called when the color map has been edited
     
    375670        vtkDataSet *ds = _dataSet->getVtkDataSet();
    376671        addRadiusArray(ds, _atomScaling, _radiusScale);
     672#ifndef MOLECULE_USE_GLYPH3D_MAPPER
    377673        if (_glypher != NULL) {
    378674            assert(ds->GetPointData() != NULL &&
     
    381677            _glypher->ScalingOn();
    382678        }
     679#else
     680        if (_atomMapper != NULL) {
     681             assert(ds->GetPointData() != NULL &&
     682                   ds->GetPointData()->GetVectors() != NULL);
     683            _atomMapper->SetScaleModeToScaleByMagnitude();
     684            _atomMapper->SetScaleArray(vtkDataSetAttributes::VECTORS);
     685            _atomMapper->ScalingOn();
     686        }
     687#endif
    383688    }
    384689}
     
    402707void Molecule::setBondRadiusScale(double scale)
    403708{
     709#ifndef MOLECULE_USE_GLYPH3D_MAPPER
    404710    if (_tuber != NULL)
    405711        _tuber->SetRadius(scale);
    406 }
     712#else
     713    if (_cylinderSource != NULL)
     714        _cylinderSource->SetRadius(scale);
     715#endif
     716}
     717
     718#ifdef MOLECULE_USE_GLYPH3D_MAPPER
     719void Molecule::setupBondPolyData()
     720{
     721    if (_dataSet == NULL)
     722        return;
     723    if (_bondPD == NULL) {
     724        _bondPD = vtkSmartPointer<vtkPolyData>::New();
     725    } else {
     726        _bondPD->Initialize();
     727    }
     728    vtkPolyData *pd = vtkPolyData::SafeDownCast(_dataSet->getVtkDataSet());
     729    if (pd == NULL)
     730        return;
     731    vtkCellArray *lines = pd->GetLines();
     732    lines->InitTraversal();
     733    vtkIdType npts, *pts;
     734    vtkSmartPointer<vtkPoints> bondPoints = vtkSmartPointer<vtkPoints>::New();
     735    vtkSmartPointer<vtkIntArray> bondElements = vtkSmartPointer<vtkIntArray>::New();
     736    vtkSmartPointer<vtkDoubleArray> bondVectors = vtkSmartPointer<vtkDoubleArray>::New();
     737    bondElements->SetName("element");
     738    bondVectors->SetName("bonds");
     739    bondVectors->SetNumberOfComponents(3);
     740    vtkDataArray *elements = NULL;
     741    if (pd->GetPointData() != NULL &&
     742        pd->GetPointData()->GetScalars() != NULL &&
     743        strcmp(pd->GetPointData()->GetScalars()->GetName(), "element") == 0) {
     744        elements = pd->GetPointData()->GetScalars();
     745    }
     746    for (int i = 0; lines->GetNextCell(npts, pts); i++) {
     747        assert(npts == 2);
     748        double pt0[3], pt1[3];
     749        double newPt0[3], newPt1[3];
     750        //double *pt0 = pd->GetPoint(pts[0]);
     751        //double *pt1 = pd->GetPoint(pts[1]);
     752        pd->GetPoint(pts[0], pt0);
     753        pd->GetPoint(pts[1], pt1);
     754        double center[3];
     755
     756        for (int j = 0; j < 3; j++)
     757            center[j] = pt0[j] + (pt1[j] - pt0[j]) * 0.5;
     758        for (int j = 0; j < 3; j++)
     759            newPt0[j] = pt0[j] + (pt1[j] - pt0[j]) * 0.25;
     760        for (int j = 0; j < 3; j++)
     761            newPt1[j] = pt0[j] + (pt1[j] - pt0[j]) * 0.75;
     762
     763        bondPoints->InsertNextPoint(newPt0);
     764        bondPoints->InsertNextPoint(newPt1);
     765
     766        TRACE("Bond %d: (%g,%g,%g)-(%g,%g,%g)-(%g,%g,%g)", i,
     767              pt0[0], pt0[1], pt0[2],
     768              center[0], center[1], center[2],
     769              pt1[0], pt1[1], pt1[2]);
     770
     771        double vec[3];
     772        for (int j = 0; j < 3; j++)
     773            vec[j] = center[j] - pt0[j];
     774
     775        bondVectors->InsertNextTupleValue(vec);
     776        bondVectors->InsertNextTupleValue(vec);
     777        TRACE("Bond %d, vec: %g,%g,%g", i, vec[0], vec[1], vec[2]);
     778
     779        if (elements != NULL) {
     780            int element = (int)elements->GetComponent(pts[0], 0);
     781            TRACE("Bond %d, elt 0: %d", i, element);
     782            bondElements->InsertNextValue(element);
     783            element = (int)elements->GetComponent(pts[1], 0);
     784            TRACE("Bond %d, elt 1: %d", i, element);
     785            bondElements->InsertNextValue(element);
     786        }
     787    }
     788    _bondPD->SetPoints(bondPoints);
     789    _bondPD->GetPointData()->SetScalars(bondElements);
     790    _bondPD->GetPointData()->SetVectors(bondVectors);
     791}
     792#endif
    407793
    408794void Molecule::addLabelArray(vtkDataSet *dataSet)
  • trunk/packages/vizservers/vtkvis/RpMolecule.h

    r3131 r3141  
    1515#include <vtkActor2D.h>
    1616#include <vtkAssembly.h>
     17#if ((VTK_MAJOR_VERSION > 5) || (VTK_MAJOR_VERSION == 5 && VTK_MINOR_VERSION >= 8))
     18#define MOLECULE_USE_GLYPH3D_MAPPER
     19#include <vtkGlyph3DMapper.h>
     20#include <vtkCylinderSource.h>
     21#else
     22#include <vtkGlyph3D.h>
    1723#include <vtkTubeFilter.h>
    18 #include <vtkGlyph3D.h>
     24#endif
    1925#include <vtkLabelPlacementMapper.h>
    2026
     
    4450        ATOMIC_RADIUS
    4551    };
     52    enum ColorMode {
     53        COLOR_BY_ELEMENTS,
     54        COLOR_BY_SCALAR,
     55        COLOR_BY_VECTOR_MAGNITUDE,
     56        COLOR_BY_VECTOR_X,
     57        COLOR_BY_VECTOR_Y,
     58        COLOR_BY_VECTOR_Z,
     59        COLOR_CONSTANT
     60    };
    4661    enum BondColorMode {
    4762        BOND_COLOR_BY_ELEMENTS,
     
    6277    }
    6378
     79    virtual void setDataSet(DataSet *dataSet,
     80                            Renderer *renderer);
     81
    6482    virtual void setClippingPlanes(vtkPlaneCollection *planes);
     83
     84    void setColorMode(ColorMode mode, DataSet::DataAttributeType type,
     85                      const char *name, double range[2] = NULL);
     86
     87    void setColorMode(ColorMode mode,
     88                      const char *name, double range[2] = NULL);
     89
     90    void setColorMode(ColorMode mode);
    6591
    6692    void setColorMap(ColorMap *colorMap);
     
    104130    virtual void update();
    105131
     132#ifdef MOLECULE_USE_GLYPH3D_MAPPER
     133    void setupBondPolyData();
     134#endif
     135
    106136    static void addLabelArray(vtkDataSet *dataSet);
    107137
     
    111141    double _radiusScale;
    112142    AtomScaling _atomScaling;
     143    bool _labelsOn;
    113144    ColorMap *_colorMap;
    114     bool _labelsOn;
     145    ColorMode _colorMode;
     146    std::string _colorFieldName;
     147    DataSet::DataAttributeType _colorFieldType;
     148    double _colorFieldRange[2];
     149    double _vectorMagnitudeRange[2];
     150    double _vectorComponentRange[3][2];
     151    Renderer *_renderer;
    115152
    116153    vtkSmartPointer<vtkLookupTable> _lut;
     
    118155    vtkSmartPointer<vtkActor> _bondProp;
    119156    vtkSmartPointer<vtkActor2D> _labelProp;
     157#ifndef MOLECULE_USE_GLYPH3D_MAPPER
    120158    vtkSmartPointer<vtkGlyph3D> _glypher;
    121159    vtkSmartPointer<vtkTubeFilter> _tuber;
    122160    vtkSmartPointer<vtkPolyDataMapper> _atomMapper;
    123161    vtkSmartPointer<vtkPolyDataMapper> _bondMapper;
     162#else
     163    vtkSmartPointer<vtkPolyData> _bondPD;
     164    vtkSmartPointer<vtkCylinderSource> _cylinderSource;
     165    vtkSmartPointer<vtkGlyph3DMapper> _atomMapper;
     166    vtkSmartPointer<vtkGlyph3DMapper> _bondMapper;
     167#endif
    124168    vtkSmartPointer<vtkLabelPlacementMapper> _labelMapper;
    125169};
Note: See TracChangeset for help on using the changeset viewer.