Geant4 Cross Reference |
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 // G4ClonedTrajectory 27 // 28 // Class description: 29 // 30 // This class is identical to G4Trajectory 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 G4Trajectory object. 35 // 36 // Makoto Asai (JLab) - Oct.2024 37 // -------------------------------------------------------------------- 38 // 39 #ifndef G4ClonedTrajectory_hh 40 #define G4ClonedTrajectory_hh 1 41 42 #include "G4Allocator.hh" 43 #include "G4ParticleDefinition.hh" // Include from 'particle+matter' 44 #include "G4Step.hh" 45 #include "G4Track.hh" 46 #include "G4ClonedTrajectoryPoint.hh" // Include from 'tracking' 47 #include "G4VTrajectory.hh" 48 #include "G4ios.hh" // Include from 'system' 49 #include "globals.hh" // Include from 'global' 50 #include "G4Threading.hh" 51 52 #include "trkgdefs.hh" 53 #include <stdlib.h> // Include from 'system' 54 55 #include <vector> 56 57 class G4Polyline; 58 class G4Trajectory; 59 60 class G4ClonedTrajectory : public G4VTrajectory 61 { 62 using G4ClonedTrajectoryPointContainer = std::vector<G4VTrajectoryPoint*>; 63 64 public: 65 // Constructors/Destructor 66 67 G4ClonedTrajectory() = default; 68 G4ClonedTrajectory(const G4Trajectory&); 69 ~G4ClonedTrajectory() override; 70 71 // Operators 72 73 inline void* operator new(size_t); 74 inline void operator delete(void*); 75 inline G4bool operator==(const G4ClonedTrajectory& r) const; 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 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 G4int GetPointEntries() const override { return G4int(positionRecord->size()); } 93 G4VTrajectoryPoint* GetPoint(G4int i) const override { return (*positionRecord)[i]; } 94 void MergeTrajectory(G4VTrajectory* secondTrajectory) override; 95 96 G4ParticleDefinition* GetParticleDefinition(); 97 98 const std::map<G4String, G4AttDef>* GetAttDefs() const override; 99 std::vector<G4AttValue>* CreateAttValues() const override; 100 101 private: 102 G4ClonedTrajectoryPointContainer* positionRecord = nullptr; 103 G4int fTrackID = 0; 104 G4int fParentID = 0; 105 G4int PDGEncoding = 0; 106 G4double PDGCharge = 0.0; 107 G4String ParticleName = "dummy"; 108 G4double initialKineticEnergy = 0.0; 109 G4ThreeVector initialMomentum; 110 }; 111 112 extern G4TRACKING_DLL G4Allocator<G4ClonedTrajectory>*& aClonedTrajectoryAllocator(); 113 114 inline void* G4ClonedTrajectory::operator new(size_t) 115 { 116 if (aClonedTrajectoryAllocator() == nullptr) { 117 aClonedTrajectoryAllocator() = new G4Allocator<G4ClonedTrajectory>; 118 } 119 return (void*)aClonedTrajectoryAllocator()->MallocSingle(); 120 } 121 122 inline void G4ClonedTrajectory::operator delete(void* aTrajectory) 123 { 124 aClonedTrajectoryAllocator()->FreeSingle((G4ClonedTrajectory*)aTrajectory); 125 } 126 127 inline G4bool G4ClonedTrajectory::operator==(const G4ClonedTrajectory& r) const { return (this == &r); } 128 129 #endif 130