Geant4 Cross Reference |
1 // 2 // ******************************************************************** 3 // * License and Disclaimer * 4 // * * 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. * 10 // * * 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitation of liability. * 17 // * * 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************************************** 25 26 #include "G4VtkPolydataInstanceTensorPipeline.hh" 27 28 #include "G4VtkViewer.hh" 29 #include "G4VtkVisContext.hh" 30 31 #include <vtkActor.h> 32 #include <vtkDataArray.h> 33 #include <vtkDoubleArray.h> 34 #include <vtkPointData.h> 35 #include <vtkPolyData.h> 36 #include <vtkPolyDataMapper.h> 37 #include <vtkProperty.h> 38 #include <vtkTensorGlyphColor.h> 39 40 std::size_t G4VtkPolydataInstanceTensorPipeline::MakeHash(const G4Polyhedron& polyhedron, 41 const G4VtkVisContext& vc) 42 { 43 // Get view parameters that the user can force through the vis attributes, thereby over-riding the 44 // current view parameter. 45 G4ViewParameters::DrawingStyle drawing_style = vc.fDrawingStyle; 46 47 std::size_t vhash = 0; 48 std::hash_combine(vhash, static_cast<int>(drawing_style)); 49 50 std::size_t phash = std::hash<G4Polyhedron>{}(polyhedron); 51 52 std::size_t hash = 0; 53 std::hash_combine(hash, phash); 54 std::hash_combine(hash, vhash); 55 56 return hash; 57 } 58 59 G4VtkPolydataInstanceTensorPipeline::G4VtkPolydataInstanceTensorPipeline(G4String nameIn, 60 const G4VtkVisContext& vcIn) 61 : G4VtkPolydataInstancePipeline(nameIn, vcIn) 62 { 63 // Set pipeline type 64 SetTypeName(G4String("G4VtkPolydataInstanceTensorPipeline")); 65 66 // make polydata from instance data 67 instancePolydata = vtkSmartPointer<vtkPolyData>::New(); 68 instancePolydata->SetPoints(instancePosition); 69 instancePolydata->GetPointData()->SetTensors(instanceTransform); 70 instancePolydata->GetPointData()->SetVectors(instanceColour); 71 instancePolydata->GetPointData()->SetScalars(instanceColour); 72 73 // tensor glyph 74 instanceTensorGlyph = vtkSmartPointer<vtkTensorGlyphColor>::New(); 75 instanceTensorGlyph->SetInputData(instancePolydata); 76 instanceTensorGlyph->SetSourceConnection(GetFinalFilter()->GetOutputPort()); 77 instanceTensorGlyph->ColorGlyphsOn(); 78 instanceTensorGlyph->ScalingOff(); 79 instanceTensorGlyph->ThreeGlyphsOff(); 80 instanceTensorGlyph->ExtractEigenvaluesOff(); 81 instanceTensorGlyph->SetColorModeToScalars(); 82 AddFilter(instanceTensorGlyph); 83 84 // set polydata mapper 85 mapper->SetInputConnection(GetFinalFilter()->GetOutputPort()); 86 87 // set actor 88 actor->SetMapper(mapper); 89 actor->SetVisibility(1); 90 91 // shading parameters 92 actor->GetProperty()->SetAmbient(0.2); 93 actor->GetProperty()->SetDiffuse(0.7); 94 actor->GetProperty()->SetSpecular(0.1); 95 actor->GetProperty()->SetSpecularPower(1); 96 97 // add to renderer 98 vc.fViewer->renderer->AddActor(GetActor()); 99 } 100 101 void G4VtkPolydataInstanceTensorPipeline::Print() 102 { 103 G4cout << "G4VtkPolydataInstanceTensorPipeline " << GetName() << G4endl; 104 G4VtkPolydataPipeline::Print(); 105 }