Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VTrajectory.cc
Go to the documentation of this file.
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// G4VTrajectory class implementation
27//
28// Contact:
29// Questions and comments to this code should be sent to
30// Katsuya Amako (e-mail: [email protected])
31// Makoto Asai (e-mail: [email protected])
32// Takashi Sasaki (e-mail: [email protected])
33// --------------------------------------------------------------------
34
35#include "G4VTrajectory.hh"
36
37#include "G4AttCheck.hh"
38#include "G4AttDef.hh"
39#include "G4AttDefStore.hh"
40#include "G4AttValue.hh"
41#include "G4Colour.hh"
42#include "G4Polyline.hh"
43#include "G4Polymarker.hh"
44#include "G4VTrajectoryPoint.hh"
45#include "G4VVisManager.hh"
46#include "G4VisAttributes.hh"
47
48G4bool G4VTrajectory::operator==(const G4VTrajectory& right) const { return (this == &right); }
49
50void G4VTrajectory::ShowTrajectory(std::ostream& os) const
51{
52 // Makes use of attribute values implemented in the concrete class.
53 // Note: the user needs to follow with new-line or end-of-string,
54 // depending on the nature of os
55
56 std::vector<G4AttValue>* attValues = CreateAttValues();
57 const std::map<G4String, G4AttDef>* attDefs = GetAttDefs();
58
59 // Ensure validity...
60 //
61 if (G4AttCheck(attValues, attDefs).Check("G4VTrajectory::ShowTrajectory")) {
62 return;
63 }
64
65 os << "Trajectory:";
66
67 for (const auto& attValue : *attValues) {
68 auto iAttDef = attDefs->find(attValue.GetName());
69 os << "\n " << iAttDef->second.GetDesc() << " (" << attValue.GetName()
70 << "): " << attValue.GetValue();
71 }
72
73 delete attValues; // AttValues must be deleted after use.
74
75 // Now do trajectory points...
76
77 for (G4int i = 0; i < GetPointEntries(); ++i) {
78 G4VTrajectoryPoint* aTrajectoryPoint = GetPoint(i);
79 attValues = aTrajectoryPoint->CreateAttValues();
80 attDefs = aTrajectoryPoint->GetAttDefs();
81
82 // Ensure validity...
83 //
84 if (G4AttCheck(attValues, attDefs).Check("G4VTrajectory::ShowTrajectory")) {
85 return;
86 }
87
88 for (const auto& attValue : *attValues) {
89 const auto iAttDef = attDefs->find(attValue.GetName());
90 os << "\n " << iAttDef->second.GetDesc() << " (" << attValue.GetName()
91 << "): " << attValue.GetValue();
92 }
93
94 delete attValues; // AttValues must be deleted after use.
95 }
96 os << std::endl;
97}
98
100{
102
103 if (pVVisManager != nullptr) {
104 pVVisManager->DispatchToModel(*this);
105 }
106}
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
virtual std::vector< G4AttValue > * CreateAttValues() const
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual G4VTrajectoryPoint * GetPoint(G4int i) const =0
virtual G4int GetPointEntries() const =0
virtual void ShowTrajectory(std::ostream &os=G4cout) const
virtual std::vector< G4AttValue > * CreateAttValues() const
G4bool operator==(const G4VTrajectory &right) const
virtual void DrawTrajectory() const
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual void DispatchToModel(const G4VTrajectory &)=0
static G4VVisManager * GetConcreteInstance()