QwAnalysis
VQwBridgingMethod.h
Go to the documentation of this file.
1 /*!
2  * \file VQwBridgingMethod.h
3  * \brief Definition of the bridging method interface
4  *
5  * \author Wouter Deconinck
6  * \date 2010-01-23
7  */
8 
9 #ifndef VQWBRIDGINGMETHOD_H
10 #define VQWBRIDGINGMETHOD_H
11 
12 // System headers
13 #include <iostream>
14 
15 // ROOT headers
16 #include "TVector3.h"
17 
18 // Qweak headers
19 #include "QwUnits.h"
20 #include "QwTrack.h"
21 #include "QwParameterFile.h"
22 
23 // Forward declarations
24 class QwPartialTrack;
25 
26 /**
27  * \class VQwBridgingMethod
28  * \ingroup QwTracking
29  * \brief Interface to the various bridging methods
30  */
32 
33  public:
34 
35  /// Default constructor
37  /// Destructor
38  virtual ~VQwBridgingMethod() { };
39 
40 
41  /// \brief Bridge from the front to back partial track (pure virtual)
42  virtual const QwTrack* Bridge(const QwPartialTrack* front, const QwPartialTrack* back) = 0;
43 
44 
45  protected:
46 
47  /// Estimate the momentum based only on the direction
48  virtual double EstimateInitialMomentum(const TVector3& direction) const;
49 
50 }; // class VQwBridgingMethod
51 
52 
53 /**
54  * Estimate the momentum for a first attempt at bridging.
55  *
56  * This takes into account only the direction of the front partial track and assumes
57  * an elastic electron-proton reaction without any energy loss. Any effects of the
58  * energy loss should be small enough that they can be recovered from the bridging
59  * algorithm in little time, and this avoids any bias of LH2 over Al events (E_loss).
60  *
61  * @param direction Direction of front partial track
62  * @return Initial momentum
63  */
64 inline double VQwBridgingMethod::EstimateInitialMomentum(const TVector3& direction) const
65 {
66  double Mp = Qw::Mp; // proton mass
67  double e0 = 1.165 * Qw::GeV; // beam energy
68  // \todo Get the beam energy from the event (?)
69 
70  // Kinematics for elastic e+p scattering
71  double cth = direction.CosTheta();
72  return e0 / (1.0 + e0 / Mp * (1 - cth));
73 }
74 
75 #endif // VQWBRIDGINGMETHOD_H
virtual double EstimateInitialMomentum(const TVector3 &direction) const
Estimate the momentum based only on the direction.
Definition of the track class.
Interface to the various bridging methods.
Contains the complete track as a concatenation of partial tracks.
Definition: QwTrack.h:30
static const double GeV
Definition: QwUnits.h:95
static const double Mp
Mass of the proton.
Definition: QwUnits.h:119
virtual const QwTrack * Bridge(const QwPartialTrack *front, const QwPartialTrack *back)=0
Bridge from the front to back partial track (pure virtual)
VQwBridgingMethod()
Default constructor.
virtual ~VQwBridgingMethod()
Destructor.
Contains the straight part of a track in one region only.