Loading...
Searching...
No Matches
WLSTrajectory.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//
27/// \file optical/wls/src/WLSTrajectory.cc
28/// \brief Implementation of the WLSTrajectory class
29//
30//
31
32#include "WLSTrajectory.hh"
33
34#include "WLSTrajectoryPoint.hh"
35
36#include "G4AttDef.hh"
37#include "G4AttDefStore.hh"
38#include "G4AttValue.hh"
39#include "G4Colour.hh"
40#include "G4ParticleTable.hh"
41#include "G4ParticleTypes.hh"
42#include "G4Polyline.hh"
43#include "G4Polymarker.hh"
44#include "G4UIcommand.hh"
45#include "G4UnitsTable.hh"
46#include "G4VisAttributes.hh"
47#include "G4VVisManager.hh"
48
50
51//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
52
54{
55 fParticleDefinition = aTrack->GetDefinition();
56 fParticleName = fParticleDefinition->GetParticleName();
57 fPDGCharge = fParticleDefinition->GetPDGCharge();
58 fPDGEncoding = fParticleDefinition->GetPDGEncoding();
59 fTrackID = aTrack->GetTrackID();
60 fParentID = aTrack->GetParentID();
61 fInitialMomentum = aTrack->GetMomentum();
63 // Following is for the first trajectory point
64 fpPointsContainer->push_back(new WLSTrajectoryPoint(aTrack));
65}
66
67//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
68
71{
72 fParticleDefinition = right.fParticleDefinition;
73 fParticleName = right.fParticleName;
74 fPDGCharge = right.fPDGCharge;
75 fPDGEncoding = right.fPDGEncoding;
76 fTrackID = right.fTrackID;
77 fParentID = right.fParentID;
78 fInitialMomentum = right.fInitialMomentum;
80
81 for(size_t i = 0; i < right.fpPointsContainer->size(); ++i)
82 {
83 auto rightPoint = (WLSTrajectoryPoint*) ((*(right.fpPointsContainer))[i]);
84 fpPointsContainer->push_back(new WLSTrajectoryPoint(*rightPoint));
85 }
86}
87
88//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
89
91{
92 for(size_t i = 0; i < fpPointsContainer->size(); ++i)
93 {
94 delete(*fpPointsContainer)[i];
95 }
96 fpPointsContainer->clear();
97 delete fpPointsContainer;
98}
99
100//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
101
102void WLSTrajectory::ShowTrajectory(std::ostream& os) const
103{
104 // Invoke the default implementation in G4VTrajectory...
105 G4VTrajectory::ShowTrajectory(os);
106 // ... or override with your own code here.
107}
108
109//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
110
112{
113 fpPointsContainer->push_back(new WLSTrajectoryPoint(aStep));
114}
115
116//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
117
119{
120 return (G4ParticleTable::GetParticleTable()->FindParticle(fParticleName));
121}
122
123//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
124
126{
127 if(!secondTrajectory)
128 return;
129
130 auto second = (WLSTrajectory*) secondTrajectory;
131 G4int ent = second->GetPointEntries();
132 // initial point of the second trajectory should not be merged
133 for(G4int i = 1; i < ent; ++i)
134 {
135 fpPointsContainer->push_back((*(second->fpPointsContainer))[i]);
136 }
137 delete(*second->fpPointsContainer)[0];
138 second->fpPointsContainer->clear();
139}
140
141//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
142
143const std::map<G4String, G4AttDef>* WLSTrajectory::GetAttDefs() const
144{
145 G4bool isNew;
146 std::map<G4String, G4AttDef>* store =
147 G4AttDefStore::GetInstance("Trajectory", isNew);
148
149 if(isNew)
150 {
151 G4String ID("ID");
152 (*store)[ID] = G4AttDef(ID, "Track ID", "Bookkeeping", "", "G4int");
153
154 G4String PID("PID");
155 (*store)[PID] = G4AttDef(PID, "Parent ID", "Bookkeeping", "", "G4int");
156
157 G4String PN("PN");
158 (*store)[PN] = G4AttDef(PN, "Particle Name", "Physics", "", "G4String");
159
160 G4String Ch("Ch");
161 (*store)[Ch] = G4AttDef(Ch, "Charge", "Physics", "e+", "G4double");
162
163 G4String PDG("PDG");
164 (*store)[PDG] = G4AttDef(PDG, "PDG Encoding", "Physics", "", "G4int");
165
166 G4String IMom("IMom");
167 (*store)[IMom] = G4AttDef(IMom, "Momentum of track at start of trajectory",
168 "Physics", "G4BestUnit", "G4ThreeVector");
169
170 G4String IMag("IMag");
171 (*store)[IMag] =
172 G4AttDef(IMag, "Magnitude of momentum of track at start of trajectory",
173 "Physics", "G4BestUnit", "G4double");
174
175 G4String NTP("NTP");
176 (*store)[NTP] = G4AttDef(NTP, "No. of points", "Bookkeeping", "", "G4int");
177 }
178 return store;
179}
180
181//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
182
183std::vector<G4AttValue>* WLSTrajectory::CreateAttValues() const
184{
185 auto values = new std::vector<G4AttValue>;
186
187 values->push_back(
188 G4AttValue("ID", G4UIcommand::ConvertToString(fTrackID), ""));
189
190 values->push_back(
191 G4AttValue("PID", G4UIcommand::ConvertToString(fParentID), ""));
192
193 values->push_back(G4AttValue("PN", fParticleName, ""));
194
195 values->push_back(
196 G4AttValue("Ch", G4UIcommand::ConvertToString(fPDGCharge), ""));
197
198 values->push_back(
199 G4AttValue("PDG", G4UIcommand::ConvertToString(fPDGEncoding), ""));
200
201 values->push_back(
202 G4AttValue("IMom", G4BestUnit(fInitialMomentum, "Energy"), ""));
203
204 values->push_back(
205 G4AttValue("IMag", G4BestUnit(fInitialMomentum.mag(), "Energy"), ""));
206
207 values->push_back(
208 G4AttValue("NTP", G4UIcommand::ConvertToString(GetPointEntries()), ""));
209
210 return values;
211}
Definition of the WLSTrajectoryPoint class.
G4ThreadLocal G4Allocator< WLSTrajectory > * WLSTrajectoryAllocator
Definition of the WLSTrajectory class.
std::vector< G4VTrajectoryPoint * > WLSTrajectoryPointContainer
int GetPointEntries() const override
void ShowTrajectory(std::ostream &os=G4cout) const override
G4ThreeVector fInitialMomentum
WLSTrajectoryPointContainer * fpPointsContainer
std::vector< G4AttValue > * CreateAttValues() const override
G4ParticleDefinition * GetParticleDefinition()
G4ParticleDefinition * fParticleDefinition
void MergeTrajectory(G4VTrajectory *secondTrajectory) override
~WLSTrajectory() override
G4String fParticleName
G4double fPDGCharge
const std::map< G4String, G4AttDef > * GetAttDefs() const override
void AppendStep(const G4Step *aStep) override
WLSTrajectory()=default

Applications | User Support | Publications | Collaboration