Loading...
Searching...
No Matches
WLSTrajectoryPoint.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/WLSTrajectoryPoint.cc
28/// \brief Implementation of the WLSTrajectoryPoint class
29//
30//
31#include "WLSTrajectoryPoint.hh"
32
33#include "G4AttDef.hh"
34#include "G4AttDefStore.hh"
35#include "G4AttValue.hh"
36#include "G4Step.hh"
37#include "G4StepStatus.hh"
38#include "G4Track.hh"
39#include "G4UIcommand.hh"
40#include "G4UnitsTable.hh"
41#include "G4VProcess.hh"
42
44
45//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
46
48 : G4TrajectoryPoint(aStep->GetPostStepPoint()->GetPosition())
49{
50 auto postStepPoint = aStep->GetPostStepPoint();
51 fTime = postStepPoint->GetGlobalTime();
52 fMomentum = postStepPoint->GetMomentum();
53 fStepStatus = postStepPoint->GetStepStatus();
54 if(postStepPoint->GetPhysicalVolume())
55 {
56 fVolumeName = postStepPoint->GetPhysicalVolume()->GetName();
57 }
58 else
59 {
60 fVolumeName = G4String(" ");
61 }
62}
63
64//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
65
67 : G4TrajectoryPoint(aTrack->GetPosition())
68{
69 fTime = aTrack->GetGlobalTime();
70 fMomentum = aTrack->GetMomentum();
71 fStepStatus = fUndefined;
72 fVolumeName = aTrack->GetVolume()->GetName();
73}
74
75//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
76
79{
80 fTime = right.fTime;
81 fMomentum = right.fMomentum;
82 fStepStatus = right.fStepStatus;
83 fVolumeName = right.fVolumeName;
84}
85
86//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
87
88const std::map<G4String, G4AttDef>* WLSTrajectoryPoint::GetAttDefs() const
89{
90 G4bool isNew;
91 std::map<G4String, G4AttDef>* store =
92 G4AttDefStore::GetInstance("TrajectoryPoint", isNew);
93 if(isNew)
94 {
95 G4String Pos("Pos");
96 (*store)[Pos] =
97 G4AttDef(Pos, "Position", "Physics", "G4BestUnit", "G4ThreeVector");
98
99 G4String Time("Time");
100 (*store)[Time] =
101 G4AttDef(Time, "Time", "Physics", "G4BestUnit", "G4double");
102
103 G4String Momentum("Momentum");
104 (*store)[Momentum] =
105 G4AttDef(Momentum, "Momentum", "Physics", "G4BestUnit", "G4ThreeVector");
106
107 G4String StepStatus("StepStatus");
108 (*store)[StepStatus] =
109 G4AttDef(StepStatus, "StepStatus", "Physics", "", "G4StepStatus");
110
111 G4String VolumeName("VolumeName");
112 (*store)[VolumeName] =
113 G4AttDef(VolumeName, "VolumeName", "Physics", "", "G4String");
114 }
115 return store;
116}
117
118//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
119
120std::vector<G4AttValue>* WLSTrajectoryPoint::CreateAttValues() const
121{
122 auto values = new std::vector<G4AttValue>;
123
124 values->push_back(G4AttValue("Time", G4BestUnit(fTime, "Time"), ""));
125 values->push_back(
126 G4AttValue("Momentum", G4BestUnit(fMomentum, "Momentum"), ""));
127 values->push_back(G4AttValue("StepStatus", G4UIcommand::ConvertToString(fStepStatus), ""));
128 values->push_back(G4AttValue("VolumeName", fVolumeName, ""));
129
130 return values;
131}
G4ThreadLocal G4Allocator< WLSTrajectoryPoint > * WLSTrajPointAllocator
Definition of the WLSTrajectoryPoint class.
const std::map< G4String, G4AttDef > * GetAttDefs() const override
WLSTrajectoryPoint()=default
std::vector< G4AttValue > * CreateAttValues() const override

Applications | User Support | Publications | Collaboration