QwAnalysis
QwEvent Class Reference

Contains a tracked event, i.e. all information from hits to tracks. More...

#include <QwEvent.h>

+ Inheritance diagram for QwEvent:
+ Collaboration diagram for QwEvent:

Public Member Functions

 QwEvent ()
 Default constructor. More...
 
virtual ~QwEvent ()
 Virtual destructor. More...
 
const QwEventHeaderGetEventHeader () const
 
void SetEventHeader (const QwEventHeader &eventheader)
 
void SetEventHeader (const QwEventHeader *eventheader)
 
void LoadBeamProperty (const TString &map)
 Load the beam properties from a map file. More...
 
double EnergyLossHydrogen (const double vertex_z)
 Calculate the energy loss in the hydrogen target. More...
 
void Clear (Option_t *option="")
 
void Reset (Option_t *option="")
 
void CalculateKinematics (const QwTrack *track)
 Calculate the kinematic variables for a given track. More...
 
void Print (Option_t *option="") const
 Print the event. More...
 
const TVector3 & GetVertexPosition () const
 
const TVector3 & GetVertexMomentum () const
 
Double_t GetScatteringAngle () const
 
Double_t GetScatteringVertexZ () const
 
 ClassDef (QwEvent, 3)
 
Hit list maintenance for output to ROOT files
QwHitCreateNewHit ()
 Create a new hit. More...
 
void AddHit (const QwHit *hit)
 Add an existing hit as a copy. More...
 
void ClearHits (Option_t *option="")
 Clear the list of hits. More...
 
void ResetHits (Option_t *option="")
 Reset the list of hits. More...
 
void AddHitContainer (const QwHitContainer *hitlist)
 Add the hits in a hit container as a copy. More...
 
QwHitContainerGetHitContainer ()
 Get the list of hits as a hit container. More...
 
Int_t GetNumberOfHits () const
 Get the number of hits. More...
 
const std::vector< QwHit * > & GetListOfHits () const
 Get the list of hits. More...
 
const QwHitGetHit (const int hit) const
 Get the specified hit. More...
 
void PrintHits (Option_t *option="") const
 Print the list of hits. More...
 
Tree line list maintenance for output to ROOT files
QwTreeLineCreateNewTreeLine ()
 Create a new tree line. More...
 
void AddTreeLine (const QwTreeLine *treeline)
 Add an existing tree line as a copy. More...
 
void AddTreeLineList (const QwTreeLine *treelinelist)
 Add a list of existing tree lines as a copy. More...
 
void AddTreeLineList (const std::vector< QwTreeLine * > &treelinelist)
 Add a list of existing tree lines as a copy. More...
 
void ClearTreeLines (Option_t *option="")
 Clear the list of tree lines. More...
 
void ResetTreeLines (Option_t *option="")
 Reset the list of tree lines. More...
 
Int_t GetNumberOfTreeLines () const
 Get the number of tree lines. More...
 
const std::vector< QwTreeLine * > & GetListOfTreeLines () const
 Get the list of tree lines. More...
 
const std::vector< QwTreeLine * > GetListOfTreeLines (EQwRegionID region, EQwDirectionID direction) const
 Get the list of tree lines in region and direction. More...
 
const QwTreeLineGetTreeLine (const int tl) const
 Get the specified tree line. More...
 
void PrintTreeLines (Option_t *option="") const
 Print the list of tree lines. More...
 
Partial track list maintenance for output to ROOT files
QwPartialTrackCreateNewPartialTrack ()
 Create a new partial track. More...
 
void AddPartialTrack (const QwPartialTrack *partialtrack)
 Add an existing partial track as a copy. More...
 
void AddPartialTrackList (const std::vector< QwPartialTrack * > &partialtracklist)
 Add a list of existing partial tracks as a copy. More...
 
void ClearPartialTracks (Option_t *option="")
 Clear the list of partial tracks. More...
 
void ResetPartialTracks (Option_t *option="")
 Reset the list of partial tracks. More...
 
Int_t GetNumberOfPartialTracks () const
 Get the number of partial tracks. More...
 
const std::vector
< QwPartialTrack * > 
GetListOfPartialTracks (EQwRegionID region, EQwDetectorPackage package) const
 Get the list of partial tracks in region and direction. More...
 
const std::vector
< QwPartialTrack * > & 
GetListOfPartialTracks () const
 Get the list of partial tracks. More...
 
const QwPartialTrackGetPartialTrack (const int pt) const
 Get the specified partial track. More...
 
void PrintPartialTracks (Option_t *option="") const
 Print the list of partial tracks. More...
 
Track list maintenance for output to ROOT files
QwTrackCreateNewTrack ()
 Create a new track. More...
 
void AddTrack (const QwTrack *track)
 Add an existing track as a copy. More...
 
void AddTrackList (const std::vector< QwTrack * > &tracklist)
 Add a list of existing partial tracks as a copy. More...
 
void ClearTracks (Option_t *option="")
 Clear the list of tracks. More...
 
void ResetTracks (Option_t *option="")
 Reset the list of tracks. More...
 
Int_t GetNumberOfTracks () const
 Get the number of tracks. More...
 
const std::vector< QwTrack * > & GetListOfTracks () const
 Get the list of tracks. More...
 
const QwTrackGetTrack (const int t) const
 Get the specified track. More...
 
void PrintTracks (Option_t *option="") const
 Print the list of tracks. More...
 
- Public Member Functions inherited from QwObjectCounter< QwEvent >
 QwObjectCounter ()
 Default constructor. More...
 
 QwObjectCounter (const QwObjectCounter &)
 Copy constructor. More...
 
virtual ~QwObjectCounter ()
 Destructor. More...
 

Data Fields

QwEventHeaderfEventHeader
 
Int_t fNQwHits
 Number of QwHits in the array. More...
 
std::vector< QwHit * > fQwHits
 Array of QwHits. More...
 
Int_t fNQwTreeLines
 Number of QwTreeLines in the array. More...
 
std::vector< QwTreeLine * > fQwTreeLines
 Array of QwTreeLines. More...
 
Int_t fNQwPartialTracks
 Number of QwPartialTracks in the array. More...
 
std::vector< QwPartialTrack * > fQwPartialTracks
 Array of QwPartialTracks. More...
 
Int_t fNQwTracks
 Number of QwTracks in the array. More...
 
std::vector< QwTrack * > fQwTracks
 Array of QwTracks. More...
 
Double_t fPrimaryQ2
 < Momentum transfer Q^2 assuming elastic scattering with hydrogen energy loss More...
 
QwTreeLinefTreeLine [kNumPackages][kNumRegions][kNumDirections]
 
Main detector light yield
std::vector< float > fMD_LeftNbOfPEs
 
std::vector< float > fMD_RightNbOfPEs
 
std::vector< float > fMD_TotalNbOfPEs
 
Generic kinematic variables
double fHydrogenEnergyLoss
 Pre-scattering target energy loss assuming LH2 target. More...
 
double fScatteringAngle
 Scattering angle. More...
 
double fScatteringVertexZ
 Scattering vertex z position. More...
 
double fScatteringVertexR
 Scattering vertex radial distance. More...
 
double fCrossSection
 
double fPreScatteringEnergy
 
double fOriginVertexEnergy
 
TVector3 fVertexPosition
 Vertex position. More...
 
TVector3 fVertexMomentum
 Vertex momentum. More...
 
Kinematics under various assumptions
QwKinematics fKin
 Inclusive scattering. More...
 
QwKinematics fKinWithLoss
 Scattering with hydrogen energy loss. More...
 
QwKinematics fKinElastic
 Scattering assuming elastic reaction. More...
 
QwKinematics fKinElasticWithLoss
 Scattering assuming elastic reaction and hydrogen energy loss. More...
 

Static Public Attributes

static double fBeamEnergy = 1.165 * Qw::GeV
 Electron beam energy. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from QwObjectCounter< QwEvent >
static size_t GetObjectsCreated ()
 Get number of objects ever created. More...
 
static size_t GetObjectsAlive ()
 Get number of objects still alive. More...
 

Detailed Description

Contains a tracked event, i.e. all information from hits to tracks.

A QwEvent contains all event information, from hits over partial track to complete tracks. It serves as the final product of the tracking code.

Definition at line 156 of file QwEvent.h.

Constructor & Destructor Documentation

QwEvent::QwEvent ( )

Default constructor.

Definition at line 20 of file QwEvent.cc.

References Clear(), fEventHeader, fNQwHits, fNQwPartialTracks, fNQwTracks, fNQwTreeLines, fTreeLine, kNumDirections, kNumPackages, and kNumRegions.

23  fPrimaryQ2(0.0)
24 {
25  // Reset the event header
26  fEventHeader = 0;
27 
28  fNQwHits = 0;
29  fNQwTreeLines = 0;
31  fNQwTracks = 0;
32 
33  // Loop over all pointer objects
34  for (int i = 0; i < kNumPackages; ++i) {
35  for (int j = 0; j < kNumRegions; ++j) {
36  // Null the treeline pointers
37  for (int k = 0; k < kNumDirections; ++k) {
38  fTreeLine[i][j][k] = 0;
39  } // end of loop over directions
40  } // end of loop over regions
41  } // end of loop over packages
42 
43  // Clear the event
44  Clear();
45 }
Double_t fPrimaryQ2
&lt; Momentum transfer Q^2 assuming elastic scattering with hydrogen energy loss
Definition: QwEvent.h:355
Int_t fNQwHits
Number of QwHits in the array.
Definition: QwEvent.h:168
double fScatteringAngle
Scattering angle.
Definition: QwEvent.h:336
double fScatteringVertexR
Scattering vertex radial distance.
Definition: QwEvent.h:338
double fCrossSection
Definition: QwEvent.h:339
Int_t fNQwTreeLines
Number of QwTreeLines in the array.
Definition: QwEvent.h:172
Int_t fNQwTracks
Number of QwTracks in the array.
Definition: QwEvent.h:180
double fScatteringVertexZ
Scattering vertex z position.
Definition: QwEvent.h:337
QwEventHeader * fEventHeader
Definition: QwEvent.h:161
double fOriginVertexEnergy
Definition: QwEvent.h:341
double fPreScatteringEnergy
Definition: QwEvent.h:340
void Clear(Option_t *option="")
Definition: QwEvent.cc:77
Int_t fNQwPartialTracks
Number of QwPartialTracks in the array.
Definition: QwEvent.h:176
QwTreeLine * fTreeLine[kNumPackages][kNumRegions][kNumDirections]
Definition: QwEvent.h:358

+ Here is the call graph for this function:

QwEvent::~QwEvent ( )
virtual

Virtual destructor.

Definition at line 48 of file QwEvent.cc.

References Clear(), fEventHeader, fTreeLine, kNumDirections, kNumPackages, kNumRegions, and QwTreeLine::next.

49 {
50  // Clear the event
51  Clear();
52 
53  // Loop over all allocated objects
54  for (int i = 0; i < kNumPackages; ++i) {
55 
56  for (int j = 0; j < kNumRegions; ++j) {
57 
58  // Delete all those treelines
59  for (int k = 0; k < kNumDirections; ++k) {
60  QwTreeLine* tl = fTreeLine[i][j][k];
61  while (tl) {
62  QwTreeLine* tl_next = tl->next;
63  delete tl;
64  tl = tl_next;
65  }
66  } // end of loop over directions
67  } // end of loop over regions
68  } // end of loop over packages
69 
70 
71  // Delete the event header
72  if (fEventHeader) delete fEventHeader;
73 }
QwEventHeader * fEventHeader
Definition: QwEvent.h:161
One-dimensional (u, v, or x) track stubs and associated hits.
Definition: QwTreeLine.h:51
void Clear(Option_t *option="")
Definition: QwEvent.cc:77
QwTreeLine * next
Definition: QwTreeLine.h:231
QwTreeLine * fTreeLine[kNumPackages][kNumRegions][kNumDirections]
Definition: QwEvent.h:358

+ Here is the call graph for this function:

Member Function Documentation

void QwEvent::AddHit ( const QwHit hit)

Add an existing hit as a copy.

Definition at line 281 of file QwEvent.cc.

References QwLog::endl(), fNQwHits, fQwHits, and QwError.

Referenced by AddHitContainer(), and CreateNewHit().

282 {
283  if (hit) fQwHits.push_back(new QwHit(hit));
284  else QwError << "trying to add null hit" << QwLog::endl;
285  fNQwHits++;
286 }
Int_t fNQwHits
Number of QwHits in the array.
Definition: QwEvent.h:168
static std::ostream & endl(std::ostream &)
End of the line.
Definition: QwLog.cc:299
Hit structure uniquely defining each hit.
Definition: QwHit.h:43
std::vector< QwHit * > fQwHits
Array of QwHits.
Definition: QwEvent.h:169
#define QwError
Predefined log drain for errors.
Definition: QwLog.h:40

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void QwEvent::AddHitContainer ( const QwHitContainer hitlist)

Add the hits in a hit container as a copy.

Definition at line 322 of file QwEvent.cc.

References AddHit().

Referenced by QwTreeEventBuffer::GetSpecificEvent().

323 {
324  for (QwHitContainer::const_iterator hit = hitlist->begin();
325  hit != hitlist->end(); hit++) {
326  const QwHit* p = &(*hit);
327  AddHit(p);
328  }
329 }
void AddHit(const QwHit *hit)
Add an existing hit as a copy.
Definition: QwEvent.cc:281
Hit structure uniquely defining each hit.
Definition: QwHit.h:43

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void QwEvent::AddPartialTrack ( const QwPartialTrack partialtrack)

Add an existing partial track as a copy.

Definition at line 442 of file QwEvent.cc.

References fNQwPartialTracks, and fQwPartialTracks.

Referenced by AddPartialTrackList(), CreateNewPartialTrack(), and QwTreeEventBuffer::GetSpecificEvent().

443 {
444  fQwPartialTracks.push_back(new QwPartialTrack(partialtrack));
446 }
std::vector< QwPartialTrack * > fQwPartialTracks
Array of QwPartialTracks.
Definition: QwEvent.h:177
Int_t fNQwPartialTracks
Number of QwPartialTracks in the array.
Definition: QwEvent.h:176
Contains the straight part of a track in one region only.

+ Here is the caller graph for this function:

void QwEvent::AddPartialTrackList ( const std::vector< QwPartialTrack * > &  partialtracklist)

Add a list of existing partial tracks as a copy.

Definition at line 449 of file QwEvent.cc.

References AddPartialTrack().

450 {
451  for (size_t i = 0; i < partialtracklist.size(); i++)
452  if (partialtracklist[i]->IsValid())
453  AddPartialTrack(partialtracklist[i]);
454 }
void AddPartialTrack(const QwPartialTrack *partialtrack)
Add an existing partial track as a copy.
Definition: QwEvent.cc:442

+ Here is the call graph for this function:

void QwEvent::AddTrack ( const QwTrack track)

Add an existing track as a copy.

Definition at line 515 of file QwEvent.cc.

References fNQwTracks, and fQwTracks.

Referenced by AddTrackList(), and CreateNewTrack().

516 {
517  fQwTracks.push_back(new QwTrack(track));
518  ++fNQwTracks;
519 }
std::vector< QwTrack * > fQwTracks
Array of QwTracks.
Definition: QwEvent.h:181
Int_t fNQwTracks
Number of QwTracks in the array.
Definition: QwEvent.h:180
Contains the complete track as a concatenation of partial tracks.
Definition: QwTrack.h:30

+ Here is the caller graph for this function:

void QwEvent::AddTrackList ( const std::vector< QwTrack * > &  tracklist)

Add a list of existing partial tracks as a copy.

Definition at line 522 of file QwEvent.cc.

References AddTrack().

523 {
524  for (size_t i = 0; i < tracklist.size(); i++)
525  //if (tracklist[i]->IsValid()) // TODO
526  AddTrack(tracklist[i]);
527 }
void AddTrack(const QwTrack *track)
Add an existing track as a copy.
Definition: QwEvent.cc:515

+ Here is the call graph for this function:

void QwEvent::AddTreeLine ( const QwTreeLine treeline)

Add an existing tree line as a copy.

Definition at line 354 of file QwEvent.cc.

References fNQwTreeLines, and fQwTreeLines.

Referenced by AddTreeLineList(), CreateNewTreeLine(), and QwTreeEventBuffer::GetSpecificEvent().

355 {
356  fQwTreeLines.push_back(new QwTreeLine(treeline));
357  fNQwTreeLines++;
358 }
Int_t fNQwTreeLines
Number of QwTreeLines in the array.
Definition: QwEvent.h:172
std::vector< QwTreeLine * > fQwTreeLines
Array of QwTreeLines.
Definition: QwEvent.h:173
One-dimensional (u, v, or x) track stubs and associated hits.
Definition: QwTreeLine.h:51

+ Here is the caller graph for this function:

void QwEvent::AddTreeLineList ( const QwTreeLine treelinelist)

Add a list of existing tree lines as a copy.

Definition at line 361 of file QwEvent.cc.

References AddTreeLine(), and QwTreeLine::next.

362 {
363  for (const QwTreeLine *treeline = treelinelist;
364  treeline; treeline = treeline->next){
365  if (treeline->IsValid()){
366  AddTreeLine(treeline);
367  }
368  }
369 }
void AddTreeLine(const QwTreeLine *treeline)
Add an existing tree line as a copy.
Definition: QwEvent.cc:354
One-dimensional (u, v, or x) track stubs and associated hits.
Definition: QwTreeLine.h:51
QwTreeLine * next
Definition: QwTreeLine.h:231

+ Here is the call graph for this function:

void QwEvent::AddTreeLineList ( const std::vector< QwTreeLine * > &  treelinelist)

Add a list of existing tree lines as a copy.

Definition at line 372 of file QwEvent.cc.

References AddTreeLine().

373 {
374  for (size_t i = 0; i < treelinelist.size(); i++)
375  if (treelinelist[i]->IsValid())
376  AddTreeLine(treelinelist[i]);
377 }
void AddTreeLine(const QwTreeLine *treeline)
Add an existing tree line as a copy.
Definition: QwEvent.cc:354

+ Here is the call graph for this function:

void QwEvent::CalculateKinematics ( const QwTrack track)

Calculate the kinematic variables for a given track.

Calculate the kinematic variables for a given track

Parameters
trackReconstructed track

Definition at line 161 of file QwEvent.cc.

References Qw::cm, Qw::deg, QwLog::endl(), EnergyLossHydrogen(), fBeamEnergy, QwTrack::fFront, fHydrogenEnergyLoss, fKin, fKinElastic, fKinElasticWithLoss, fKinWithLoss, QwTrack::fMomentum, QwKinematics::fNu, QwKinematics::fP0, QwKinematics::fPp, fPrimaryQ2, QwKinematics::fQ2, fScatteringAngle, fScatteringVertexR, fScatteringVertexZ, fVertexMomentum, fVertexPosition, QwKinematics::fW2, QwKinematics::fX, QwKinematics::fY, QwPartialTrack::GetMomentumDirection(), QwPartialTrack::GetMomentumDirectionTheta(), QwPartialTrack::GetPosition(), QwPartialTrack::GetVertexZ(), Qw::GeV, Qw::GeV2, Qw::MeV, Qw::Mp, and QwWarning.

162 {
163  if (! track) {
164  QwWarning << "QwEvent::CalculateKinematics: "
165  << "called with null track pointer." << QwLog::endl;
166  return;
167  }
168 
169  if (! (track->fFront)) {
170  QwWarning << "QwEvent::CalculateKinematics: "
171  << "full track with null front track pointer." << QwLog::endl;
172  return;
173  }
174 
175  // Beam energy
176  Double_t energy = fBeamEnergy;
177 
178  // Longitudinal vertex in target from front partial track
179  Double_t vertex_z = track->fFront->GetVertexZ();
180  fVertexPosition = track->fFront->GetPosition(vertex_z);
184 
185  // Scattering angle from front partial track
186  Double_t theta = track->fFront->GetMomentumDirectionTheta();
187  Double_t cos_theta = cos(theta);
188  fScatteringAngle = theta / Qw::deg;
189 
190  // Prescattering energy loss
191  Double_t pre_loss = EnergyLossHydrogen(vertex_z);
192  fHydrogenEnergyLoss = pre_loss;
193 
194  // Generic scattering without energy loss
195  Double_t P0 = energy;
196  Double_t PP = track->fMomentum;
197  Double_t Q2 = 2.0 * P0 * PP * (1 - cos_theta);
198  fKin.fP0 = P0 / Qw::GeV;
199  fKin.fPp = PP / Qw::GeV;
200  fKin.fQ2 = Q2 / Qw::GeV2;
201  fKin.fNu = (P0 - PP) / Qw::GeV;
202  fKin.fW2 = (Qw::Mp * Qw::Mp + 2.0 * Qw::Mp * (P0 - PP) - Q2) / Qw::GeV2;
203  fKin.fX = Q2 / (2.0 * Qw::Mp * (P0 - PP));
204  fKin.fY = (P0 - PP) / P0;
205 
206  // Generic scattering with energy loss
207  P0 = energy - pre_loss;
208  PP = track->fMomentum;
209  Q2 = 2.0 * P0 * PP * (1 - cos_theta);
210  fKinWithLoss.fP0 = P0 / Qw::GeV;
211  fKinWithLoss.fPp = PP / Qw::GeV;
212  fKinWithLoss.fQ2 = Q2 / Qw::GeV2;
213  fKinWithLoss.fNu = (P0 - PP) / Qw::GeV;
214  fKinWithLoss.fW2 = (Qw::Mp * Qw::Mp + 2.0 * Qw::Mp * (P0 - PP) - Q2) / Qw::GeV2;
215  fKinWithLoss.fX = Q2 / (2.0 * Qw::Mp * (P0 - PP));
216  fKinWithLoss.fY = (P0 - PP) / P0;
217 
218  // Elastic scattering without energy loss
219  P0 = energy * Qw::MeV;
220  PP = Qw::Mp * P0 / (Qw::Mp + P0 * (1 - cos_theta));
221  Q2 = 2.0 * P0 * PP * (1 - cos_theta);
222  fKinElastic.fP0 = P0 / Qw::GeV;
223  fKinElastic.fPp = PP / Qw::GeV;
224  fKinElastic.fQ2 = Q2 / Qw::GeV2;
225  fKinElastic.fNu = (P0 - PP) / Qw::GeV;
226  fKinElastic.fW2 = (Qw::Mp * Qw::Mp + 2.0 * Qw::Mp * (P0 - PP) - Q2) / Qw::GeV2;
227  fKinElastic.fX = Q2 / (2.0 * Qw::Mp * (P0 - PP));
228  fKinElastic.fY = (P0 - PP) / P0;
229 
230  // Elastic scattering with energy loss
231  P0 = (energy - pre_loss) * Qw::MeV;
232  PP = Qw::Mp * P0 / (Qw::Mp + P0 * (1 - cos_theta));
233  Q2 = 2.0 * P0 * PP * (1 - cos_theta);
237  fKinElasticWithLoss.fNu = (P0 - PP) / Qw::GeV;
238  fKinElasticWithLoss.fW2 = (Qw::Mp * Qw::Mp + 2.0 * Qw::Mp * (P0 - PP) - Q2) / Qw::GeV2;
239  fKinElasticWithLoss.fX = Q2 / (2.0 * Qw::Mp * (P0 - PP));
240  fKinElasticWithLoss.fY = (P0 - PP) / P0;
241 
243 }
Double_t fPrimaryQ2
&lt; Momentum transfer Q^2 assuming elastic scattering with hydrogen energy loss
Definition: QwEvent.h:355
TVector3 fVertexPosition
Vertex position.
Definition: QwEvent.h:342
double EnergyLossHydrogen(const double vertex_z)
Calculate the energy loss in the hydrogen target.
Definition: QwEvent.cc:127
double fScatteringAngle
Scattering angle.
Definition: QwEvent.h:336
static const double deg
Definition: QwUnits.h:106
double fW2
Invariant mass squared of the recoiling target system.
Definition: QwEvent.h:138
double fScatteringVertexR
Scattering vertex radial distance.
Definition: QwEvent.h:338
static double fBeamEnergy
Electron beam energy.
Definition: QwEvent.h:320
Double_t GetMomentumDirectionTheta() const
Return the theta angle.
QwKinematics fKin
Inclusive scattering.
Definition: QwEvent.h:348
TVector3 fVertexMomentum
Vertex momentum.
Definition: QwEvent.h:343
const TVector3 GetPosition(const double z) const
Return the vertex at position z.
QwKinematics fKinElasticWithLoss
Scattering assuming elastic reaction and hydrogen energy loss.
Definition: QwEvent.h:351
double fP0
Incoming momentum .
Definition: QwEvent.h:133
double fQ2
Four-momentum transfer squared .
Definition: QwEvent.h:137
double fNu
Energy loss of the electron .
Definition: QwEvent.h:139
const TVector3 GetMomentumDirection() const
Return the direction.
double fPp
Outgoing momentum .
Definition: QwEvent.h:136
QwKinematics fKinWithLoss
Scattering with hydrogen energy loss.
Definition: QwEvent.h:349
static const double GeV2
Definition: QwUnits.h:97
double fScatteringVertexZ
Scattering vertex z position.
Definition: QwEvent.h:337
Double_t GetVertexZ() const
Get the vertex position.
double fY
Fractional energy loss .
Definition: QwEvent.h:141
static const double GeV
Definition: QwUnits.h:95
static const double Mp
Mass of the proton.
Definition: QwUnits.h:119
double fHydrogenEnergyLoss
Pre-scattering target energy loss assuming LH2 target.
Definition: QwEvent.h:335
static std::ostream & endl(std::ostream &)
End of the line.
Definition: QwLog.cc:299
QwPartialTrack * fFront
Front partial track.
Definition: QwTrack.h:147
double fX
Bjorken-x scaling variable .
Definition: QwEvent.h:140
static const double MeV
Definition: QwUnits.h:94
#define QwWarning
Predefined log drain for warnings.
Definition: QwLog.h:45
double fMomentum
Spectrometer momentum.
Definition: QwTrack.h:100
QwKinematics fKinElastic
Scattering assuming elastic reaction.
Definition: QwEvent.h:350
static const double cm
Length units: base unit is mm.
Definition: QwUnits.h:61

+ Here is the call graph for this function:

QwEvent::ClassDef ( QwEvent  ,
 
)
void QwEvent::Clear ( Option_t *  option = "")

Definition at line 77 of file QwEvent.cc.

References ClearHits(), ClearPartialTracks(), ClearTracks(), and ClearTreeLines().

Referenced by QwEvent(), and ~QwEvent().

78 {
79  ClearHits(option);
80  ClearTreeLines(option);
81  ClearPartialTracks(option);
82  ClearTracks(option);
83 }
void ClearHits(Option_t *option="")
Clear the list of hits.
Definition: QwEvent.cc:289
void ClearTracks(Option_t *option="")
Clear the list of tracks.
Definition: QwEvent.cc:530
void ClearPartialTracks(Option_t *option="")
Clear the list of partial tracks.
Definition: QwEvent.cc:457
void ClearTreeLines(Option_t *option="")
Clear the list of tree lines.
Definition: QwEvent.cc:380

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void QwEvent::ClearHits ( Option_t *  option = "")

Clear the list of hits.

Definition at line 289 of file QwEvent.cc.

References fNQwHits, and fQwHits.

Referenced by Clear(), and ResetHits().

290 {
291  for (size_t i = 0; i < fQwHits.size(); i++)
292  delete fQwHits.at(i);
293  fQwHits.clear();
294  fNQwHits = 0;
295 }
Int_t fNQwHits
Number of QwHits in the array.
Definition: QwEvent.h:168
std::vector< QwHit * > fQwHits
Array of QwHits.
Definition: QwEvent.h:169

+ Here is the caller graph for this function:

void QwEvent::ClearPartialTracks ( Option_t *  option = "")

Clear the list of partial tracks.

Definition at line 457 of file QwEvent.cc.

References fNQwPartialTracks, and fQwPartialTracks.

Referenced by Clear(), and ResetPartialTracks().

458 {
459  for (size_t i = 0; i < fQwPartialTracks.size(); i++)
460  delete fQwPartialTracks.at(i);
461  fQwPartialTracks.clear();
462  fNQwPartialTracks = 0;
463 }
std::vector< QwPartialTrack * > fQwPartialTracks
Array of QwPartialTracks.
Definition: QwEvent.h:177
Int_t fNQwPartialTracks
Number of QwPartialTracks in the array.
Definition: QwEvent.h:176

+ Here is the caller graph for this function:

void QwEvent::ClearTracks ( Option_t *  option = "")

Clear the list of tracks.

Definition at line 530 of file QwEvent.cc.

References fNQwTracks, and fQwTracks.

Referenced by Clear(), and ResetTracks().

531 {
532  for (size_t i = 0; i < fQwTracks.size(); i++)
533  delete fQwTracks.at(i);
534  fQwTracks.clear();
535  fNQwTracks = 0;
536 }
std::vector< QwTrack * > fQwTracks
Array of QwTracks.
Definition: QwEvent.h:181
Int_t fNQwTracks
Number of QwTracks in the array.
Definition: QwEvent.h:180

+ Here is the caller graph for this function:

void QwEvent::ClearTreeLines ( Option_t *  option = "")

Clear the list of tree lines.

Definition at line 380 of file QwEvent.cc.

References fNQwTreeLines, and fQwTreeLines.

Referenced by Clear(), and ResetTreeLines().

381 {
382  for (size_t i = 0; i < fQwTreeLines.size(); i++) {
383  QwTreeLine* tl = fQwTreeLines.at(i);
384  delete tl;
385  }
386  fQwTreeLines.clear();
387  fNQwTreeLines = 0;
388 }
Int_t fNQwTreeLines
Number of QwTreeLines in the array.
Definition: QwEvent.h:172
std::vector< QwTreeLine * > fQwTreeLines
Array of QwTreeLines.
Definition: QwEvent.h:173
One-dimensional (u, v, or x) track stubs and associated hits.
Definition: QwTreeLine.h:51

+ Here is the caller graph for this function:

QwHit * QwEvent::CreateNewHit ( )

Create a new hit.

Definition at line 273 of file QwEvent.cc.

References AddHit().

274 {
275  QwHit* hit = new QwHit();
276  AddHit(hit);
277  return hit;
278 }
void AddHit(const QwHit *hit)
Add an existing hit as a copy.
Definition: QwEvent.cc:281
Hit structure uniquely defining each hit.
Definition: QwHit.h:43

+ Here is the call graph for this function:

QwPartialTrack * QwEvent::CreateNewPartialTrack ( )

Create a new partial track.

Definition at line 434 of file QwEvent.cc.

References AddPartialTrack().

435 {
436  QwPartialTrack* partialtrack = new QwPartialTrack();
437  AddPartialTrack(partialtrack);
438  return partialtrack;
439 }
void AddPartialTrack(const QwPartialTrack *partialtrack)
Add an existing partial track as a copy.
Definition: QwEvent.cc:442
Contains the straight part of a track in one region only.

+ Here is the call graph for this function:

QwTrack * QwEvent::CreateNewTrack ( )

Create a new track.

Definition at line 507 of file QwEvent.cc.

References AddTrack().

508 {
509  QwTrack* track = new QwTrack();
510  AddTrack(track);
511  return track;
512 }
void AddTrack(const QwTrack *track)
Add an existing track as a copy.
Definition: QwEvent.cc:515
Contains the complete track as a concatenation of partial tracks.
Definition: QwTrack.h:30

+ Here is the call graph for this function:

QwTreeLine * QwEvent::CreateNewTreeLine ( )

Create a new tree line.

Definition at line 346 of file QwEvent.cc.

References AddTreeLine().

347 {
348  QwTreeLine* treeline = new QwTreeLine();
349  AddTreeLine(treeline);
350  return treeline;
351 }
void AddTreeLine(const QwTreeLine *treeline)
Add an existing tree line as a copy.
Definition: QwEvent.cc:354
One-dimensional (u, v, or x) track stubs and associated hits.
Definition: QwTreeLine.h:51

+ Here is the call graph for this function:

double QwEvent::EnergyLossHydrogen ( const double  vertex_z)

Calculate the energy loss in the hydrogen target.

Calculate the energy loss in the hydrogen target

Parameters
vertex_zLongitudinal position of vertex (absolute coordinates)
Returns
Energy loss

Definition at line 127 of file QwEvent.cc.

References Qw::MeV.

Referenced by CalculateKinematics().

128 {
129  // Target parameters
130  double target_z_length = 34.35; // Target Length (cm)
131  double target_z_position= -652.67; // Target center position (cm) in Z
132  double target_z_front = target_z_position - 0.5 * target_z_length;
133  double target_z_back = target_z_position + 0.5 * target_z_length;
134 
135  // Energy loss in hydrogen target
136  double pre_loss = 0.0;
137  if (vertex_z < target_z_front) {
138  // Before target
139  pre_loss = 0.05 * Qw::MeV;
140 
141  } else if (vertex_z >= target_z_front && vertex_z <= target_z_back) {
142  // Inside target
143  double depth = vertex_z - target_z_front;
144  //pre_loss = 20.08 * Qw::MeV * depth/target_z_length; // linear approximation
145  pre_loss = 0.05 + depth*0.6618 - 0.003462*depth*depth; // quadratic fit approximation
146 
147  } else {
148  // After target
149  pre_loss = 18.7 * Qw::MeV; // averaged total Eloss, full target (excluding downstream Al window)
150  }
151 
152  return pre_loss;
153 }
static const double MeV
Definition: QwUnits.h:94

+ Here is the caller graph for this function:

const QwEventHeader* QwEvent::GetEventHeader ( ) const
inline

Definition at line 192 of file QwEvent.h.

References fEventHeader.

192 { return fEventHeader; };
QwEventHeader * fEventHeader
Definition: QwEvent.h:161
const QwHit * QwEvent::GetHit ( const int  hit) const

Get the specified hit.

Definition at line 304 of file QwEvent.cc.

References fQwHits, and GetNumberOfHits().

305 {
306  if (hit < 0 || hit > GetNumberOfHits()) return 0;
307  return fQwHits.at(hit);
308 }
Int_t GetNumberOfHits() const
Get the number of hits.
Definition: QwEvent.h:230
std::vector< QwHit * > fQwHits
Array of QwHits.
Definition: QwEvent.h:169

+ Here is the call graph for this function:

QwHitContainer * QwEvent::GetHitContainer ( )

Get the list of hits as a hit container.

Definition at line 332 of file QwEvent.cc.

References QwLog::endl(), fQwHits, and QwError.

Referenced by QwTreeEventBuffer::GetHitContainer().

333 {
334  QwHitContainer* hitlist = new QwHitContainer();
335  for (std::vector<QwHit*>::const_iterator hit = fQwHits.begin();
336  hit != fQwHits.end(); hit++)
337  if (*hit)
338  hitlist->push_back(*hit);
339  else
340  QwError << "null hit in hit list" << QwLog::endl;
341  return hitlist;
342 }
static std::ostream & endl(std::ostream &)
End of the line.
Definition: QwLog.cc:299
std::vector< QwHit * > fQwHits
Array of QwHits.
Definition: QwEvent.h:169
#define QwError
Predefined log drain for errors.
Definition: QwLog.h:40

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

const std::vector<QwHit*>& QwEvent::GetListOfHits ( ) const
inline

Get the list of hits.

Definition at line 232 of file QwEvent.h.

References fQwHits.

232 { return fQwHits; };
std::vector< QwHit * > fQwHits
Array of QwHits.
Definition: QwEvent.h:169
const std::vector< QwPartialTrack * > QwEvent::GetListOfPartialTracks ( EQwRegionID  region,
EQwDetectorPackage  package 
) const

Get the list of partial tracks in region and direction.

Definition at line 488 of file QwEvent.cc.

References fQwPartialTracks.

491 {
492  std::vector<QwPartialTrack*> partialtracks;
493  for (std::vector<QwPartialTrack*>::const_iterator
494  pt = fQwPartialTracks.begin();
495  pt != fQwPartialTracks.end(); pt++) {
496  if ((*pt)->GetRegion() == region) {
497  if ((*pt)->GetPackage() == package) {
498  partialtracks.push_back(*pt);
499  }
500  }
501  }
502  return partialtracks;
503 }
std::vector< QwPartialTrack * > fQwPartialTracks
Array of QwPartialTracks.
Definition: QwEvent.h:177
const std::vector<QwPartialTrack*>& QwEvent::GetListOfPartialTracks ( ) const
inline

Get the list of partial tracks.

Definition at line 282 of file QwEvent.h.

References fQwPartialTracks.

282 { return fQwPartialTracks; };
std::vector< QwPartialTrack * > fQwPartialTracks
Array of QwPartialTracks.
Definition: QwEvent.h:177
const std::vector<QwTrack*>& QwEvent::GetListOfTracks ( ) const
inline

Get the list of tracks.

Definition at line 304 of file QwEvent.h.

References fQwTracks.

304 { return fQwTracks; };
std::vector< QwTrack * > fQwTracks
Array of QwTracks.
Definition: QwEvent.h:181
const std::vector<QwTreeLine*>& QwEvent::GetListOfTreeLines ( ) const
inline

Get the list of tree lines.

Definition at line 256 of file QwEvent.h.

References fQwTreeLines.

256 { return fQwTreeLines; };
std::vector< QwTreeLine * > fQwTreeLines
Array of QwTreeLines.
Definition: QwEvent.h:173
const std::vector< QwTreeLine * > QwEvent::GetListOfTreeLines ( EQwRegionID  region,
EQwDirectionID  direction 
) const

Get the list of tree lines in region and direction.

Definition at line 417 of file QwEvent.cc.

References fQwTreeLines.

418 {
419  std::vector<QwTreeLine*> treelines;
420  for (std::vector<QwTreeLine*>::const_iterator
421  tl = fQwTreeLines.begin();
422  tl != fQwTreeLines.end(); tl++) {
423  if ((*tl)->GetRegion() == region) {
424  if ((*tl)->GetDirection() == direction) {
425  treelines.push_back(*tl);
426  }
427  }
428  }
429  return treelines;
430 }
std::vector< QwTreeLine * > fQwTreeLines
Array of QwTreeLines.
Definition: QwEvent.h:173
Int_t QwEvent::GetNumberOfHits ( ) const
inline

Get the number of hits.

Definition at line 230 of file QwEvent.h.

References fNQwHits.

Referenced by GetHit().

230 { return fNQwHits; };
Int_t fNQwHits
Number of QwHits in the array.
Definition: QwEvent.h:168

+ Here is the caller graph for this function:

Int_t QwEvent::GetNumberOfPartialTracks ( ) const
inline

Get the number of partial tracks.

Definition at line 278 of file QwEvent.h.

References fNQwPartialTracks.

Referenced by GetPartialTrack().

278 { return fNQwPartialTracks; };
Int_t fNQwPartialTracks
Number of QwPartialTracks in the array.
Definition: QwEvent.h:176

+ Here is the caller graph for this function:

Int_t QwEvent::GetNumberOfTracks ( ) const
inline

Get the number of tracks.

Definition at line 302 of file QwEvent.h.

References fNQwTracks.

Referenced by GetTrack().

302 { return fNQwTracks; };
Int_t fNQwTracks
Number of QwTracks in the array.
Definition: QwEvent.h:180

+ Here is the caller graph for this function:

Int_t QwEvent::GetNumberOfTreeLines ( ) const
inline

Get the number of tree lines.

Definition at line 254 of file QwEvent.h.

References fNQwTreeLines.

Referenced by GetTreeLine(), and QwTrackingWorker::ProcessEvent().

254 { return fNQwTreeLines; };
Int_t fNQwTreeLines
Number of QwTreeLines in the array.
Definition: QwEvent.h:172

+ Here is the caller graph for this function:

const QwPartialTrack * QwEvent::GetPartialTrack ( const int  pt) const

Get the specified partial track.

Definition at line 472 of file QwEvent.cc.

References fQwPartialTracks, and GetNumberOfPartialTracks().

473 {
474  if (pt < 0 || pt > GetNumberOfPartialTracks()) return 0;
475  return fQwPartialTracks.at(pt);
476 }
std::vector< QwPartialTrack * > fQwPartialTracks
Array of QwPartialTracks.
Definition: QwEvent.h:177
Int_t GetNumberOfPartialTracks() const
Get the number of partial tracks.
Definition: QwEvent.h:278

+ Here is the call graph for this function:

Double_t QwEvent::GetScatteringAngle ( ) const
inline

Definition at line 319 of file QwEvent.h.

References fScatteringAngle.

319 { return fScatteringAngle;};
double fScatteringAngle
Scattering angle.
Definition: QwEvent.h:336
Double_t QwEvent::GetScatteringVertexZ ( ) const
inline

Definition at line 320 of file QwEvent.h.

References fScatteringVertexZ.

320 { return fScatteringVertexZ;};
double fScatteringVertexZ
Scattering vertex z position.
Definition: QwEvent.h:337
const QwTrack * QwEvent::GetTrack ( const int  t) const

Get the specified track.

Definition at line 545 of file QwEvent.cc.

References fQwTracks, and GetNumberOfTracks().

Referenced by DrawTreeLineVDC(), and QwTrackingWorker::ProcessEvent().

546 {
547  if (t < 0 || t > GetNumberOfTracks()) return 0;
548  return fQwTracks.at(t);
549 }
std::vector< QwTrack * > fQwTracks
Array of QwTracks.
Definition: QwEvent.h:181
Int_t GetNumberOfTracks() const
Get the number of tracks.
Definition: QwEvent.h:302

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

const QwTreeLine * QwEvent::GetTreeLine ( const int  tl) const

Get the specified tree line.

Definition at line 397 of file QwEvent.cc.

References fQwTreeLines, and GetNumberOfTreeLines().

398 {
399  if (tl < 0 || tl > GetNumberOfTreeLines()) return 0;
400  return fQwTreeLines.at(tl);
401 }
Int_t GetNumberOfTreeLines() const
Get the number of tree lines.
Definition: QwEvent.h:254
std::vector< QwTreeLine * > fQwTreeLines
Array of QwTreeLines.
Definition: QwEvent.h:173

+ Here is the call graph for this function:

const TVector3& QwEvent::GetVertexMomentum ( ) const
inline

Definition at line 318 of file QwEvent.h.

References fVertexMomentum.

318 { return fVertexMomentum;};
TVector3 fVertexMomentum
Vertex momentum.
Definition: QwEvent.h:343
const TVector3& QwEvent::GetVertexPosition ( ) const
inline

Definition at line 317 of file QwEvent.h.

References fVertexPosition.

317 { return fVertexPosition;};
TVector3 fVertexPosition
Vertex position.
Definition: QwEvent.h:342
void QwEvent::LoadBeamProperty ( const TString &  map)

Load the beam properties from a map file.

Load the beam properties from a map file.

Parameters
mapName of map file

Definition at line 99 of file QwEvent.cc.

References QwLog::endl(), fBeamEnergy, Qw::GeV, Qw::MeV, QwMessage, and QwParameterFile::TrimComment().

100 {
101  QwParameterFile mapstr(map.Data());
102  while (mapstr.ReadNextLine())
103  {
104  mapstr.TrimComment(); // Remove everything after a comment character.
105  mapstr.TrimWhitespace(); // Get rid of leading and trailing spaces.
106  if (mapstr.LineIsEmpty()) continue;
107 
108  TString varname, varvalue;
109  if (mapstr.HasVariablePair("=",varname,varvalue)) {
110  // This is a declaration line. Decode it.
111  varname.ToLower();
112  if (varname == "energy") {
113  fBeamEnergy = atof(varvalue.Data()) * Qw::MeV;
114  QwMessage << "Beam energy set to " << fBeamEnergy / Qw::GeV << " GeV" << QwLog::endl;
115  }
116  }
117  }
118 }
#define QwMessage
Predefined log drain for regular messages.
Definition: QwLog.h:50
void TrimComment(const char commentchar)
static double fBeamEnergy
Electron beam energy.
Definition: QwEvent.h:320
static const double GeV
Definition: QwUnits.h:95
static std::ostream & endl(std::ostream &)
End of the line.
Definition: QwLog.cc:299
static const double MeV
Definition: QwUnits.h:94

+ Here is the call graph for this function:

void QwEvent::Print ( Option_t *  option = "") const

Print the event.

Definition at line 246 of file QwEvent.cc.

References QwLog::endl(), fEventHeader, fKin, QwKinematics::fP0, QwKinematics::fPp, QwKinematics::fQ2, PrintHits(), PrintPartialTracks(), PrintTracks(), PrintTreeLines(), and QwOut.

247 {
248  // Event header
250  // Event kinematics
251  QwOut << "P0 = " << fKin.fP0 << " GeV/c" << QwLog::endl;
252  QwOut << "PP = " << fKin.fPp << " GeV/c" << QwLog::endl;
253  QwOut << "Q^2 = " << fKin.fQ2 << " (GeV/c)^2" << QwLog::endl;
254  //std::cout << "weight = " << fCrossSectionWeight << std::endl;
255  //std::cout << "energy = " << fTotalEnergy/Qw::MeV << " MeV" << std::endl;
256  //std::cout << "momentum = " << fMomentum / Qw::MeV << " MeV" << std::endl;
257  //std::cout << "vertex position = " << fVertexPosition.Z()/Qw::cm << " cm" << std::endl;
258  //std::cout << "vertex momentum = " << fVertexMomentum.Z()/Qw::MeV << " MeV" << std::endl;
259 
260  // Event content
261  QwOut << "Hits in this event:" << QwLog::endl;
262  PrintHits(option);
263  QwOut << "Tree lines in this event:" << QwLog::endl;
264  PrintTreeLines(option);
265  QwOut << "Partial tracks in this event:" << QwLog::endl;
266  PrintPartialTracks(option);
267  QwOut << "Tracks in this event:" << QwLog::endl;
268  PrintTracks(option);
269 }
#define QwOut
Predefined log drain for explicit output.
Definition: QwLog.h:35
QwKinematics fKin
Inclusive scattering.
Definition: QwEvent.h:348
double fP0
Incoming momentum .
Definition: QwEvent.h:133
double fQ2
Four-momentum transfer squared .
Definition: QwEvent.h:137
double fPp
Outgoing momentum .
Definition: QwEvent.h:136
void PrintTracks(Option_t *option="") const
Print the list of tracks.
Definition: QwEvent.cc:552
QwEventHeader * fEventHeader
Definition: QwEvent.h:161
static std::ostream & endl(std::ostream &)
End of the line.
Definition: QwLog.cc:299
void PrintHits(Option_t *option="") const
Print the list of hits.
Definition: QwEvent.cc:311
void PrintTreeLines(Option_t *option="") const
Print the list of tree lines.
Definition: QwEvent.cc:404
void PrintPartialTracks(Option_t *option="") const
Print the list of partial tracks.
Definition: QwEvent.cc:479

+ Here is the call graph for this function:

void QwEvent::PrintHits ( Option_t *  option = "") const

Print the list of hits.

Definition at line 311 of file QwEvent.cc.

References fQwHits.

Referenced by Print().

312 {
313  for (std::vector<QwHit*>::const_iterator hit = fQwHits.begin();
314  hit != fQwHits.end(); ++hit) {
315  if ((*hit)->GetDriftDistance() != -5)
316  std::cout << **hit << std::endl;
317  }
318 }
std::vector< QwHit * > fQwHits
Array of QwHits.
Definition: QwEvent.h:169

+ Here is the caller graph for this function:

void QwEvent::PrintPartialTracks ( Option_t *  option = "") const

Print the list of partial tracks.

Definition at line 479 of file QwEvent.cc.

References fQwPartialTracks.

Referenced by Print().

480 {
481  for (std::vector<QwPartialTrack*>::const_iterator parttrack =
482  fQwPartialTracks.begin(); parttrack != fQwPartialTracks.end();
483  ++parttrack)
484  std::cout << **parttrack << std::endl;
485 }
std::vector< QwPartialTrack * > fQwPartialTracks
Array of QwPartialTracks.
Definition: QwEvent.h:177

+ Here is the caller graph for this function:

void QwEvent::PrintTracks ( Option_t *  option = "") const

Print the list of tracks.

Definition at line 552 of file QwEvent.cc.

References fQwTracks.

Referenced by Print().

553 {
554  for (std::vector<QwTrack*>::const_iterator track = fQwTracks.begin();
555  track != fQwTracks.end(); track++)
556  std::cout << **track << std::endl;
557 }
std::vector< QwTrack * > fQwTracks
Array of QwTracks.
Definition: QwEvent.h:181

+ Here is the caller graph for this function:

void QwEvent::PrintTreeLines ( Option_t *  option = "") const

Print the list of tree lines.

Definition at line 404 of file QwEvent.cc.

References fQwTreeLines, and QwTreeLine::next.

Referenced by Print().

405 {
406  for (std::vector<QwTreeLine*>::const_iterator treeline = fQwTreeLines.begin();
407  treeline != fQwTreeLines.end(); treeline++) {
408  std::cout << **treeline << std::endl;
409  QwTreeLine* tl = (*treeline)->next;
410  while (tl) {
411  std::cout << *tl << std::endl;
412  tl = tl->next;
413  }
414  }
415 }
std::vector< QwTreeLine * > fQwTreeLines
Array of QwTreeLines.
Definition: QwEvent.h:173
One-dimensional (u, v, or x) track stubs and associated hits.
Definition: QwTreeLine.h:51
QwTreeLine * next
Definition: QwTreeLine.h:231

+ Here is the caller graph for this function:

void QwEvent::Reset ( Option_t *  option = "")

Definition at line 86 of file QwEvent.cc.

References ResetHits(), ResetPartialTracks(), ResetTracks(), and ResetTreeLines().

Referenced by main().

87 {
88  ResetHits(option);
89  ResetTreeLines(option);
90  ResetPartialTracks(option);
91  ResetTracks(option);
92 }
void ResetTracks(Option_t *option="")
Reset the list of tracks.
Definition: QwEvent.cc:539
void ResetPartialTracks(Option_t *option="")
Reset the list of partial tracks.
Definition: QwEvent.cc:466
void ResetHits(Option_t *option="")
Reset the list of hits.
Definition: QwEvent.cc:298
void ResetTreeLines(Option_t *option="")
Reset the list of tree lines.
Definition: QwEvent.cc:391

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void QwEvent::ResetHits ( Option_t *  option = "")

Reset the list of hits.

Definition at line 298 of file QwEvent.cc.

References ClearHits().

Referenced by Reset().

299 {
300  ClearHits();
301 }
void ClearHits(Option_t *option="")
Clear the list of hits.
Definition: QwEvent.cc:289

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void QwEvent::ResetPartialTracks ( Option_t *  option = "")

Reset the list of partial tracks.

Definition at line 466 of file QwEvent.cc.

References ClearPartialTracks().

Referenced by Reset().

467 {
469 }
void ClearPartialTracks(Option_t *option="")
Clear the list of partial tracks.
Definition: QwEvent.cc:457

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void QwEvent::ResetTracks ( Option_t *  option = "")

Reset the list of tracks.

Definition at line 539 of file QwEvent.cc.

References ClearTracks().

Referenced by Reset().

540 {
541  ClearTracks();
542 }
void ClearTracks(Option_t *option="")
Clear the list of tracks.
Definition: QwEvent.cc:530

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void QwEvent::ResetTreeLines ( Option_t *  option = "")

Reset the list of tree lines.

Definition at line 391 of file QwEvent.cc.

References ClearTreeLines().

Referenced by Reset().

392 {
393  ClearTreeLines();
394 }
void ClearTreeLines(Option_t *option="")
Clear the list of tree lines.
Definition: QwEvent.cc:380

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void QwEvent::SetEventHeader ( const QwEventHeader eventheader)
inline

Definition at line 193 of file QwEvent.h.

References fEventHeader.

Referenced by QwTreeEventBuffer::GetSpecificEvent().

193  {
194  if (fEventHeader) delete fEventHeader;
195  fEventHeader = new QwEventHeader(eventheader);
196  };
Contains header information of a tracked event.
Definition: QwEvent.h:37
QwEventHeader * fEventHeader
Definition: QwEvent.h:161

+ Here is the caller graph for this function:

void QwEvent::SetEventHeader ( const QwEventHeader eventheader)
inline

Definition at line 197 of file QwEvent.h.

References fEventHeader.

197  {
198  if (fEventHeader) delete fEventHeader;
199  fEventHeader = new QwEventHeader(*eventheader);
200  };
Contains header information of a tracked event.
Definition: QwEvent.h:37
QwEventHeader * fEventHeader
Definition: QwEvent.h:161

Field Documentation

Double_t QwEvent::fBeamEnergy = 1.165 * Qw::GeV
static

Electron beam energy.

Definition at line 320 of file QwEvent.h.

Referenced by CalculateKinematics(), and LoadBeamProperty().

double QwEvent::fCrossSection

Definition at line 339 of file QwEvent.h.

Referenced by QwTreeEventBuffer::GetSpecificEvent().

QwEventHeader* QwEvent::fEventHeader

Definition at line 161 of file QwEvent.h.

Referenced by GetEventHeader(), Print(), QwEvent(), SetEventHeader(), and ~QwEvent().

double QwEvent::fHydrogenEnergyLoss

Pre-scattering target energy loss assuming LH2 target.

Definition at line 335 of file QwEvent.h.

Referenced by CalculateKinematics().

QwKinematics QwEvent::fKin

Inclusive scattering.

Definition at line 348 of file QwEvent.h.

Referenced by CalculateKinematics(), QwTreeEventBuffer::GetSpecificEvent(), and Print().

QwKinematics QwEvent::fKinElastic

Scattering assuming elastic reaction.

Definition at line 350 of file QwEvent.h.

Referenced by CalculateKinematics().

QwKinematics QwEvent::fKinElasticWithLoss

Scattering assuming elastic reaction and hydrogen energy loss.

Definition at line 351 of file QwEvent.h.

Referenced by CalculateKinematics().

QwKinematics QwEvent::fKinWithLoss

Scattering with hydrogen energy loss.

Definition at line 349 of file QwEvent.h.

Referenced by CalculateKinematics().

std::vector<float> QwEvent::fMD_LeftNbOfPEs

Definition at line 328 of file QwEvent.h.

Referenced by QwTreeEventBuffer::GetSpecificEvent().

std::vector<float> QwEvent::fMD_RightNbOfPEs

Definition at line 329 of file QwEvent.h.

Referenced by QwTreeEventBuffer::GetSpecificEvent().

std::vector<float> QwEvent::fMD_TotalNbOfPEs

Definition at line 330 of file QwEvent.h.

Referenced by QwTreeEventBuffer::GetSpecificEvent().

Int_t QwEvent::fNQwHits

Number of QwHits in the array.

Definition at line 168 of file QwEvent.h.

Referenced by AddHit(), ClearHits(), GetNumberOfHits(), and QwEvent().

Int_t QwEvent::fNQwPartialTracks

Number of QwPartialTracks in the array.

Definition at line 176 of file QwEvent.h.

Referenced by AddPartialTrack(), ClearPartialTracks(), GetNumberOfPartialTracks(), and QwEvent().

Int_t QwEvent::fNQwTracks

Number of QwTracks in the array.

Definition at line 180 of file QwEvent.h.

Referenced by AddTrack(), ClearTracks(), GetNumberOfTracks(), and QwEvent().

Int_t QwEvent::fNQwTreeLines

Number of QwTreeLines in the array.

Definition at line 172 of file QwEvent.h.

Referenced by AddTreeLine(), ClearTreeLines(), GetNumberOfTreeLines(), and QwEvent().

double QwEvent::fOriginVertexEnergy

Definition at line 341 of file QwEvent.h.

Referenced by QwTreeEventBuffer::GetSpecificEvent().

double QwEvent::fPreScatteringEnergy

Definition at line 340 of file QwEvent.h.

Referenced by QwTreeEventBuffer::GetSpecificEvent().

Double_t QwEvent::fPrimaryQ2

< Momentum transfer Q^2 assuming elastic scattering with hydrogen energy loss

Definition at line 355 of file QwEvent.h.

Referenced by CalculateKinematics(), and QwTreeEventBuffer::GetSpecificEvent().

std::vector<QwHit*> QwEvent::fQwHits

Array of QwHits.

Definition at line 169 of file QwEvent.h.

Referenced by AddHit(), ClearHits(), GetHit(), GetHitContainer(), GetListOfHits(), and PrintHits().

std::vector<QwPartialTrack*> QwEvent::fQwPartialTracks

Array of QwPartialTracks.

Definition at line 177 of file QwEvent.h.

Referenced by AddPartialTrack(), ClearPartialTracks(), GetListOfPartialTracks(), GetPartialTrack(), and PrintPartialTracks().

std::vector<QwTrack*> QwEvent::fQwTracks

Array of QwTracks.

Definition at line 181 of file QwEvent.h.

Referenced by AddTrack(), ClearTracks(), GetListOfTracks(), GetTrack(), and PrintTracks().

std::vector<QwTreeLine*> QwEvent::fQwTreeLines

Array of QwTreeLines.

Definition at line 173 of file QwEvent.h.

Referenced by AddTreeLine(), ClearTreeLines(), GetListOfTreeLines(), GetTreeLine(), and PrintTreeLines().

double QwEvent::fScatteringAngle

Scattering angle.

Definition at line 336 of file QwEvent.h.

Referenced by CalculateKinematics(), GetScatteringAngle(), and QwTreeEventBuffer::GetSpecificEvent().

double QwEvent::fScatteringVertexR

Scattering vertex radial distance.

Definition at line 338 of file QwEvent.h.

Referenced by CalculateKinematics(), and QwTreeEventBuffer::GetSpecificEvent().

double QwEvent::fScatteringVertexZ

Scattering vertex z position.

Definition at line 337 of file QwEvent.h.

Referenced by CalculateKinematics(), GetScatteringVertexZ(), and QwTreeEventBuffer::GetSpecificEvent().

list of tree lines [upper/lower][region][u/v/x/y]

Definition at line 358 of file QwEvent.h.

Referenced by QwEvent(), and ~QwEvent().

TVector3 QwEvent::fVertexMomentum

Vertex momentum.

Definition at line 343 of file QwEvent.h.

Referenced by CalculateKinematics(), QwTreeEventBuffer::GetSpecificEvent(), and GetVertexMomentum().

TVector3 QwEvent::fVertexPosition

Vertex position.

Definition at line 342 of file QwEvent.h.

Referenced by CalculateKinematics(), QwTreeEventBuffer::GetSpecificEvent(), and GetVertexPosition().


The documentation for this class was generated from the following files: