Hall A ROOT/C++ Analyzer (podd)
Loading...
Searching...
No Matches
THaVDCPlane.h
Go to the documentation of this file.
1#ifndef Podd_THaVDCPlane_h_
2#define Podd_THaVDCPlane_h_
3
5// //
6// THaVDCPlane //
7// //
9
10#include "THaSubDetector.h"
11#include "THaVDCWire.h"
12#include "THaVDCCluster.h"
13#include "TClonesArray.h"
14#include "THaVDCHit.h"
15#include <cassert>
16#include <vector>
17
18namespace VDC {
19 class TimeToDistConv;
20}
21
22class THaEvData;
23class THaVDC;
24
26
27public:
28
29 explicit THaVDCPlane( const char* name="", const char* description="",
30 THaDetectorBase* parent = nullptr );
31 virtual ~THaVDCPlane();
32
33 virtual void Clear( Option_t* opt="" );
34 virtual Int_t Decode( const THaEvData& ); // Raw data -> hits
35 virtual Int_t ApplyTimeCorrection(); // Drift time correction
36 virtual Int_t FindClusters(); // Hits -> clusters
37 virtual Int_t FitTracks(); // Clusters -> tracks
38
39 //Get and Set functions
40 Int_t GetNClusters() const { return fClusters->GetLast()+1; }
41 TClonesArray* GetClusters() const { return fClusters; }
43 { assert( i>=0 && i<GetNClusters() );
44 return static_cast<THaVDCCluster*>( fClusters->UncheckedAt(i) ); }
45
46 Int_t GetNWires() const { return fWires->GetLast()+1; }
47 TClonesArray* GetWires() const { return fWires; }
49 { assert( i>=0 && i<GetNWires() );
50 return static_cast<THaVDCWire*>( fWires->UncheckedAt(i) ); }
51
52 Int_t GetNHits() const { return fHits->GetLast()+1; }
53 TClonesArray* GetHits() const { return fHits; }
55 { assert( i>=0 && i<GetNHits() );
56 return static_cast<THaVDCHit*>( fHits->UncheckedAt(i) ); }
57
58 Int_t GetNWiresHit() const { return fNWiresHit; }
59
60 const TVector3& GetCenter() const { return fCenter; }
61 Double_t GetZ() const { return fCenter.Z(); }
62 Double_t GetWBeg() const { return fWBeg; }
63 Double_t GetWSpac() const { return fWSpac; }
64 Double_t GetWAngle() const { return fWAngle; }
65 Double_t GetCosWAngle() const { return fCosWAngle; }
66 Double_t GetSinWAngle() const { return fSinWAngle; }
67 Double_t GetTDCRes() const { return fTDCRes; }
68 Double_t GetDriftVel() const { return fDriftVel; }
69 Double_t GetMinTime() const { return fMinTime; }
70 Double_t GetMaxTime() const { return fMaxTime; }
71 Double_t GetMaxTdiff() const { return fMaxTdiff; }
73
74// Double_t GetT0() const { return fT0; }
75// Int_t GetNumBins() const { return fNumBins; }
76// Float_t *GetTable() const { return fTable; }
77
78 // Override function from THaDetectorBase to check positions
79 // relative to fCenter
80 virtual Bool_t IsInActiveArea( Double_t x, Double_t y ) const;
81
82 void UpdateGeometry( Double_t x, Double_t y, bool force = false );
83
84protected:
85
86 // Event data
87 //Use TClonesArray::GetLast()+1 to get the number of wires, hits, & clusters
89 TClonesArray* fHits; // Fired wires
90 TClonesArray* fClusters; // Clusters
91
92 Int_t fNHits; // Total number of hits (including multihits)
93 Int_t fNWiresHit; // Number of wires with one or more hits
94 Int_t fNpass; // Number of passes over hits in FindClusters()
95
96 // Configuration
97 Int_t fMinClustSize; // Minimum number of wires needed for a cluster
98 Int_t fMaxClustSpan; // Maximum size of cluster in wire spacings
99 Int_t fNMaxGap; // Max gap in wire numbers in a cluster
100 Int_t fMinTime; // Min and Max limits of TDC times for clusters
102 UInt_t fMaxThits; // Max allowed number of hits per wire per event
103 Double_t fMinTdiff; // Min and Max limits of times between wires in cluster
105 Double_t fTDCRes; // TDC Resolution ( s / channel)
106 Double_t fDriftVel; // Drift velocity in the wire plane (m/s)
107 Double_t fT0Resolution; // (Average) resolution of cluster time offset fit
108 Bool_t fOnlyFastestHit; // Only record earliest hit for each wire
109 Bool_t fNoNegativeTime; // Disallow negative drift times
110
111 // Geometry
112 TVector3 fCenter; // Plane center in VDC coordinate system (m)
113 Double_t fWBeg; // Position of 1st wire along wire coord sys (m)
114 Double_t fWSpac; // Wire spacing and direction (m)
115 Double_t fWAngle; // Angle (rad) between dispersive direction (x) and
116 // normal to wires in dir of increasing wire position
117 Double_t fSinWAngle; // sine of fWAngle, for efficiency
118 Double_t fCosWAngle; // cosine of fWAngle, for efficiency
119
120 // Lookup table parameters
121// Double_t fT0; // calculated zero time
122// Int_t fNumBins; // size of lookup table
123// Float_t *fTable; // time-to-distance lookup table
124
125 VDC::TimeToDistConv* fTTDConv; // Time-to-distance converter for this plane's wires
126
127 THaVDC* fVDC; // VDC detector to which this plane belongs
128
129 // Temporary storage shared between member functions
133
134 virtual void MakePrefix();
135 virtual Int_t ReadDatabase( const TDatime& date );
136 virtual Int_t DefineVariables( EMode mode = kDefine );
137 virtual Int_t ReadGeometry( FILE* file, const TDatime& date,
138 Bool_t required = false );
139
140 virtual Int_t StoreHit( const DigitizerHitInfo_t& hitinfo, UInt_t data );
141 virtual void PrintDecodedData( const THaEvData& evdata ) const;
142
143private:
144 Int_t ReadDatabaseErrcheck( const std::vector<Float_t>& tdc_offsets,
145 const char* here );
146 Int_t ReadGeometryErrcheck( const std::vector<Double_t>& position,
147 const std::vector<Double_t>& size,
148 const char* here );
149 Int_t CreateTTDConv( const char* classname,
150 const std::vector<Double_t>& ttd_param,
151 const char* here );
152
153ClassDef(THaVDCPlane,0) // VDCPlane class
154};
155
157
158#endif
int Int_t
unsigned int UInt_t
bool Bool_t
double Double_t
const char Option_t
#define ClassDef(name, id)
char name[80]
static const char *const here
Definition THaVar.cxx:64
Int_t CreateTTDConv(const char *classname, const std::vector< Double_t > &ttd_param, const char *here)
Double_t fTDCRes
Double_t GetMaxTdiff() const
Definition THaVDCPlane.h:71
Double_t GetT0Resolution() const
Definition THaVDCPlane.h:72
Double_t fWAngle
virtual Int_t FitTracks()
Double_t GetMaxTime() const
Definition THaVDCPlane.h:70
Int_t GetNWiresHit() const
Definition THaVDCPlane.h:58
Int_t ReadDatabaseErrcheck(const std::vector< Float_t > &tdc_offsets, const char *here)
THaVDCWire * GetWire(Int_t i) const
Definition THaVDCPlane.h:48
Bool_t fNoNegativeTime
THaVDCCluster * GetCluster(Int_t i) const
Definition THaVDCPlane.h:42
Int_t fMaxClustSpan
Definition THaVDCPlane.h:98
Int_t GetNWires() const
Definition THaVDCPlane.h:46
Double_t GetSinWAngle() const
Definition THaVDCPlane.h:66
Double_t GetTDCRes() const
Definition THaVDCPlane.h:67
THaVDCHit * GetHit(Int_t i) const
Definition THaVDCPlane.h:54
TClonesArray * fWires
Definition THaVDCPlane.h:88
const TVector3 & GetCenter() const
Definition THaVDCPlane.h:60
Double_t fSinWAngle
TClonesArray * GetHits() const
Definition THaVDCPlane.h:53
Int_t fNWiresHit
Definition THaVDCPlane.h:93
virtual void Clear(Option_t *opt="")
Int_t fNMaxGap
Definition THaVDCPlane.h:99
Double_t fWSpac
Double_t GetWAngle() const
Definition THaVDCPlane.h:64
THaVDC * fVDC
virtual Int_t ReadDatabase(const TDatime &date)
TClonesArray * fHits
Definition THaVDCPlane.h:89
Double_t GetWSpac() const
Definition THaVDCPlane.h:63
TClonesArray * fClusters
Definition THaVDCPlane.h:90
Int_t fMinClustSize
Definition THaVDCPlane.h:97
TClonesArray * GetClusters() const
Definition THaVDCPlane.h:41
Double_t fMaxTdiff
THaVDCWire * fPrevWire
Double_t fWBeg
virtual Int_t FindClusters()
TClonesArray * GetWires() const
Definition THaVDCPlane.h:47
Int_t GetNClusters() const
Definition THaVDCPlane.h:40
UInt_t fMaxData
virtual Int_t Decode(const THaEvData &)
virtual Bool_t IsInActiveArea(Double_t x, Double_t y) const
virtual Int_t ApplyTimeCorrection()
void UpdateGeometry(Double_t x, Double_t y, bool force=false)
virtual Int_t StoreHit(const DigitizerHitInfo_t &hitinfo, UInt_t data)
Double_t GetMinTime() const
Definition THaVDCPlane.h:69
virtual ~THaVDCPlane()
Double_t fMinTdiff
Double_t GetZ() const
Definition THaVDCPlane.h:61
Double_t fCosWAngle
Double_t fT0Resolution
virtual void MakePrefix()
VDC::TimeToDistConv * fTTDConv
Bool_t fOnlyFastestHit
Double_t GetWBeg() const
Definition THaVDCPlane.h:62
virtual Int_t DefineVariables(EMode mode=kDefine)
TVector3 fCenter
Double_t fDriftVel
UInt_t fMaxThits
Int_t ReadGeometryErrcheck(const std::vector< Double_t > &position, const std::vector< Double_t > &size, const char *here)
virtual void PrintDecodedData(const THaEvData &evdata) const
Double_t GetCosWAngle() const
Definition THaVDCPlane.h:65
Int_t GetNHits() const
Definition THaVDCPlane.h:52
virtual Int_t ReadGeometry(FILE *file, const TDatime &date, Bool_t required=false)
Double_t GetDriftVel() const
Definition THaVDCPlane.h:68
TObject * UncheckedAt(Int_t i) const
Int_t GetLast() const override
Double_t Z() const