Hall A ROOT/C++ Analyzer (podd)
Loading...
Searching...
No Matches
THaSpectrometer.h
Go to the documentation of this file.
1#ifndef Podd_THaSpectrometer_h_
2#define Podd_THaSpectrometer_h_
3
5//
6// THaSpectrometer
7//
9
10#include "THaApparatus.h"
11#include "THaVertexModule.h"
12#include "THaTrackingModule.h"
13#include "TClonesArray.h"
14#include "TVector3.h"
15#include "TRotation.h"
16#include "THaParticleInfo.h"
17#include "THaPidDetector.h"
18#include <cassert>
19
20class THaTrack;
21class TList;
22class THaCut;
23
25 public THaVertexModule {
26
27public:
28 virtual ~THaSpectrometer();
29
30 // Main functions
31 virtual void Clear( Option_t* opt="");
32 virtual Int_t CoarseTrack();
33 virtual Int_t CoarseReconstruct();
34 virtual EStatus Init( const TDatime& run_time );
35 virtual Int_t Track();
36 virtual Int_t Reconstruct();
37 virtual Int_t CalcPID();
38 virtual Int_t FindVertices( TClonesArray& tracks ) = 0;
39 virtual Int_t TrackCalc() = 0;
40
41 // Auxiliary functions
42 virtual Int_t AddDetector( THaDetector* det, Bool_t quiet = false,
43 Bool_t first = false );
44 virtual Int_t AddPidParticle( const char* shortname,
45 const char* name,
46 Double_t mass, Int_t charge = 0 );
47 virtual void DefinePidParticles();
49 Int_t GetNpidParticles() const;
50 Int_t GetNpidDetectors() const;
53 Int_t GetNTracks() const { return fTracks->GetLast()+1; }
54 TClonesArray* GetTracks() const { return fTracks; }
55 TClonesArray* GetTrackPID() const { return fTrackPID; }
56 virtual const TVector3& GetVertex() const;
57 virtual Bool_t HasVertex() const;
58
59 Bool_t IsDone( UInt_t stage ) const;
60 Bool_t IsPID() const { return fPID; }
62 void SetPID( Bool_t b = true ) { fPID = b; }
63
64 // The following is specific to small-acceptance pointing spectrometers
65 // using spectrometer-specific coordinates such as TRANSPORT
66 const TRotation& GetToLabRot() const { return fToLabRot; }
67 const TRotation& GetToTraRot() const { return fToTraRot; }
68 const TVector3& GetPointingOffset() const { return fPointingOffset; }
69 Double_t GetThetaGeo() const { return fThetaGeo; }
70 Double_t GetPhiGeo() const { return fPhiGeo; }
71 Double_t GetThetaSph() const { return fThetaSph; }
72 Double_t GetPhiSph() const { return fPhiSph; }
73 Double_t GetPcentral() const { return fPcentral; }
74 Double_t GetCollDist() const { return fCollDist; }
75
77 Bool_t bend_down );
78
79 virtual void TrackToLab( THaTrack& track, TVector3& pvect ) const;
80 virtual void TransportToLab( Double_t p, Double_t th,
81 Double_t ph, TVector3& pvect ) const;
82 virtual void LabToTransport( const TVector3& vertex,
83 const TVector3& pvect,
84 TVector3& tvertex,
85 Double_t* ray ) const;
86 void LabToTransport( const TVector3& vertex,
87 const TVector3& pvect,
88 Double_t* ray ) const;
95
96protected:
97
98 static const Int_t kInitTrackMultiplicity = 5;
99
101 TClonesArray* fTrackPID; //PID info for the tracks
102 TList* fTrackingDetectors; //Tracking detectors
103 TList* fNonTrackingDetectors; //Non-tracking detectors
104 TObjArray* fPidDetectors; //PID detectors
105 TObjArray* fPidParticles; //Particles for which we want PID
106 THaTrack* fGoldenTrack; //Golden track within fTracks
107
108 // The following is specific to small-acceptance pointing spectrometers
109 TRotation fToLabRot; //Rotation matrix from TRANSPORT to lab
110 TRotation fToTraRot; //Rotation matrix from lab to TRANSPORT
111 TVector3 fPointingOffset; //Optical point in lab coordinate system
112 Double_t fThetaGeo; //In-plane geographic central angle (rad)
113 Double_t fPhiGeo; //Out-of-plane geographic central angle (rad)
114 Double_t fThetaSph, fPhiSph; //Central angles in spherical coords. (rad)
115 Double_t fSinThGeo, fCosThGeo; //Sine and cosine of central angles
116 Double_t fSinPhGeo, fCosPhGeo; // in geographical coordinates
117 Double_t fSinThSph, fCosThSph; //Sine and cosine of central angles in
118 Double_t fSinPhSph, fCosPhSph; // spherical coordinates
119 Double_t fPcentral; //Central momentum (GeV)
120 Double_t fCollDist; //Distance from collimator to target center (m)
121
122 // Status flags
123 UInt_t fStagesDone; //Bitfield of completed analysis stages
124 Bool_t fPID; //PID enabled
125
126 // only derived classes can construct me
127 THaSpectrometer( const char* name, const char* description );
128
129 virtual Int_t DefineVariables( EMode mode = kDefine );
130 virtual Int_t ReadRunDatabase( const TDatime& date );
131 virtual void ListInit(); // Initialize lists of detector types
132 virtual void PidInit(); // Initialize PID structures
133
134 ClassDef(THaSpectrometer,1) // A generic spectrometer
135};
136
137
138//---------------- inlines ----------------------------------------------------
140{
141 return fPidParticles->GetLast()+1;
142}
143
144//_____________________________________________________________________________
146{
147 return fPidDetectors->GetLast()+1;
148}
149
150//_____________________________________________________________________________
152{
153 assert(dynamic_cast<THaParticleInfo*>(fPidParticles->At(i)));
154 return static_cast<THaParticleInfo*>( fPidParticles->At(i) );
155}
156
157//_____________________________________________________________________________
159{
160 assert(dynamic_cast<THaPidDetector*>(fPidDetectors->At(i)));
161 return static_cast<THaPidDetector*>( fPidDetectors->At(i) );
162}
163
164//_____________________________________________________________________________
165inline
167{
168 return ((fStagesDone & stage) != 0);
169}
170
171//_____________________________________________________________________________
172inline
174 const TVector3& pvect,
175 Double_t* ray ) const
176{
177 TVector3 dummy;
178 LabToTransport( vertex, pvect, dummy, ray );
179}
180
181#endif
182
int Int_t
unsigned int UInt_t
bool Bool_t
double Double_t
const char Option_t
#define ClassDef(name, id)
#define BIT(n)
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t b
Bool_t IsDone(UInt_t stage) const
virtual Int_t AddDetector(THaDetector *det, Bool_t quiet=false, Bool_t first=false)
Double_t GetPcentral() const
THaTrack * fGoldenTrack
virtual Bool_t HasVertex() const
virtual Int_t CoarseReconstruct()
THaParticleInfo * GetPidParticleInfo(Int_t i) const
void SetPID(Bool_t b=true)
THaTrack * GetGoldenTrack() const
Int_t GetNTracks() const
TClonesArray * GetTracks() const
virtual void DefinePidParticles()
Int_t GetNpidDetectors() const
void SetCentralAngles(Double_t th, Double_t ph, Bool_t bend_down)
Int_t GetNpidParticles() const
const TRotation & GetToTraRot() const
virtual void LabToTransport(const TVector3 &vertex, const TVector3 &pvect, TVector3 &tvertex, Double_t *ray) const
Double_t GetPhiGeo() const
Double_t GetCollDist() const
TClonesArray * fTracks
const TRotation & GetToLabRot() const
virtual Int_t TrackCalc()=0
virtual Int_t AddPidParticle(const char *shortname, const char *name, Double_t mass, Int_t charge=0)
virtual Int_t ReadRunDatabase(const TDatime &date)
const TVector3 & GetPointingOffset() const
TClonesArray * GetTrackPID() const
THaPidDetector * GetPidDetector(Int_t i) const
virtual void TransportToLab(Double_t p, Double_t th, Double_t ph, TVector3 &pvect) const
Double_t GetThetaSph() const
virtual void Clear(Option_t *opt="")
virtual Int_t CoarseTrack()
static const Int_t kInitTrackMultiplicity
TObjArray * fPidDetectors
virtual Int_t Reconstruct()
virtual Int_t FindVertices(TClonesArray &tracks)=0
virtual void ListInit()
TClonesArray * fTrackPID
virtual void PidInit()
virtual Int_t CalcPID()
void SetGoldenTrack(THaTrack *t)
virtual const TVector3 & GetVertex() const
Double_t GetThetaGeo() const
virtual Int_t DefineVariables(EMode mode=kDefine)
TList * fNonTrackingDetectors
Double_t GetPhiSph() const
virtual void TrackToLab(THaTrack &track, TVector3 &pvect) const
virtual Int_t Track()
TObjArray * fPidParticles
virtual ~THaSpectrometer()
Bool_t IsPID() const
TObject * At(Int_t idx) const override
Int_t GetLast() const override
double * vertex