Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4AdjointTrackingAction.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// G4AdjointTrackingAction
27//
28// Class description:
29//
30// This class represents actions taken place at the start/end point
31// of processing one track during an adjoint simulation
32
33// Author: L. Desorgher, SpaceIT GmbH
34// Contract: ESA contract 21435/08/NL/AT
35// Customer: ESA/ESTEC
36// --------------------------------------------------------------------
37#ifndef G4AdjointTrackingAction_hh
38#define G4AdjointTrackingAction_hh 1
39
40#include "G4ThreeVector.hh"
42#include "globals.hh"
43
44#include <vector>
45
47class G4Track;
49
51{
52 public:
54 ~G4AdjointTrackingAction() override = default;
55
56 void PreUserTrackingAction(const G4Track*) override;
57 void PostUserTrackingAction(const G4Track*) override;
60
62 {
63 theUserFwdTrackingAction = anAction;
64 }
66 {
67 return last_pos_vec[i];
68 }
70 {
71 return last_direction_vec[i];
72 }
73 inline G4double GetEkinAtEndOfLastAdjointTrack(std::size_t i = 0) { return last_ekin_vec[i]; }
75 {
76 return last_ekin_nuc_vec[i];
77 }
78 inline G4double GetWeightAtEndOfLastAdjointTrack(std::size_t i = 0) { return last_weight_vec[i]; }
79 inline G4double GetCosthAtEndOfLastAdjointTrack(std::size_t i = 0) { return last_cos_th_vec[i]; }
80 inline const G4String& GetFwdParticleNameAtEndOfLastAdjointTrack() { return last_fwd_part_name; }
82 {
83 return last_fwd_part_PDGEncoding_vec[i];
84 }
85 inline G4bool GetIsAdjointTrackingMode() { return is_adjoint_tracking_mode; }
86 inline G4int GetLastFwdParticleIndex(std::size_t i = 0) { return last_fwd_part_index_vec[i]; }
87 inline std::size_t GetNbOfAdointTracksReachingTheExternalSurface() { return last_pos_vec.size(); }
88 inline void SetListOfPrimaryFwdParticles(std::vector<G4ParticleDefinition*>* aListOfParticles)
89 {
90 pListOfPrimaryFwdParticles = aListOfParticles;
91 }
92
93 private:
94 G4AdjointSteppingAction* theAdjointSteppingAction = nullptr;
95 G4UserTrackingAction* theUserFwdTrackingAction = nullptr;
96 G4bool is_adjoint_tracking_mode = false;
97
98 // Adjoint particle information on the external surface
99 // ----------------------------------------------------
100 G4ThreeVector last_pos;
101 G4ThreeVector last_direction;
102 G4double last_ekin = 0.0, last_ekin_nuc = 0.0;
103 // last_ekin_nuc=last_ekin/nuc, nuc is 1 if not a nucleus
104 G4double last_cos_th = 0.0;
105 G4String last_fwd_part_name;
106 G4int last_fwd_part_PDGEncoding = 0;
107 G4double last_weight = 0.0;
108 G4int last_fwd_part_index = 0;
109 std::vector<G4ParticleDefinition*>* pListOfPrimaryFwdParticles = nullptr;
110
111 std::vector<G4ThreeVector> last_pos_vec;
112 std::vector<G4ThreeVector> last_direction_vec;
113 std::vector<G4double> last_ekin_vec;
114 std::vector<G4double> last_ekin_nuc_vec;
115 std::vector<G4double> last_cos_th_vec;
116 std::vector<G4double> last_weight_vec;
117 std::vector<G4int> last_fwd_part_PDGEncoding_vec;
118 std::vector<G4int> last_fwd_part_index_vec;
119};
120
121#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4double GetEkinAtEndOfLastAdjointTrack(std::size_t i=0)
void SetListOfPrimaryFwdParticles(std::vector< G4ParticleDefinition * > *aListOfParticles)
G4double GetEkinNucAtEndOfLastAdjointTrack(std::size_t i=0)
G4ThreeVector GetPositionAtEndOfLastAdjointTrack(std::size_t i=0)
const G4String & GetFwdParticleNameAtEndOfLastAdjointTrack()
G4int GetFwdParticlePDGEncodingAtEndOfLastAdjointTrack(std::size_t i=0)
G4double GetWeightAtEndOfLastAdjointTrack(std::size_t i=0)
G4AdjointTrackingAction(G4AdjointSteppingAction *anAction)
std::size_t GetNbOfAdointTracksReachingTheExternalSurface()
G4ThreeVector GetDirectionAtEndOfLastAdjointTrack(std::size_t i=0)
void PostUserTrackingAction(const G4Track *) override
G4int GetLastFwdParticleIndex(std::size_t i=0)
~G4AdjointTrackingAction() override=default
G4double GetCosthAtEndOfLastAdjointTrack(std::size_t i=0)
void SetUserForwardTrackingAction(G4UserTrackingAction *anAction)
void PreUserTrackingAction(const G4Track *) override