Geant4-11
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: Katsuya.Amako@kek.jp)
31// Makoto Asai (e-mail: asai@slac.stanford.edu)
32// Takashi Sasaki (e-mail: Takashi.Sasaki@kek.jp)
33// --------------------------------------------------------------------
34
35#include "G4VTrajectory.hh"
36#include "G4VTrajectoryPoint.hh"
37#include "G4AttDefStore.hh"
38#include "G4AttDef.hh"
39#include "G4AttValue.hh"
40#include "G4AttCheck.hh"
41#include "G4VVisManager.hh"
42#include "G4VisAttributes.hh"
43#include "G4Polyline.hh"
44#include "G4Polymarker.hh"
45#include "G4Colour.hh"
46
48{
49}
50
52{
53}
54
56{
57 return (this==&right);
58}
59
60void G4VTrajectory::ShowTrajectory(std::ostream& os) const
61{
62 // Makes use of attribute values implemented in the concrete class.
63 // Note: the user needs to follow with new-line or end-of-string,
64 // depending on the nature of os
65
66 std::vector<G4AttValue>* attValues = CreateAttValues();
67 const std::map<G4String,G4AttDef>* attDefs = GetAttDefs();
68
69 // Ensure validity...
70 //
71 if (G4AttCheck(attValues,attDefs).Check("G4VTrajectory::ShowTrajectory"))
72 {
73 return;
74 }
75
76 os << "Trajectory:";
77
78 for (auto iAttVal = attValues->cbegin();
79 iAttVal != attValues->cend(); ++iAttVal)
80 {
81 std::map<G4String,G4AttDef>::const_iterator iAttDef =
82 attDefs->find(iAttVal->GetName());
83 os << "\n " << iAttDef->second.GetDesc()
84 << " (" << iAttVal->GetName()
85 << "): " << iAttVal->GetValue();
86 }
87
88 delete attValues; // AttValues must be deleted after use.
89
90 // Now do trajectory points...
91
92 for (G4int i=0; i<GetPointEntries(); ++i)
93 {
94 G4VTrajectoryPoint* aTrajectoryPoint = GetPoint(i);
95 attValues = aTrajectoryPoint->CreateAttValues();
96 attDefs = aTrajectoryPoint->GetAttDefs();
97
98 // Ensure validity...
99 //
100 if (G4AttCheck(attValues,attDefs).Check("G4VTrajectory::ShowTrajectory"))
101 {
102 return;
103 }
104
105 for (auto iAttVal = attValues->cbegin();
106 iAttVal != attValues->cend(); ++iAttVal)
107 {
108 std::map<G4String,G4AttDef>::const_iterator iAttDef =
109 attDefs->find(iAttVal->GetName());
110 os << "\n " << iAttDef->second.GetDesc()
111 << " (" << iAttVal->GetName()
112 << "): " << iAttVal->GetValue();
113 }
114
115 delete attValues; // AttValues must be deleted after use.
116 }
117 os << std::endl;
118}
119
121{
123
124 if (pVVisManager != nullptr)
125 {
126 pVVisManager->DispatchToModel(*this);
127 }
128}
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 ~G4VTrajectory()
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual void DispatchToModel(const G4VTrajectory &)=0
static G4VVisManager * GetConcreteInstance()