Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ClonedSmoothTrajectory.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// G4ClonedSmoothTrajectory
27//
28// Class description:
29//
30// This class is identical to G4SmoothTrajectory 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 G4SmoothTrajectory object.
35//
36// Makoto Asai (JLab) - Oct.2024
37// --------------------------------------------------------------------
38#ifndef G4ClonedSmoothTrajectory_hh
39#define G4ClonedSmoothTrajectory_hh 1
40
41#include "G4Allocator.hh"
42#include "G4ParticleDefinition.hh" // Include from 'particle+matter'
43#include "G4ClonedSmoothTrajectoryPoint.hh" // Include from 'tracking'
44#include "G4Step.hh"
45#include "G4Track.hh"
46#include "G4VTrajectory.hh"
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 //
68
69 // Operators
70 //
71 inline G4bool operator==(const G4ClonedSmoothTrajectory& r) const;
72 inline void* operator new(size_t);
73 inline void operator delete(void*);
74
75 // Get/Set functions
76 //
77 inline G4int GetTrackID() const override { return fTrackID; }
78 inline G4int GetParentID() const override { return fParentID; }
79 inline G4String GetParticleName() const override { return ParticleName; }
80 inline G4double GetCharge() const override { return PDGCharge; }
81 inline G4int GetPDGEncoding() const override { return PDGEncoding; }
82 inline G4double GetInitialKineticEnergy() const { return initialKineticEnergy; }
83 inline G4ThreeVector GetInitialMomentum() const override { return initialMomentum; }
84
85 // Other member functions
86 //
87 void ShowTrajectory(std::ostream& os = G4cout) const override;
88 void DrawTrajectory() const override;
89 void AppendStep(const G4Step* aStep) override;
90 G4int GetPointEntries() const override { return (G4int)positionRecord->size(); }
91 G4VTrajectoryPoint* GetPoint(G4int i) const override { return (*positionRecord)[i]; }
92 void MergeTrajectory(G4VTrajectory* secondTrajectory) override;
93
95
96 // Get method for HEPRep style attributes
97 //
98 const std::map<G4String, G4AttDef>* GetAttDefs() const override;
99 std::vector<G4AttValue>* CreateAttValues() const override;
100
101 private:
102 G4TrajectoryPointContainer* positionRecord = nullptr;
103 G4int fTrackID = 0;
104 G4int fParentID = 0;
105 G4int PDGEncoding = 0;
106 G4double PDGCharge = 0.0;
107 G4String ParticleName = "";
108 G4double initialKineticEnergy = 0.0;
109 G4ThreeVector initialMomentum;
110};
111
113
114inline void* G4ClonedSmoothTrajectory::operator new(size_t)
115{
116 if (aClonedSmoothTrajectoryAllocator() == nullptr) {
118 }
119 return (void*)aClonedSmoothTrajectoryAllocator()->MallocSingle();
120}
121
122inline void G4ClonedSmoothTrajectory::operator delete(void* aTrajectory)
123{
125}
126
128{
129 return (this == &r);
130}
131
132#endif
G4TRACKING_DLL G4Allocator< G4ClonedSmoothTrajectory > *& aClonedSmoothTrajectoryAllocator()
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4GLOB_DLL std::ostream G4cout
void DrawTrajectory() const override
G4bool operator==(const G4ClonedSmoothTrajectory &r) const
G4int GetPDGEncoding() const override
G4VTrajectoryPoint * GetPoint(G4int i) const override
G4int GetPointEntries() const override
G4int GetParentID() const override
const std::map< G4String, G4AttDef > * GetAttDefs() const override
void ShowTrajectory(std::ostream &os=G4cout) const override
G4ThreeVector GetInitialMomentum() const override
std::vector< G4AttValue > * CreateAttValues() const override
G4int GetTrackID() const override
G4ParticleDefinition * GetParticleDefinition()
G4String GetParticleName() const override
void MergeTrajectory(G4VTrajectory *secondTrajectory) override
void AppendStep(const G4Step *aStep) override
G4double GetCharge() const override
G4ClonedSmoothTrajectory()=default
G4VTrajectory()=default
#define G4TRACKING_DLL
Definition trkgdefs.hh:45