Geant4-11
G4AdjointTrackingAction.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// G4AdjointTrackingAction class implementation
27//
28// Author: L. Desorgher, SpaceIT GmbH
29// Contract: ESA contract 21435/08/NL/AT
30// Customer: ESA/ESTEC
31// --------------------------------------------------------------------
32
34#include "G4ParticleTable.hh"
35#include "G4Track.hh"
37
38// --------------------------------------------------------------------
39//
42 : theAdjointSteppingAction(anAction)
43{
44}
45
46// --------------------------------------------------------------------
47//
49{
50}
51
52// --------------------------------------------------------------------
53//
55{
56 G4String partType = aTrack->GetParticleDefinition()->GetParticleType();
57 if (G4StrUtil::contains(partType, "adjoint"))
58 {
61 }
62 else
63 {
66 {
68 }
69 }
71}
72
73// --------------------------------------------------------------------
74//
76{
77 // important to have it here !
78 //
81
83 {
84 if (theUserFwdTrackingAction != nullptr)
85 {
87 }
88 }
90 {
97 last_fwd_part_name.erase(0,4);
102 if (aPartDef->GetParticleType() == "adjoint_nucleus")
103 {
104 G4double nb_nuc = G4double(aPartDef->GetBaryonNumber());
105 last_ekin_nuc /= nb_nuc;
106 }
107
109 std::size_t i = 0;
110 while(i< pListOfPrimaryFwdParticles->size() && last_fwd_part_index<0)
111 {
112 if ((*pListOfPrimaryFwdParticles)[i]->GetParticleName()
114 {
116 }
117 ++i;
118 }
119
120 // Fill the vectors
121 //
122 last_pos_vec.push_back(last_pos);
124 last_ekin_vec.push_back(last_ekin);
126 last_cos_th_vec.push_back(last_cos_th);
127 last_weight_vec.push_back(last_weight);
130 }
131 else
132 {
133 }
134}
135
136// --------------------------------------------------------------------
137//
139{
140 last_pos_vec.clear();
141 last_direction_vec.clear();
142 last_ekin_vec.clear();
143 last_ekin_nuc_vec.clear();
144 last_cos_th_vec.clear();
145 last_weight_vec.clear();
148}
double G4double
Definition: G4Types.hh:83
double z() const
double mag() const
void SetAdjointTrackingMode(G4bool aBool)
G4ParticleDefinition * GetLastPartDef()
void SetPrimWeight(G4double weight)
std::vector< G4ThreeVector > last_direction_vec
std::vector< G4double > last_ekin_nuc_vec
virtual void PreUserTrackingAction(const G4Track *)
virtual void PostUserTrackingAction(const G4Track *)
G4UserTrackingAction * theUserFwdTrackingAction
std::vector< G4double > last_cos_th_vec
G4AdjointTrackingAction(G4AdjointSteppingAction *anAction)
G4AdjointSteppingAction * theAdjointSteppingAction
std::vector< G4ThreeVector > last_pos_vec
std::vector< G4ParticleDefinition * > * pListOfPrimaryFwdParticles
std::vector< G4double > last_weight_vec
std::vector< G4int > last_fwd_part_index_vec
std::vector< G4int > last_fwd_part_PDGEncoding_vec
std::vector< G4double > last_ekin_vec
const G4String & GetParticleType() const
const G4String & GetParticleName() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetWeight() const
virtual void PostUserTrackingAction(const G4Track *)
virtual void PreUserTrackingAction(const G4Track *)
G4bool contains(const G4String &str, std::string_view ss)
Check if a string contains a given substring.