Hall C ROOT/C++ Analyzer (hcana)
Loading...
Searching...
No Matches
THcPrimaryKine.cxx
Go to the documentation of this file.
1
8#include "THcPrimaryKine.h"
10#include "THcGlobals.h"
11#include "THcParmList.h"
12#include "THaRunBase.h"
13#include "THaRunParameters.h"
14#include "THaBeam.h"
15#include "VarDef.h"
16#include "TMath.h"
17#include <cstring>
18#include <cstdio>
19#include <cstdlib>
20#include <iostream>
21
22using namespace std;
23
25
26//_____________________________________________________________________________
27THcPrimaryKine::THcPrimaryKine( const char* name, const char* description,
28 const char* spectro, Double_t particle_mass,
29 Double_t target_mass ) :
30 THaPhysicsModule(name,description), fM(particle_mass),
31 fMA(target_mass), fSpectroName(spectro), fSpectro(NULL), fBeam(NULL)
32{
33 // Standard constructor. Must specify particle mass. Incident particles
34 // are assumed to be along z_lab.
35
36}
37
38//_____________________________________________________________________________
39THcPrimaryKine::THcPrimaryKine( const char* name, const char* description,
40 const char* spectro, const char* beam,
41 Double_t target_mass )
42 : THaPhysicsModule(name,description), fM(-1.0), fMA(target_mass),
43 fSpectroName(spectro), fBeamName(beam), fSpectro(NULL), fBeam(NULL)
44{
45 // Constructor with specification of optional beam module.
46 // Particle mass will normally come from the beam module.
47
48}
49//_____________________________________________________________________________
50
51//_____________________________________________________________________________
53{
54 // Destructor
55
57}
58
59//_____________________________________________________________________________
61{
62 // Clear event-by-event variables.
63
66 = fThetaQ = fPhiQ = kBig;
67 // Clear the 4-vectors
69 fP1 = fA = fA1 = fQ = fP0;
70
71}
72
73//_____________________________________________________________________________
75{
76 // Define/delete global variables.
77
78 if( mode == kDefine && fIsSetup ) return kOK;
79 fIsSetup = ( mode == kDefine );
80
81 RVarDef vars[] = {
82 { "Q2", "4-momentum transfer squared (GeV^2)", "fQ2" },
83 { "omega", "Energy transfer (GeV)", "fOmega" },
84 { "W2", "Invariant mass (GeV^2) for Mp ", "fW2" },
85 { "W", "Sqrt(Invariant mass) for Mp ", "fW" },
86 { "x_bj", "Bjorken x", "fXbj" },
87 { "scat_ang_rad", "Scattering angle (rad)", "fScatAngle" },
88 { "scat_ang_deg", "Scattering angle (deg)", "fScatAngle_deg" },
89 { "epsilon", "Virtual photon polarization factor", "fEpsilon" },
90 { "q3m", "Magnitude of 3-momentum transfer", "fQ3mag" },
91 { "th_q", "Theta of 3-momentum vector (rad)", "fThetaQ" },
92 { "ph_q", "Phi of 3-momentum vector (rad)", "fPhiQ" },
93 { "nu", "Energy transfer (GeV)", "fOmega" },
94 { "q_x", "x-cmp of Photon vector in the lab", "fQ.X()" },
95 { "q_y", "y-cmp of Photon vector in the lab", "fQ.Y()" },
96 { "q_z", "z-cmp of Photon vector in the lab", "fQ.Z()" },
97 { 0 }
98 };
99 return DefineVarsFromList( vars, mode );
100}
101
102//_____________________________________________________________________________
104{
105 // Initialize the module.
106 // Locate the spectrometer apparatus named in fSpectroName and save
107 // pointer to it.
108
109 fSpectro = dynamic_cast<THcHallCSpectrometer*>
110 ( FindModule( fSpectroName.Data(), "THcHallCSpectrometer"));
111 if( !fSpectro ) {
113 return fStatus;
114 }
115
116 // Optional beam apparatus
117 if( fBeamName.Length() > 0 ) {
118 fBeam = dynamic_cast<THaBeamModule*>
119 ( FindModule( fBeamName.Data(), "THaBeamModule") );
120 if( !fBeam ) {
122 return fStatus;
123 }
124 if( fM <= 0.0 )
125 fM = fBeam->GetBeamInfo()->GetM();
126 }
127
128 // Standard initialization. Calls this object's DefineVariables().
129 if( THaPhysicsModule::Init( run_time ) != kOK )
130 return fStatus;
131
132 if( fM <= 0.0 ) {
133 Error( Here("Init"), "Particle mass not defined. Module "
134 "initialization failed" );
136 }
137 return fStatus;
138}
139
140//_____________________________________________________________________________
142{
143 // Calculate electron kinematics for the Golden Track of the spectrometer
144 if( !IsOK() || !gHaRun ) return -1;
145
147 if( !trkifo || !trkifo->IsOK() ) return 1;
148 // Determine 4-momentum of incident particle.
149 // If a beam module given, use it to get the beam momentum. This
150 // module may apply corrections for beam energy loss, variations, etc.
151
152 Double_t xptar = trkifo->GetTheta() + fOopCentralOffset;
153 TVector3 pvect;
154 fSpectro->TransportToLab(trkifo->GetP(), xptar, trkifo->GetPhi(), pvect);
155
156 if( fBeam ) {
158 } else {
159 // If no beam given, assume beam along z_lab
161 fP0.SetXYZM( 0.0, 0.0, p_in, fM );
162 }
163
164 fP1.SetVectM( pvect, fM );
165 fA.SetXYZM( 0.0, 0.0, 0.0, fMA ); // Assume target at rest
166
167 // proton mass (for x_bj)
168 const Double_t Mp = 0.93827;
169 fMp.SetXYZM( 0.0, 0.0, 0.0, Mp ); // Assume target at rest
170
171
172 // Standard electron kinematics
173 fQ = fP0 - fP1; // cqx, cqy, cqz, omega
174 fQ2 = -fQ.M2();
175 fQ3mag = fQ.P();
176 fOmega = fQ.E();
177 // cqxzabs = TMath::Sqrt(fQ.X()*fQ.X() + fQ.Y()*fQ.Y());
178 fA1 = fA + fQ;
179 // fW2 = fA1.M2();
180 fMp1 = fMp + fQ;
181 fW2 = fMp1.M2();
182 if (fW2>0) fW = TMath::Sqrt(fW2);
183 fScatAngle = fP0.Angle( fP1.Vect() );
184 fEpsilon = 1.0 / ( 1.0 + 2.0*fQ3mag*fQ3mag/fQ2*
185 TMath::Power( TMath::Tan(fScatAngle/2.0), 2.0 ));
187 fThetaQ = fQ.Theta();
188 fPhiQ = fQ.Phi();
189 fXbj = fQ2/(2.0*Mp*fOmega);
190
191 fDataValid = true;
192 return 0;
193}
194
196{
197
198#ifdef WITH_DEBUG
199 cout << "In THcPrimaryKine::ReadDatabase()" << endl;
200#endif
201
202 char prefix[2];
203
204 prefix[0] = tolower(GetName()[0]);
205 prefix[1] = '\0';
206
207 fOopCentralOffset = 0.0;
208 DBRequest list[]={
209 {"gtargmass_amu", &fMA_amu, kDouble, 0, 1},
210 {Form("%s_oopcentral_offset",prefix), &fOopCentralOffset,kDouble, 0, 1},
211 {0}
212 };
213 gHcParms->LoadParmValues((DBRequest*)&list);
214 //
215 fMA= fMA_amu*931.5/1000.;
216 return kOK;
217}
218
219//_____________________________________________________________________________
221{
222 if( !IsInit())
223 fM = m;
224 else
225 PrintInitError("SetMass()");
226}
227
228//_____________________________________________________________________________
230{
231 if( !IsInit())
232 fMA = m;
233 else
234 PrintInitError("SetTargetMass()");
235}
236
237//_____________________________________________________________________________
238void THcPrimaryKine::SetSpectrometer( const char* name )
239{
240 if( !IsInit())
241 fSpectroName = name;
242 else
243 PrintInitError("SetSpectrometer()");
244}
245
246//_____________________________________________________________________________
247void THcPrimaryKine::SetBeam( const char* name )
248{
249 if( !IsInit())
250 fBeamName = name;
251 else
252 PrintInitError("SetBeam()");
253}
int Int_t
const Data_t kBig
double Double_t
const char Option_t
Option_t Option_t TPoint TPoint const char mode
R__EXTERN class THaRunBase * gHaRun
ClassImp(VDC::AnalyticTTDConv) using namespace std
R__EXTERN class THcParmList * gHcParms
Definition THcGlobals.h:11
char * Form(const char *fmt,...)
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
Bool_t IsInit() const
Double_t GetM() const
const TVector3 & GetPvect() const
THaBeamInfo * GetBeamInfo()
void PrintInitError(const char *here)
virtual void Clear(Option_t *opt="")
THaRunParameters * GetParameters() const
Double_t GetBeamP() const
virtual void TransportToLab(Double_t p, Double_t th, Double_t ph, TVector3 &pvect) const
Double_t GetTheta() const
Double_t GetP() const
Double_t GetPhi() const
Bool_t IsOK() const
THaTrackInfo * GetTrackInfo()
A standard Hall C spectrometer apparatus.
Int_t LoadParmValues(const DBRequest *list, const char *prefix="")
Retrieve parameter values from the parameter cache.
Class for the Calculate kinematics of scattering of the primary (beam) particle. These are usually th...
virtual Int_t DefineVariables(EMode mode=kDefine)
void SetMass(Double_t m)
THaBeamModule * fBeam
THcPrimaryKine(const char *name, const char *description, const char *spectro="", Double_t particle_mass=0.0, Double_t target_mass=0.0)
virtual Int_t ReadDatabase(const TDatime &date)
THcHallCSpectrometer * fSpectro
Double_t fScatAngle_deg
Double_t fScatAngle
void SetSpectrometer(const char *name)
virtual void Clear(Option_t *opt="")
void SetBeam(const char *name)
virtual ~THcPrimaryKine()
void SetTargetMass(Double_t m)
Double_t fOopCentralOffset
virtual Int_t Process(const THaEvData &)
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)
const char * GetName() const override
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 Sqrt(Double_t x)
Double_t Tan(Double_t)
constexpr Double_t RadToDeg()
STL namespace.
TMarker m