Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VtkUnstructuredGridPipeline Class Reference

#include <G4VtkUnstructuredGridPipeline.hh>

+ Inheritance diagram for G4VtkUnstructuredGridPipeline:

Classes

class  PseudoSceneForCubicalCells
 
class  PseudoSceneForTetCells
 
class  PseudoSceneVtkBase
 

Public Member Functions

 G4VtkUnstructuredGridPipeline (G4String name, const G4VtkVisContext &vc)
 
 ~G4VtkUnstructuredGridPipeline ()=default
 
void AddFilter (vtkSmartPointer< vtkUnstructuredGridAlgorithm > f)
 
vtkSmartPointer< vtkUnstructuredGridAlgorithm > GetFilter (G4int iFilter)
 
std::size_t GetNumberOfFilters ()
 
vtkSmartPointer< vtkUnstructuredGridAlgorithm > GetFinalFilter ()
 
void Modified ()
 
void Clear ()
 
void Print ()
 
void Enable ()
 
void Disable ()
 
void SetUnstructuredGridData (const G4Mesh &mesh)
 
- Public Member Functions inherited from G4VVtkPipeline
 G4VVtkPipeline ()
 
 G4VVtkPipeline (G4String nameIn, G4String typeIn, const G4VtkVisContext &vcIn, G4bool disableParentIn=false, vtkSmartPointer< vtkRenderer > rendererIn=nullptr)
 
virtual ~G4VVtkPipeline ()
 
void SetDisableParent (G4bool disableParentIn)
 
bool GetDisableParent ()
 
void SetName (G4String nameIn)
 
const G4String GetName ()
 
void SetTypeName (G4String typeNameIn)
 
G4String GetTypeName ()
 
G4VtkVisContextGetVtkVisContext ()
 
void AddChildPipeline (G4VVtkPipeline *child)
 
G4VVtkPipelineGetChildPipeline (G4int iPipeline)
 
G4int GetNumberOfChildPipelines ()
 
std::vector< G4VVtkPipeline * > GetChildPipelines ()
 
void ClearChildPipeline ()
 

Additional Inherited Members

- Protected Attributes inherited from G4VVtkPipeline
G4String name
 
G4String type
 
G4bool disableParent
 
std::vector< G4VVtkPipeline * > childPipelines
 
vtkSmartPointer< vtkRenderer > renderer
 
G4VtkVisContext vc
 

Detailed Description

Definition at line 56 of file G4VtkUnstructuredGridPipeline.hh.

Constructor & Destructor Documentation

◆ G4VtkUnstructuredGridPipeline()

G4VtkUnstructuredGridPipeline::G4VtkUnstructuredGridPipeline ( G4String name,
const G4VtkVisContext & vc )

Definition at line 63 of file G4VtkUnstructuredGridPipeline.cc.

63 :
64 G4VVtkPipeline(nameIn, "G4VtkUnstructuredPipeline", vcIn, false, vcIn.fViewer->renderer) {
65
67
68 pointColourValues = vtkSmartPointer<vtkDoubleArray>::New();
69 cellColourValues = vtkSmartPointer<vtkDoubleArray>::New();
70 pointColourValues->SetNumberOfComponents(4);
71 cellColourValues->SetNumberOfComponents(4);
72
73 pointColourIndices = vtkSmartPointer<vtkDoubleArray>::New();
74 cellColourIndices = vtkSmartPointer<vtkDoubleArray>::New();
75 pointColourIndices->SetNumberOfComponents(1);
76 cellColourIndices->SetNumberOfComponents(1);
77
79 colourLUT->DiscretizeOn();
80
82 unstructuredGrid->SetPoints(points);
83 unstructuredGrid->GetPointData()->SetScalars(pointColourValues);
84 unstructuredGrid->GetCellData()->SetScalars(cellColourValues);
85
86 // clean filter
87#if 0
89 clean->SetInputData(unstructuredGrid);
90 clean->ToleranceIsAbsoluteOff();
91 clean->SetTolerance(1e-6);
92#endif
93
94 // Clip filter
96 vtkNew<vtkPlane> plane;
97 clip->SetClipFunction(plane);
98 clip->SetInputData(unstructuredGrid);
99
100 // Triangle filter
102 tri->SetInputData(unstructuredGrid);
103
104 // create dataset mapper
106 mapper->SetScalarModeToUseCellData();
107 mapper->SetColorModeToDirectScalars();
108 mapper->SetInputData(unstructuredGrid);
109 //mapper->SetInputConnection(clip->GetOutputPort());
110 //mapper->SetInputConnection(tri->GetOutputPort());
111
112 // create volume mapper
114 volumeMapper->SetScalarModeToUseCellData();
115 volumeMapper->SetInputConnection(tri->GetOutputPort());
116
117 // create volume properties
118 auto volumeProp = vtkSmartPointer<vtkVolumeProperty>::New();
119 volumeProp->SetColor(colourLUT);
120
121 // create actor
123 actor->SetMapper(mapper);
124 actor->SetVisibility(1);
125
126 // create volume
128 volume->SetMapper(volumeMapper);
129 //volume->SetProperty(volumeProp);
130 volume->SetVisibility(1);
131
132 // add to renderer
133 vc.fViewer->renderer->AddActor(actor);
134 //vc.fViewer->renderer->AddVolume(volume);
135}
G4VtkVisContext vc
vtkNew< vtkRenderer > renderer
const G4VtkViewer * fViewer

◆ ~G4VtkUnstructuredGridPipeline()

G4VtkUnstructuredGridPipeline::~G4VtkUnstructuredGridPipeline ( )
default

Member Function Documentation

◆ AddFilter()

void G4VtkUnstructuredGridPipeline::AddFilter ( vtkSmartPointer< vtkUnstructuredGridAlgorithm > f)
inline

Definition at line 63 of file G4VtkUnstructuredGridPipeline.hh.

63{ filters.push_back(f); }

◆ Clear()

void G4VtkUnstructuredGridPipeline::Clear ( )
inlinevirtual

Reimplemented from G4VVtkPipeline.

Definition at line 69 of file G4VtkUnstructuredGridPipeline.hh.

69{};

◆ Disable()

void G4VtkUnstructuredGridPipeline::Disable ( )
inlinevirtual

Implements G4VVtkPipeline.

Definition at line 73 of file G4VtkUnstructuredGridPipeline.hh.

73{};

◆ Enable()

void G4VtkUnstructuredGridPipeline::Enable ( )
inlinevirtual

Implements G4VVtkPipeline.

Definition at line 72 of file G4VtkUnstructuredGridPipeline.hh.

72{};

◆ GetFilter()

vtkSmartPointer< vtkUnstructuredGridAlgorithm > G4VtkUnstructuredGridPipeline::GetFilter ( G4int iFilter)
inline

Definition at line 64 of file G4VtkUnstructuredGridPipeline.hh.

64{ return filters[iFilter]; }

◆ GetFinalFilter()

vtkSmartPointer< vtkUnstructuredGridAlgorithm > G4VtkUnstructuredGridPipeline::GetFinalFilter ( )
inline

Definition at line 66 of file G4VtkUnstructuredGridPipeline.hh.

66{ return filters[filters.size() - 1]; }

◆ GetNumberOfFilters()

std::size_t G4VtkUnstructuredGridPipeline::GetNumberOfFilters ( )
inline

Definition at line 65 of file G4VtkUnstructuredGridPipeline.hh.

65{ return filters.size(); }

◆ Modified()

void G4VtkUnstructuredGridPipeline::Modified ( )
inlinevirtual

Reimplemented from G4VVtkPipeline.

Definition at line 68 of file G4VtkUnstructuredGridPipeline.hh.

68{};

◆ Print()

void G4VtkUnstructuredGridPipeline::Print ( )
inlinevirtual

Reimplemented from G4VVtkPipeline.

Definition at line 70 of file G4VtkUnstructuredGridPipeline.hh.

70{};

◆ SetUnstructuredGridData()

void G4VtkUnstructuredGridPipeline::SetUnstructuredGridData ( const G4Mesh & mesh)

Definition at line 264 of file G4VtkUnstructuredGridPipeline.cc.

264 {
265
266 // Modelling parameters
268 tmpMP.SetCulling(true); // This avoids drawing transparent...
269 tmpMP.SetCullingInvisible(true); // ... or invisble volumes.
270
271 // Physical volume model
272 const auto& container = mesh.GetContainerVolume();
273 const G4bool useFullExtent = true; // To avoid calculating the extent
274 G4PhysicalVolumeModel tmpPVModel(container,
276 G4Transform3D(), // so that positions are in local coordinates
277 &tmpMP,
278 useFullExtent);
279
280 // Graphics scene
281 if(mesh.GetMeshType() == G4Mesh::tetrahedron) {
282 PseudoSceneForTetCells pseudoScene(&tmpPVModel, mesh.GetMeshDepth(), points,
283 pointColourValues, cellColourValues,
284 pointColourIndices, cellColourIndices,
285 colourLUT, unstructuredGrid);
286 tmpPVModel.DescribeYourselfTo(pseudoScene);
287 //G4cout << pseudoScene.GetNumberOfPoints() << " " << pseudoScene.GetNumberOfCells() << " " << pseudoScene.GetNumberOfAddedCells() << G4endl;
288 //pseudoScene.DumpColourMap();
289 }
290 else if(mesh.GetMeshType() == G4Mesh::nested3DRectangular) {
291 PseudoSceneForCubicalCells pseudoScene(&tmpPVModel, mesh.GetMeshDepth(), points,
292 pointColourValues, cellColourValues,
293 pointColourIndices, cellColourIndices,
294 colourLUT, unstructuredGrid);
295 tmpPVModel.DescribeYourselfTo(pseudoScene);
296 }
297 else if(mesh.GetMeshType() == G4Mesh::rectangle) {
298 PseudoSceneForCubicalCells pseudoScene(&tmpPVModel, mesh.GetMeshDepth(), points,
299 pointColourValues, cellColourValues,
300 pointColourIndices, cellColourIndices,
301 colourLUT, unstructuredGrid);
302 tmpPVModel.DescribeYourselfTo(pseudoScene);
303 }
304
305}
HepGeom::Transform3D G4Transform3D
bool G4bool
Definition G4Types.hh:86
MeshType GetMeshType() const
Definition G4Mesh.hh:75
G4VPhysicalVolume * GetContainerVolume() const
Definition G4Mesh.hh:73
@ rectangle
Definition G4Mesh.hh:53
@ nested3DRectangular
Definition G4Mesh.hh:54
@ tetrahedron
Definition G4Mesh.hh:57
G4int GetMeshDepth() const
Definition G4Mesh.hh:76
void SetCulling(G4bool)
void SetCullingInvisible(G4bool)

The documentation for this class was generated from the following files: