Hall A ROOT/C++ Analyzer (podd)
Loading...
Searching...
No Matches
THaSecondaryKine.cxx
Go to the documentation of this file.
1//*-- Author : Ole Hansen 17-Mar-03
2
4//
5// THaSecondaryKine
6//
7// Secondary particle kinematics module.
8// This module calculates kinematic quantities for the hadron side of
9// two-arm coincidence reactions. The reaction is
10//
11// e + A -> e' + X + B
12//
13// where e' and X are measured, and B is undetected.
14// As input, the module requires information about the electron kinematics
15// (from a THaPrimaryKine module, e.g. THaElectronKine), the momentum of
16// the secondary particle X (from a THaTrackingModule, e.g. a THaSpectrometer),
17// and the mass of X.
18//
20
21#include "THaSecondaryKine.h"
22#include "THaPrimaryKine.h"
23#include "THaTrackingModule.h"
24#include "VarDef.h"
25#include "TLorentzVector.h"
26#include "TVector3.h"
27#include "TMath.h"
28#include "TRotation.h"
29
30using namespace std;
31using namespace Podd;
32
34
35//_____________________________________________________________________________
36THaSecondaryKine::THaSecondaryKine( const char* name, const char* description,
37 const char* secondary_spectro,
38 const char* primary_kine,
39 Double_t secondary_mass ) :
40 THaPhysicsModule(name,description),
41 fTheta_xq(kBig), fPhi_xq(kBig), fTheta_bq(kBig), fPhi_bq(kBig),
42 fXangle(kBig), fPmiss(kBig), fPmiss_x(kBig), fPmiss_y(kBig), fPmiss_z(kBig),
43 fEmiss(kBig), fMrecoil(kBig), fErecoil(kBig), fTX(kBig), fTB(kBig),
44 fPX_cm(kBig), fTheta_x_cm(kBig), fPhi_x_cm(kBig), fTheta_b_cm(kBig),
45 fPhi_b_cm(kBig), fTX_cm(kBig), fTB_cm(kBig), fTtot_cm(kBig), fMandelS(kBig),
46 fMandelT(kBig), fMandelU(kBig),
47 fMX(secondary_mass), fSpectroName(secondary_spectro), fSpectro(nullptr),
48 fPrimaryName(primary_kine), fPrimary(nullptr)
49{
50 // Constructor
51
52}
53
54//_____________________________________________________________________________
56{
57 // Destructor
58
60}
61
62//_____________________________________________________________________________
76
77//_____________________________________________________________________________
79{
80 // Define/delete global variables.
81
82 RVarDef vars[] = {
83 { "th_xq", "Polar angle of detected particle with q (rad)",
84 "fTheta_xq" },
85 { "ph_xq", "Azimuth of detected particle with scattering plane (rad)",
86 "fPhi_xq" },
87 { "th_bq", "Polar angle of recoil system with q (rad)", "fTheta_bq" },
88 { "ph_bq", "Azimuth of recoil system with scattering plane (rad)",
89 "fPhi_bq" },
90 { "xangle", "Angle of detected particle with scattered electron (rad)",
91 "fXangle" },
92 { "pmiss", "Missing momentum magnitude (GeV), nuclear physics "
93 "definition (-pB)", "fPmiss" },
94 { "pmiss_x", "x-component of p_miss wrt q (GeV)", "fPmiss_x" },
95 { "pmiss_y", "y-component of p_miss wrt q (GeV)", "fPmiss_y" },
96 { "pmiss_z", "z-component of p_miss, along q (GeV)", "fPmiss_z" },
97 { "emiss", "Missing energy (GeV), nuclear physics definition "
98 "omega-Tx-Tb", "fEmiss" },
99 { "Mrecoil", "Invariant mass of recoil system (GeV)", "fMrecoil" },
100 { "Erecoil", "Total energy of recoil system (GeV)", "fErecoil" },
101 { "Prec_x", "x-component of recoil system in lab (GeV/c)", "fB.X()" },
102 { "Prec_y", "y-component of recoil system in lab (GeV/c)", "fB.Y()" },
103 { "Prec_z", "z-component of recoil system in lab (GeV/c)", "fB.Z()" },
104 { "tx", "Kinetic energy of detected particle (GeV)", "fTX" },
105 { "tb", "Kinetic energy of recoil system (GeV)", "fTB" },
106 { "px_cm", "Magnitude of X momentum in CM system (GeV)", "fPX_cm" },
107 { "thx_cm", "Polar angle of X in CM system wrt q (rad)", "fTheta_x_cm" },
108 { "phx_cm", "Azimuth of X in CM system wrt q (rad)", "fPhi_x_cm" },
109 { "thb_cm", "Polar angle of recoil systm in CM wrt q (rad)",
110 "fTheta_b_cm" },
111 { "phb_cm", "Azimuth of recoil system in CM wrt q (rad)", "fPhi_b_cm" },
112 { "tx_cm", "Kinetic energy of X in CM (GeV)", "fTX_cm" },
113 { "tb_cm", "Kinetic energy of B in CM (GeV)", "fTB_cm" },
114 { "t_tot_cm", "Total CM kinetic energy", "fTtot_cm" },
115 { "MandelS", "Mandelstam s for secondary vertex (GeV^2)", "fMandelS" },
116 { "MandelT", "Mandelstam t for secondary vertex (GeV^2)", "fMandelT" },
117 { "MandelU", "Mandelstam u for secondary vertex (GeV^2)", "fMandelU" },
118 { nullptr }
119 };
120 return DefineVarsFromList( vars, mode );
121}
122
123//_____________________________________________________________________________
125{
126 // Initialize the module.
127 // Locate the kinematics module named in fPrimaryName, which describes
128 // the primary interaction kinematics, and save pointer to it.
129
130 // Standard initialization. Calls this object's DefineVariables().
131 if( THaPhysicsModule::Init( run_time ) != kOK )
132 return fStatus;
133
134 fSpectro = dynamic_cast<THaTrackingModule*>
135 ( FindModule( fSpectroName.Data(), "THaTrackingModule"));
136 if( fStatus ) return fStatus;
137
138 fPrimary = dynamic_cast<THaPrimaryKine*>
139 ( FindModule( fPrimaryName.Data(), "THaPrimaryKine"));
140 return fStatus;
141}
142
143//_____________________________________________________________________________
145{
146 // Calculate the kinematics.
147
148 if( !IsOK() ) return -1;
149
150 // Tracking information from the secondary spectrometer
152 if( !trkifo || !trkifo->IsOK() ) return 1;
153
154 // Require valid input data
155 if( !fPrimary->DataValid() ) return 2;
156
157 // Measured momentum of secondary particle X in lab
158 const TVector3& pX3 = trkifo->GetPvect();
159 // 4-momentum of X
160 fX.SetVectM( pX3, fMX );
161
162 // 4-momenta of the the primary interaction
163 const TLorentzVector* pA = fPrimary->GetA(); // Initial target
164 const TLorentzVector* pA1 = fPrimary->GetA1(); // Final target
165 const TLorentzVector* pQ = fPrimary->GetQ(); // Momentum xfer
166 const TLorentzVector* pP1 = fPrimary->GetP1(); // Final electron
167 Double_t omega = fPrimary->GetOmega(); // Energy xfer
168
169 // 4-momentum of undetected recoil system B
170 fB = *pA1 - fX;
171
172 // Angle of X with scattered primary particle
173 fXangle = fX.Angle( pP1->Vect());
174
175 // Angles of X and B wrt q-vector
176 // xq and bq are the 3-momentum vectors of X and B expressed in
177 // the coordinate system where q is the z-axis and the x-axis
178 // lies in the scattering plane (defined by q and e') and points
179 // in the direction of e', so the out-of-plane angle lies within
180 // -90<phi_xq<90deg if X is detected on the downstream/forward side of q.
181 TRotation rot_to_q;
182 rot_to_q.SetZAxis( pQ->Vect(), pP1->Vect()).Invert();
183 TVector3 xq = fX.Vect();
184 TVector3 bq = fB.Vect();
185 xq *= rot_to_q;
186 bq *= rot_to_q;
187 fTheta_xq = xq.Theta(); //"theta_xq"
188 fPhi_xq = xq.Phi(); //"out-of-plane angle", "phi"
189 fTheta_bq = bq.Theta();
190 fPhi_bq = bq.Phi();
191
192 // Missing momentum and components wrt q-vector
193 // The definition of p_miss as the negative of the undetected recoil
194 // momentum is the standard nuclear physics convention.
195 TVector3 p_miss = -bq;
196 fPmiss = p_miss.Mag(); //=fB.P()
197 //The missing momentum components in the q coordinate system.
198 //FIXME: just store p_miss in member variable
199 fPmiss_x = p_miss.X();
200 fPmiss_y = p_miss.Y();
201 fPmiss_z = p_miss.Z();
202
203 // Invariant mass of the recoil system, a.k.a. "missing mass".
204 // This invariant mass equals MB(ground state) plus any excitation energy.
205 fMrecoil = fB.M();
206
207 // Kinetic energies of X and B
208 fTX = fX.E() - fMX;
209 fTB = fB.E() - fMrecoil;
210
211 // Standard nuclear physics definition of "missing energy":
212 // binding energy of X in the target (= removal energy of X).
213 // NB: If X is knocked out of a lower shell, the recoil system carries
214 // a significant excitation energy. This excitation is included in Emiss
215 // here, as it should, since it results from the binding of X.
216 fEmiss = omega - fTX - fTB;
217
218 // In production reactions, the "missing energy" is defined
219 // as the total energy of the undetected recoil system.
220 // This is the "missing mass", Mrecoil, plus any kinetic energy.
221 fErecoil = fB.E();
222
223 // Calculate some interesting quantities in the CM system of A'.
224 // NB: If the target is initially at rest, the A'-vector (spatial part)
225 // is the same as the q-vector, so we could just reuse the q rotation
226 // matrix. The following is completely general, i.e. allows for a moving
227 // target.
228
229 // Boost of the A' system.
230 // NB: ROOT's Boost() boosts particle to lab if the boost vector points
231 // along the particle's lab velocity, so we need the negative here to
232 // boost from the lab to the particle frame.
233 Double_t beta = pA1->P()/(pA1->E());
234 TVector3 boost(0.,0.,-beta);
235
236 // CM vectors of X and B.
237 // Express X and B in the frame where q is along the z-axis
238 // - the typical head-on collision picture.
239 TRotation rot_to_A1;
240 rot_to_A1.SetZAxis( pA1->Vect(), pP1->Vect()).Invert();
241 TVector3 x_cm_vect = fX.Vect();
242 x_cm_vect *= rot_to_A1;
243 TLorentzVector x_cm(x_cm_vect,fX.E());
244 x_cm.Boost(boost);
245 TLorentzVector b_cm(-x_cm.Vect(), pA1->E()-x_cm.E());
246 fPX_cm = x_cm.P();
247 // pB_cm, by construction, is the same as pX_cm.
248
249 // CM angles of X and B in the A' frame
250 fTheta_x_cm = x_cm.Theta();
251 fPhi_x_cm = x_cm.Phi();
252 fTheta_b_cm = b_cm.Theta();
253 fPhi_b_cm = b_cm.Phi();
254
255 // CM kinetic energies of X and B and total kinetic energy
256 fTX_cm = x_cm.E() - fMX;
257 fTB_cm = b_cm.E() - fMrecoil;
259
260 // Mandelstam variables for the secondary vertex.
261 // These variables are defined for the reaction gA->XB,
262 // where g is the virtual photon (with momentum q), and A, X, and B
263 // are as before.
264
265 fMandelS = (*pQ+*pA).M2();
266 fMandelT = (*pQ-fX).M2();
267 fMandelU = (*pQ-fB).M2();
268
269 fDataValid = true;
270 return 0;
271}
272
273//_____________________________________________________________________________
275{
276 // Query the run database for any parameters of the module that were not
277 // set by the constructor. This module has only one parameter,
278 // the mass of the detected secondary particle X.
279 //
280 // First searches for "<prefix>.MX", then, if not found, for "MX".
281 // If still not found, use proton mass.
282
284 if( err ) return err;
285
286 FILE* f = OpenRunDBFile( date );
287 if( !f ) return kFileError;
288
289 if ( fMX <= 0.0 ) {
290 TString name(fPrefix), tag("MX"); name += tag;
291 Int_t st = LoadDBvalue( f, date, name.Data(), fMX );
292 if( st )
293 LoadDBvalue( f, date, tag.Data(), fMX );
294 if( fMX <= 0.0 ) fMX = 0.938;
295 }
296 fclose(f);
297 return 0;
298}
299
300//_____________________________________________________________________________
302{
303 if( !IsInit())
304 fMX = m;
305 else
306 PrintInitError("SetMX()");
307}
308
309//_____________________________________________________________________________
310void THaSecondaryKine::SetSpectrometer( const char* name )
311{
312 if( !IsInit())
314 else
315 PrintInitError("SetSpectrometer()");
316}
317
318//_____________________________________________________________________________
319void THaSecondaryKine::SetPrimary( const char* name )
320{
321 if( !IsInit())
323 else
324 PrintInitError("SetPrimary()");
325}
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]
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="")
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
bool DataValid() const
void PrintInitError(const char *here)
virtual void Clear(Option_t *opt="")
const FourVect * GetQ() const
Double_t GetOmega() const
const FourVect * GetA1() const
const FourVect * GetA() const
const FourVect * GetP1() const
void SetPrimary(const char *name)
THaTrackingModule * fSpectro
void SetSpectrometer(const char *name)
void SetMX(Double_t m)
THaSecondaryKine(const char *name, const char *description, const char *secondary_spectro="", const char *primary_kine="", Double_t secondary_mass=0.0)
virtual Int_t Process(const THaEvData &)
virtual void Clear(Option_t *opt="")
THaPrimaryKine * fPrimary
virtual Int_t DefineVariables(EMode mode=kDefine)
virtual Int_t ReadRunDatabase(const TDatime &date)
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 M() const
Double_t P() const
void Boost(const TVector3 &)
Double_t Theta() const
Double_t Phi() const
Double_t E() const
void SetVectM(const TVector3 &spatial, Double_t mass)
TRotation & SetZAxis(const TVector3 &axis)
TRotation & Invert()
const char * Data() const
Double_t Z() const
Double_t Phi() const
Double_t Y() const
Double_t X() const
Double_t Mag() const
Double_t Theta() const
double beta(double x, double y)
LVector boost(const LVector &v, const BoostVector &b)
STL namespace.
TMarker m
ClassImp(TPyArg)