Hall A ROOT/C++ Analyzer (podd)
Loading...
Searching...
No Matches
THaPrimaryKine.cxx
Go to the documentation of this file.
1//*-- Author : Ole Hansen 25-Feb-03
2
4//
5// THaPrimaryKine
6//
7// Calculate kinematics of scattering of the primary (beam) particle.
8// These are usually the electron kinematics.
9// For backwards-compatibility, primary electron kinematics are also
10// available through the derived class THaElectronKine.
11//
13
14#include "THaPrimaryKine.h"
15#include "THaTrackingModule.h"
16#include "THaRunBase.h"
17#include "THaRunParameters.h"
18#include "THaBeam.h"
19#include "VarDef.h"
20#include "TMath.h"
21
22using namespace std;
23using namespace Podd;
24
26
27//_____________________________________________________________________________
28THaPrimaryKine::THaPrimaryKine( const char* name, const char* description,
29 const char* spectro, Double_t particle_mass,
31 THaPhysicsModule(name,description),
32 fQ2(kBig), fOmega(kBig), fW2(kBig), fXbj(kBig), fScatAngle(kBig),
33 fEpsilon(kBig), fQ3mag(kBig), fThetaQ(kBig), fPhiQ(kBig),
34 fM(particle_mass), fMA(target_mass),
35 fSpectroName(spectro), fSpectro(nullptr), fBeam(nullptr)
36{
37 // Standard constructor. Must specify particle mass. Incident particles
38 // are assumed to be along z_lab.
39
40}
41
42//_____________________________________________________________________________
43THaPrimaryKine::THaPrimaryKine( const char* name, const char* description,
44 const char* spectro, const char* beam,
46 : THaPhysicsModule(name,description),
47 fQ2(kBig), fOmega(kBig), fW2(kBig), fXbj(kBig), fScatAngle(kBig),
48 fEpsilon(kBig), fQ3mag(kBig), fThetaQ(kBig), fPhiQ(kBig),
49 fM(-1.0), fMA(target_mass),
50 fSpectroName(spectro), fBeamName(beam), fSpectro(nullptr), fBeam(nullptr)
51{
52 // Constructor with specification of optional beam module.
53 // Particle mass will normally come from the beam module.
54
55}
56
57//_____________________________________________________________________________
59{
60 // Destructor
61
63}
64
65//_____________________________________________________________________________
67{
68 // Clear event-by-event variables.
69
72 = fThetaQ = fPhiQ = kBig;
73 // Clear the 4-vectors
75 fP1 = fA = fA1 = fQ = fP0;
76
77}
78
79//_____________________________________________________________________________
81{
82 // Define/delete global variables.
83
84 RVarDef vars[] = {
85 { "Q2", "4-momentum transfer squared (GeV^2)", "fQ2" },
86 { "omega", "Energy transfer (GeV)", "fOmega" },
87 { "W2", "Invariant mass of recoil system (GeV^2)", "fW2" },
88 { "x_bj", "Bjorken x", "fXbj" },
89 { "angle", "Scattering angle (rad)", "fScatAngle" },
90 { "epsilon", "Virtual photon polarization factor", "fEpsilon" },
91 { "q3m", "Magnitude of 3-momentum transfer", "fQ3mag" },
92 { "th_q", "Theta of 3-momentum vector (rad)", "fThetaQ" },
93 { "ph_q", "Phi of 3-momentum vector (rad)", "fPhiQ" },
94 { "nu", "Energy transfer (GeV)", "fOmega" },
95 { "q_x", "x-cmp of Photon vector in the lab", "fQ.X()" },
96 { "q_y", "y-cmp of Photon vector in the lab", "fQ.Y()" },
97 { "q_z", "z-cmp of Photon vector in the lab", "fQ.Z()" },
98 { nullptr }
99 };
100 return DefineVarsFromList( vars, mode );
101}
102
103//_____________________________________________________________________________
105{
106 // Initialize the module.
107 // Locate the spectrometer apparatus named in fSpectroName and save
108 // pointer to it.
109
110 fSpectro = dynamic_cast<THaTrackingModule*>
111 ( FindModule( fSpectroName.Data(), "THaTrackingModule"));
112 if( !fSpectro )
113 return fStatus;
114
115 // Optional beam apparatus
116 if( fBeamName.Length() > 0 ) {
117 fBeam = dynamic_cast<THaBeamModule*>
118 ( FindModule( fBeamName.Data(), "THaBeamModule") );
119 if( !fBeam )
120 return fStatus;
121 if( fM <= 0.0 )
122 fM = fBeam->GetBeamInfo()->GetM();
123 }
124
125 // Standard initialization. Calls this object's DefineVariables().
126 if( THaPhysicsModule::Init( run_time ) != kOK )
127 return fStatus;
128
129 if( fM <= 0.0 ) {
130 Error( Here("Init"), "Particle mass not defined. Module "
131 "initialization failed" );
133 }
134 return fStatus;
135}
136
137//_____________________________________________________________________________
139{
140 // Calculate electron kinematics for the Golden Track of the spectrometer
141
142 if( !IsOK() || !gHaRun ) return -1;
143
145 if( !trkifo || !trkifo->IsOK() ) return 1;
146
147 // Determine 4-momentum of incident particle.
148 // If a beam module given, use it to get the beam momentum. This
149 // module may apply corrections for beam energy loss, variations, etc.
150 if( fBeam ) {
152 } else {
153 // If no beam given, assume beam along z_lab
155 fP0.SetXYZM( 0.0, 0.0, p_in, fM );
156 }
157
158 fP1.SetVectM( trkifo->GetPvect(), fM );
159 fA.SetXYZM( 0.0, 0.0, 0.0, fMA ); // Assume target at rest
160
161 // proton mass (for x_bj)
162 const Double_t Mp = 0.938;
163
164 // Standard electron kinematics
165 fQ = fP0 - fP1;
166 fQ2 = -fQ.M2();
167 fQ3mag = fQ.P();
168 fOmega = fQ.E();
169 fA1 = fA + fQ;
170 fW2 = fA1.M2();
171 fScatAngle = fP0.Angle( fP1.Vect() );
172 fEpsilon = 1.0 / ( 1.0 + 2.0*fQ3mag*fQ3mag/fQ2*
173 TMath::Power( TMath::Tan(fScatAngle/2.0), 2.0 ));
174 fThetaQ = fQ.Theta();
175 fPhiQ = fQ.Phi();
176 fXbj = fQ2/(2.0*Mp*fOmega);
177
178 fDataValid = true;
179 return 0;
180}
181
182//_____________________________________________________________________________
184{
185 // Query the run database for particle and target masses.
186 //
187 // If the particle or target masses weren't given to the constructor,
188 // then search for them in the run database.
189 //
190 // For the particle mass, first search for "<prefix>.M", then,
191 // if not found, for "M". If still not found, use electron mass mass.
192 //
193 // For the target mass, first search for "<prefix>.MA", then, if not found,
194 // for "MA". If still not found, use proton mass.
195
197 if( err ) return err;
198
199 if ( fM > 0.0 && fMA > 0.0 )
200 return 0;
201
202 FILE* f = OpenRunDBFile( date );
203 if( !f ) return kFileError;
204
205 if( fM <= 0.0 ) {
206 TString name(fPrefix), tag("M"); name += tag;
207 Int_t st = LoadDBvalue( f, date, name.Data(), fM );
208 if( st )
209 LoadDBvalue( f, date, tag.Data(), fM );
210 if( fM <= 0.0 ) fM = 0.511e-3;
211 }
212 if( fMA <= 0.0 ) {
213 TString name(fPrefix), tag("MA"); name += tag;
214 Int_t st = LoadDBvalue( f, date, name.Data(), fMA );
215 if( st )
216 LoadDBvalue( f, date, tag.Data(), fMA );
217 if( fMA <= 0.0 ) fMA = 0.938;
218 }
219
220 fclose(f);
221 return 0;
222}
223
224//_____________________________________________________________________________
226{
227 if( !IsInit())
228 fM = m;
229 else
230 PrintInitError("SetMass()");
231}
232
233//_____________________________________________________________________________
235{
236 if( !IsInit())
237 fMA = m;
238 else
239 PrintInitError("SetTargetMass()");
240}
241
242//_____________________________________________________________________________
243void THaPrimaryKine::SetSpectrometer( const char* name )
244{
245 if( !IsInit())
247 else
248 PrintInitError("SetSpectrometer()");
249}
250
251//_____________________________________________________________________________
252void THaPrimaryKine::SetBeam( const char* name )
253{
254 if( !IsInit())
255 fBeamName = name;
256 else
257 PrintInitError("SetBeam()");
258}
int Int_t
const Data_t kBig
Definition DataType.h:15
#define f(i)
double Double_t
const char Option_t
Option_t Option_t TPoint TPoint const char mode
char name[80]
R__EXTERN class THaRunBase * gHaRun
Definition THaGlobals.h:16
static const Double_t target_mass
static const Double_t Mp
virtual Int_t ReadRunDatabase(const TDatime &date)
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)
virtual FILE * OpenRunDBFile(const TDatime &date)
Bool_t IsOK() const
Bool_t IsInit() const
Double_t GetM() const
const TVector3 & GetPvect() const
Definition THaBeamInfo.h:56
THaBeamInfo * GetBeamInfo()
void PrintInitError(const char *here)
virtual void Clear(Option_t *opt="")
void SetSpectrometer(const char *name)
Double_t fScatAngle
void SetTargetMass(Double_t m)
virtual Int_t DefineVariables(EMode mode=kDefine)
virtual Int_t ReadRunDatabase(const TDatime &date)
THaTrackingModule * fSpectro
virtual void Clear(Option_t *opt="")
virtual ~THaPrimaryKine()
void SetMass(Double_t m)
THaBeamModule * fBeam
THaPrimaryKine(const char *name, const char *description, const char *spectro="", Double_t particle_mass=0.0, Double_t target_mass=0.0)
void SetBeam(const char *name)
virtual Int_t Process(const THaEvData &)
THaRunParameters * GetParameters() const
Definition THaRunBase.h:61
Double_t GetBeamP() const
const TVector3 & GetPvect() const
Bool_t IsOK() const
THaTrackInfo * GetTrackInfo()
Double_t Angle(const TVector3 &v) const
void SetXYZT(Double_t x, Double_t y, Double_t z, Double_t t)
TVector3 Vect() const
Double_t M2() const
Double_t P() const
Double_t Theta() const
Double_t Phi() const
Double_t E() const
void SetVectM(const TVector3 &spatial, Double_t mass)
void SetXYZM(Double_t x, Double_t y, Double_t z, Double_t m)
virtual void Error(const char *method, const char *msgfmt,...) const
Ssiz_t Length() const
const char * Data() const
Double_t Power(Double_t x, Double_t y)
Double_t Tan(Double_t)
STL namespace.
TMarker m
ClassImp(TPyArg)