Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAEventSet.cc
Go to the documentation of this file.
1// ********************************************************************
2// * License and Disclaimer *
3// * *
4// * The Geant4 software is copyright of the Copyright Holders of *
5// * the Geant4 Collaboration. It is provided under the terms and *
6// * conditions of the Geant4 Software License, included in the file *
7// * LICENSE and available at http://cern.ch/geant4/license . These *
8// * include a list of copyright holders. *
9// * *
10// * Neither the authors of this software system, nor their employing *
11// * institutes,nor the agencies providing financial support for this *
12// * work make any representation or warranty, express or implied, *
13// * regarding this software system or assume any liability for its *
14// * use. Please see the license in the file LICENSE and URL above *
15// * for the full disclaimer and the limitation of liability. *
16// * *
17// * This code implementation is the result of the scientific and *
18// * technical work of the GEANT4 collaboration. *
19// * By using, copying, modifying or distributing the software (or *
20// * any work based on the software) you agree to acknowledge its *
21// * use in resulting scientific publications, and indicate your *
22// * acceptance of all terms of the Geant4 Software license. *
23// ********************************************************************
24//
25#include <memory>
26#include "G4DNAEventSet.hh"
27#include "globals.hh"
28#include "G4UnitsTable.hh"
30
31Event::Event(const G4double& time, const Index& index, ReactionData* pReactionData)
32 : fTimeStep(time)
33 , fIndex(index)
34 , fData(std::pair<std::unique_ptr<JumpingData>, ReactionData*>(nullptr,
35 pReactionData))
36{}
37
38Event::Event(const G4double& time, const Index& index,
39 std::unique_ptr<JumpingData>&& jumping)
40 : fTimeStep(time)
41 , fIndex(index)
42 , fData(std::pair<std::unique_ptr<JumpingData>, ReactionData*>(
43 std::move(jumping), nullptr))
44{}
45
46Event::~Event() = default;
47
49{
50 G4cout << "****PrintEvent::TimeStep : " << G4BestUnit(fTimeStep, "Time")
51 << " index : " << fIndex << " action : ";
52 if(std::get<0>(fData) == nullptr)
53 {
54 G4cout << std::get<1>(fData)->GetReactant1()->GetName() << " + "
55 << std::get<1>(fData)->GetReactant2()->GetName() << " -> "
56 << std::get<1>(fData)->GetProducts()->size() << G4endl;
57 }
58 else
59 {
60 G4cout << std::get<0>(fData)->first->GetName() << " jumping to "
61 << std::get<0>(fData)->second << G4endl;
62 }
63}
64
65G4bool comparatorEventSet::operator()(std::unique_ptr<Event> const& rhs,
66 std::unique_ptr<Event> const& lhs) const
67{
68 return rhs->GetTime() < lhs->GetTime();
69}
70
72 : fEventSet(comparatorEventSet())
73{}
74
75void G4DNAEventSet::CreateEvent(const G4double& time, const Index& index,
76 Event::ReactionData* pReactionData)
77{
78 auto pEvent = std::make_unique<Event>(time, index, pReactionData);
79 AddEvent(std::move(pEvent));
80}
81
82void G4DNAEventSet::CreateEvent(const G4double& time, const Index& index,
83 std::unique_ptr<Event::JumpingData> jum)
84{
85 auto pEvent = std::make_unique<Event>(time, index, std::move(jum));
86 AddEvent(std::move(pEvent));
87}
88
90{
91 auto it = fEventMap.find(key);
92 if(it != fEventMap.end())
93 {
94 fEventSet.erase(it->second);
95 fEventMap.erase(it);
96 }
97}
98
99void G4DNAEventSet::RemoveEvent(EventSet::iterator iter)
100{
101 auto index = (*iter)->GetIndex();
102 RemoveEventOfVoxel(index);
103}
104
105void G4DNAEventSet::AddEvent(std::unique_ptr<Event> pEvent)
106{
107 // idea is no 2 events in one key (or index)
108 auto key = pEvent->GetIndex();
110 auto it = fEventSet.emplace(std::move(pEvent));
111 fEventMap[key] = std::get<0>(it);
112}
113
114//_____________________________________________________________________________________
115
117
118[[maybe_unused]] void G4DNAEventSet::PrintEventSet()
119{
120 G4cout<<G4endl;
121 G4cout << "*****************************************************" << G4endl;
122 G4cout << "G4DNAEventSet::PrintEventSet() of : "<< this->size()<<" events "<< G4endl;
123 for(const auto& it : fEventSet)
124 {
125 (*it).PrintEvent();
126 }
127 G4cout << "End PrintEventSet()" << G4endl;
128 G4cout << "*****************************************************" << G4endl;
129 G4cout << G4endl;
130}
#define G4BestUnit(a, b)
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Event(const G4double &time, const Index &index, ReactionData *)
std::pair< MolType, Index > JumpingData
void PrintEvent() const
virtual ~Event()
virtual ~G4DNAEventSet()
void PrintEventSet()
void RemoveEventOfVoxel(const Index &key)
void AddEvent(std::unique_ptr< Event > pEvent)
void RemoveEventSet()
void CreateEvent(const G4double &time, const Index &index, Event::ReactionData *pReactionData)
void RemoveEvent(EventSet::iterator iter)
G4bool operator()(std::unique_ptr< Event > const &rhs, std::unique_ptr< Event > const &lhs) const