Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ClonedRichTrajectory.hh
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// G4ClonedRichTrajectory
27//
28// Class description:
29//
30// This class is identical to G4RichTrajectory class, but uses the
31// singleton G4Allocator so that the master thread can safely
32// delete the object instantiated by worker threads in sub-event
33// parallel mode. This class object should be instantiated as
34// a clone of G4RichTrajectory object.
35//
36// Makoto Asai (JLab) - Oct.2024
37// --------------------------------------------------------------------
38#ifndef G4ClonedRichTrajectory_HH
39#define G4ClonedRichTrajectory_HH 1
40
41#include "G4TouchableHandle.hh"
42#include "G4VTrajectory.hh"
43#include "G4ParticleDefinition.hh" // Include from 'particle+matter'
44#include "G4Step.hh"
45#include "G4Track.hh"
46#include "G4ClonedRichTrajectoryPoint.hh" // Include from 'tracking'
47#include "G4ios.hh" // Include from 'system'
48#include "globals.hh" // Include from 'global'
49
50#include "trkgdefs.hh"
51#include <stdlib.h> // Include from 'system'
52
53#include <vector>
54
55class G4Polyline;
57
59{
60 using G4TrajectoryPointContainer = std::vector<G4VTrajectoryPoint*>;
61
62 public:
63 // Constructors/destructor
64 //
66 ~G4ClonedRichTrajectory() override;
69
70 // Operators
71 //
72 G4bool operator==(const G4ClonedRichTrajectory& r) const { return (this == &r); }
73
74 inline void* operator new(size_t);
75 inline void operator delete(void*);
76
77 // Get/Set functions
78
79 inline G4int GetTrackID() const override { return fTrackID; }
80 inline G4int GetParentID() const override { return fParentID; }
81 inline G4String GetParticleName() const override { return ParticleName; }
82 inline G4double GetCharge() const override { return PDGCharge; }
83 inline G4int GetPDGEncoding() const override { return PDGEncoding; }
84 inline G4double GetInitialKineticEnergy() const { return initialKineticEnergy; }
85 inline G4ThreeVector GetInitialMomentum() const override { return initialMomentum; }
86
87 // Other (virtual) member functions
88 //
89 void ShowTrajectory(std::ostream& os = G4cout) const override;
90 void DrawTrajectory() const override;
91 void AppendStep(const G4Step* aStep) override;
92 void MergeTrajectory(G4VTrajectory* secondTrajectory) override;
93 inline G4int GetPointEntries() const override;
94 inline G4VTrajectoryPoint* GetPoint(G4int i) const override;
95
97
98 // Get methods for HepRep style attributes
99 //
100 const std::map<G4String, G4AttDef>* GetAttDefs() const override;
101 std::vector<G4AttValue>* CreateAttValues() const override;
102
103 private:
104 //G4TrajectoryPointContainer* positionRecord = nullptr;
105 G4int fTrackID = 0;
106 G4int fParentID = 0;
107 G4int PDGEncoding = 0;
108 G4double PDGCharge = 0.0;
109 G4String ParticleName = "dummy";
110 G4double initialKineticEnergy = 0.0;
111 G4ThreeVector initialMomentum;
112
113 // Extended information (only publicly accessible through AttValues)...
114 //
115 G4TrajectoryPointContainer* fpRichPointContainer = nullptr;
116 G4TouchableHandle fpInitialVolume;
117 G4TouchableHandle fpInitialNextVolume;
118 const G4VProcess* fpCreatorProcess = nullptr;
119 G4int fCreatorModelID = 0;
120 G4TouchableHandle fpFinalVolume;
121 G4TouchableHandle fpFinalNextVolume;
122 const G4VProcess* fpEndingProcess = nullptr;
123 G4double fFinalKineticEnergy = 0.0;
124};
125
127
128inline void* G4ClonedRichTrajectory::operator new(size_t)
129{
130 if (aClonedRichTrajectoryAllocator() == nullptr) {
132 }
133 return (void*)aClonedRichTrajectoryAllocator()->MallocSingle();
134}
135
136inline void G4ClonedRichTrajectory::operator delete(void* aRichTrajectory)
137{
138 aClonedRichTrajectoryAllocator()->FreeSingle((G4ClonedRichTrajectory*)aRichTrajectory);
139}
140
142{
143 return G4int(fpRichPointContainer->size());
144}
145
147{
148 return (*fpRichPointContainer)[i];
149}
150
151#endif
G4TRACKING_DLL G4Allocator< G4ClonedRichTrajectory > *& aClonedRichTrajectoryAllocator()
CLHEP::Hep3Vector G4ThreeVector
G4ReferenceCountedHandle< G4VTouchable > G4TouchableHandle
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4GLOB_DLL std::ostream G4cout
G4VTrajectoryPoint * GetPoint(G4int i) const override
G4int GetPDGEncoding() const override
void AppendStep(const G4Step *aStep) override
G4bool operator==(const G4ClonedRichTrajectory &r) const
void DrawTrajectory() const override
G4int GetTrackID() const override
std::vector< G4AttValue > * CreateAttValues() const override
G4double GetCharge() const override
G4int GetParentID() const override
G4ClonedRichTrajectory()=default
G4double GetInitialKineticEnergy() const
void ShowTrajectory(std::ostream &os=G4cout) const override
G4ParticleDefinition * GetParticleDefinition()
G4ThreeVector GetInitialMomentum() const override
G4String GetParticleName() const override
void MergeTrajectory(G4VTrajectory *secondTrajectory) override
const std::map< G4String, G4AttDef > * GetAttDefs() const override
G4ClonedRichTrajectory & operator=(const G4ClonedRichTrajectory &)=delete
G4int GetPointEntries() const override
G4VTrajectory()=default
#define G4TRACKING_DLL
Definition trkgdefs.hh:45