Hall A ROOT/C++ Analyzer (podd)
Loading...
Searching...
No Matches
THaTrackEloss.cxx
Go to the documentation of this file.
1//*-- Author : Ole Hansen 04-May-04
2
4//
5// THaTrackEloss
6//
7// Correct energy loss of a track by a constant. This is basically trivial,
8// but useful to have to provide input for kinematics modules.
9// More advanced corrections can be implemented by inheritance from
10// this module.
11//
13
14#include "THaTrackEloss.h"
15#include "THaSpectrometer.h"
16#include "THaTrack.h"
17#include "THaTrackInfo.h"
18#include "TMath.h"
19
20using namespace std;
21
22//_____________________________________________________________________________
24 const char* description,
25 const char* input_tracks,
26 Double_t particle_mass,
27 Int_t hadron_charge ) :
28 THaElossCorrection(name,description,input_tracks,particle_mass,
29 hadron_charge), fTrackModule(nullptr)
30{
31 // Normal constructor.
32
33 Clear();
34}
35
36//_____________________________________________________________________________
38{
39 // Destructor
40
42}
43
44//_____________________________________________________________________________
46{
47 // Compute the expected energy loss for the track given in trkifo.
48 //
49 // Currently, we use a very simple algorithm that computes
50 // the energy loss based on a fixed material thickness.
51 // Only the beta dependence of the energy loss is used,
52 // but there are no corrections for angle etc.
53 //
54 // May be overridden by derived classes as necessary.
55
56 Double_t p0 = trkifo->GetP();
57 Double_t beta = p0 / TMath::Sqrt(p0*p0 + fM*fM);
58 if( fElectronMode ) {
61 } else {
64 }
65}
66
67//_____________________________________________________________________________
69{
70 // Clear all event-by-event variables.
71
74}
75
76//_____________________________________________________________________________
78{
79 // Initialize the module.
80 // Locate the input tracking module named in fInputName and save
81 // pointer to it.
82
83 static const char* const here = "Init()";
84
85 // Find the input tracking module
86 fTrackModule = dynamic_cast<THaTrackingModule*>
87 ( FindModule( fInputName.Data(), "THaTrackingModule"));
88 if( !fTrackModule )
89 return fStatus;
90
91 // Get the parent spectrometer apparatus from the input module
93 if( !spectro ) {
94 Error( Here(here), "Oops. Input tracking module has no pointer "
95 "to a spectrometer?!?" );
96 return fStatus = kInitError;
97 }
98 // Needed for initialization of dependent modules in a chain
99 fTrkIfo.SetSpectrometer( spectro );
100
101 // Standard initialization. Calls this object's DefineVariables() and
102 // reads material properties from the run database.
103 THaElossCorrection::Init( run_time );
104
105 return fStatus;
106}
107
108//_____________________________________________________________________________
110{
111 // Define/delete global variables.
112
114 if( ret )
115 return ret;
116
118}
119
120//_____________________________________________________________________________
122{
123 // Calculate corrections and adjust the track parameters.
124
125 if( !IsOK() ) return -1;
126
128 if( !trkifo->IsOK() ) return 2;
129
130 // Copy the input track info
131 fTrkIfo = *trkifo;
132
133 // Compute the correction
134 Double_t p_out = fTrkIfo.GetP();
135 if( p_out <= 0.0 ) return 4; //oops
136 Double_t E_out = TMath::Sqrt(p_out*p_out + fM*fM);
137 if( !fTestMode ) {
138 // calculate pathlength for this event if we have a vertex module
139 if( fExtPathMode ) {
140 if( !fVertexModule->HasVertex() )
141 return 1;
142 fPathlength =
144 }
145 //update fEloss
146 CalcEloss(trkifo);
147 }
148 Double_t p_in = TMath::Sqrt(p_out*p_out + fEloss*fEloss + 2.0*E_out*fEloss);
149
150 // Apply the correction
151 fTrkIfo.SetP(p_in);
152
153 fDataValid = true;
154 return 0;
155}
156
157//_____________________________________________________________________________
159
int Int_t
double Double_t
const char Option_t
Option_t Option_t TPoint TPoint const char mode
char name[80]
static const char *const here
Definition THaVar.cxx:64
static Int_t DefineVarsFromList(const void *list, EType type, EMode mode, const char *def_prefix, const TObject *obj, const char *prefix, const char *here, const char *comment_subst="")
virtual const char * Here(const char *) const
THaAnalysisObject * FindModule(const char *name, const char *classname, bool do_error=true)
Bool_t IsOK() const
virtual Int_t DefineVariables(EMode mode=kDefine)
virtual void Clear(Option_t *opt="")
static Double_t ElossHadron(Int_t Z_hadron, Double_t beta, Double_t z_med, Double_t a_med, Double_t d_med, Double_t pathlength)
static Double_t ElossElectron(Double_t beta, Double_t z_med, Double_t a_med, Double_t d_med, Double_t pathlength)
THaVertexModule * fVertexModule
virtual void CalcEloss(THaTrackInfo *trkifo)
THaTrackingModule * fTrackModule
virtual Int_t Process(const THaEvData &)
virtual void Clear(Option_t *opt="")
THaTrackEloss(const char *name, const char *description, const char *input_tracks="", Double_t particle_mass=0.511e-3, Int_t hadron_charge=1)
virtual Int_t DefineVariables(EMode mode=kDefine)
virtual ~THaTrackEloss()
THaSpectrometer * GetSpectrometer() const
void SetSpectrometer(THaSpectrometer *s)
void SetP(Double_t p)
Double_t GetP() const
Bool_t IsOK() const
static const RVarDef * GetRVarDef()
THaTrackInfo * GetTrackInfo()
virtual Bool_t HasVertex() const
virtual const TVector3 & GetVertex() const
virtual void Error(const char *method, const char *msgfmt,...) const
const char * Data() const
Double_t Z() const
double beta(double x, double y)
Double_t Sqrt(Double_t x)
Double_t Abs(Double_t d)
STL namespace.
ClassImp(TPyArg)