Hall A ROOT/C++ Analyzer (podd)
Loading...
Searching...
No Matches
THaVDCChamber.cxx
Go to the documentation of this file.
1
2// //
3// THaVDCChamber //
4// //
5// Representation of a VDC chamber, consisting of one U plane and one //
6// V plane in close proximity. With the Hall A VDCs, the two planes of a //
7// chamber share the same mechanical frame and gas system. //
8// //
9// This particular implementation takes the U plane as the reference for //
10// specifying coordinates. It was also written assuming that the U plane is //
11// below the V plane (smaller z). It should work the other way around as //
12// well but that configuration has not been tested. If in doubt, one can //
13// always swap the meaning of "U" and "V" such that "U" always means the //
14// first (bottom) plane. //
15// //
17
18#include "THaVDCChamber.h"
19#include "THaVDCPlane.h"
20#include "THaVDCPoint.h"
21#include "THaVDCCluster.h"
22#include "TMath.h"
23
24//_____________________________________________________________________________
25THaVDCChamber::THaVDCChamber( const char* name, const char* description,
26 THaDetectorBase* parent ) :
27 THaSubDetector(name,description,parent),
28 // Create the U and V planes
29 fU{new THaVDCPlane( "u", "U plane", this )},
30 fV{new THaVDCPlane( "v", "V plane", this )},
31 // Create array for cluster pairs (points) representing hits
32 fPoints{new TClonesArray("THaVDCPoint", 10)}, // 10 is arbitrary
33 fSpacing(0), fSin_u(0), fCos_u(1), fSin_v(1), fCos_v(0), fInv_sin_vu(0)
34{
35 // Constructor
36
37}
38
39//_____________________________________________________________________________
41{
42 // Destructor. Delete plane objects and point array.
43
45 delete fU;
46 delete fV;
47 delete fPoints;
48}
49
50//_____________________________________________________________________________
52{
53 // Initialize the chamber class. Calls its own Init(), then initializes
54 // subdetectors, then calculates local geometry data.
55
56 if( IsZombie() || !fV || !fU )
57 return fStatus = kInitError;
58
59 if( (fStatus = THaSubDetector::Init(date )) ||
60 (fStatus = fU->Init(date )) ||
61 (fStatus = fV->Init(date )))
62 return fStatus;
63
64 fSpacing = fV->GetZ() - fU->GetZ(); // Space between U & V wire planes
65
66 // Precompute and store values for efficiency
72
73 // Use the same coordinate system as the planes
74 fXax = fU->GetXax();
75 fYax = fU->GetYax();
76 fZax = fU->GetZax();
77
78 // Calculate plane center x/y coordinates, assuming that the plane's
79 // wires are arranged symmetrically around the center
80 Double_t ubeg = fU->GetWBeg();
81 Double_t vbeg = fV->GetWBeg();
82 Double_t uend = fU->GetWBeg() + (fU->GetNWires()-1)*fU->GetWSpac();
83 Double_t vend = fV->GetWBeg() + (fV->GetNWires()-1)*fV->GetWSpac();
84 Double_t xc = 0.5 * ( UVtoX(ubeg,vbeg) + UVtoX(uend,vend) );
85 Double_t yc = 0.5 * ( UVtoY(ubeg,vbeg) + UVtoY(uend,vend) );
86 fU->UpdateGeometry(xc,yc);
87 fV->UpdateGeometry(xc,yc);
88
89 // Take the u plane as our reference plane
90 fOrigin = fU->GetOrigin();
91
92 // Estimate our size
93 fSize[0] = 0.5*TMath::Max( fU->GetXSize(), fU->GetXSize() );
94 fSize[1] = 0.5*TMath::Max( fU->GetYSize(), fU->GetYSize() );
95 fSize[2] = fSpacing + 0.5*fU->GetZSize() + 0.5*fV->GetZSize();
96
97 return fStatus = kOK;
98}
99
100//_____________________________________________________________________________
102{
103 // initialize global variables
104
105 // Register variables in global list
107 if( ret )
108 return ret;
109
110 RVarDef vars[] = {
111 { "npt", "Number of space points", "GetNPoints()" },
112 { "pt.x", "Point center x (m)", "fPoints.THaVDCPoint.X()" },
113 { "pt.y", "Point center y (m)", "fPoints.THaVDCPoint.Y()" },
114 { "pt.z", "Point center z (m)", "fPoints.THaVDCPoint.Z()" },
115 { "pt.d_x", "Point VDC x (m)", "fPoints.THaVDCPoint.fX" },
116 { "pt.d_y", "Point VDC y (m)", "fPoints.THaVDCPoint.fY" },
117 { "pt.d_th", "Point VDC tan(theta)", "fPoints.THaVDCPoint.fTheta" },
118 { "pt.d_ph", "Point VDC tan(phi)", "fPoints.THaVDCPoint.fPhi" },
119 { "pt.paired", "Point is paired", "fPoints.THaVDCPoint.HasPartner()" },
120 { nullptr }
121 };
122 return DefineVarsFromList( vars, mode );
123}
124
125//_____________________________________________________________________________
127{
128 // Set debug level of this chamber and its plane subdetectors
129
131 fU->SetDebug(level);
132 fV->SetDebug(level);
133}
134
135//_____________________________________________________________________________
137 const THaVDCCluster* vcl ) const
138{
139 // Convert U,V coordinates of the given uv cluster pair to the detector
140 // coordinate system of this chamber. Takes u as the reference plane.
141
142 Double_t u = ucl->GetIntercept(); // Intercept in U plane
143 Double_t v0 = vcl->GetIntercept(); // Intercept in V plane
144 Double_t mu = ucl->GetSlope(); // Slope of U cluster
145 Double_t mv = vcl->GetSlope(); // Slope of V cluster
146
147 // Project v0 into the u plane
148 Double_t v = v0 - mv * GetSpacing();
149
150 PointCoords_t c{ UVtoX(u,v), UVtoY(u,v), UVtoX(mu,mv), UVtoY(mu,mv) };
151
152 return c;
153}
154
155//_____________________________________________________________________________
157{
158 // Match clusters in the U plane with cluster in the V plane by t0 and
159 // geometry. Fills fPoints with good pairs.
160
161 Int_t nu = fU->GetNClusters();
162 Int_t nv = fV->GetNClusters();
163
164 // Match best in-time clusters (t_0 close to zero).
165 // Discard any out-of-time clusters (|t_0| > N sigma, N = 3, configurable).
166
167 Double_t max_u_t0 = fU->GetT0Resolution();
168 Double_t max_v_t0 = fV->GetT0Resolution();
169
170 Int_t nuv = 0;
171
172 // Consider all possible uv cluster combinations
173 for( Int_t iu = 0; iu < nu; ++iu ) {
174 auto* uClust = fU->GetCluster(iu);
175 if( TMath::Abs(uClust->GetT0()) > max_u_t0 )
176 continue;
177
178 for( Int_t iv = 0; iv < nv; ++iv ) {
179 auto* vClust = fV->GetCluster(iv);
180 if( TMath::Abs(vClust->GetT0()) > max_v_t0 )
181 continue;
182
183 // Test position to be within drift chambers
184 const PointCoords_t c = CalcDetCoords(uClust,vClust);
185 if( !fU->IsInActiveArea( c.x, c.y ) ) {
186 continue;
187 }
188
189 // This cluster pair passes preliminary tests, so we save it, regardless
190 // of possible ambiguities, which will be sorted out later
191 new( (*fPoints)[nuv++] ) THaVDCPoint( uClust, vClust, this );
192 }
193 }
194
195 // return the number of UV pairs found
196 return nuv;
197}
198
199//_____________________________________________________________________________
201{
202 // Compute track info (x, y, theta, phi) for the all matched points.
203 // Uses TRANSPORT coordinates.
204
205 Int_t nPoints = GetNPoints();
206 for (int i = 0; i < nPoints; i++) {
207 THaVDCPoint* point = GetPoint(i);
208 if( point )
209 point->CalcDetCoords();
210 }
211 return nPoints;
212}
213
214//_____________________________________________________________________________
216{
217 // Clear event-by-event data
219 fU->Clear(opt);
220 fV->Clear(opt);
221 fPoints->Clear();
222}
223
224//_____________________________________________________________________________
226{
227 // Decode both wire planes
228
229 fU->Decode(evData);
230 fV->Decode(evData);
231 return 0;
232}
233
234//_____________________________________________________________________________
236{
237 // Apply time offset correction to all our hits
238
241}
242
243//_____________________________________________________________________________
245{
246 // Find clusters in U & V planes
247
248 fU->FindClusters();
249 fV->FindClusters();
250}
251
252//_____________________________________________________________________________
254{
255 // Fit data to recalculate cluster position
256
257 fU->FitTracks();
258 fV->FitTracks();
259}
260
261//_____________________________________________________________________________
263{
264 // Coarse computation of tracks
265
266 // Apply drift time offset correction obtained in prior Decode or
267 // InterStage(Decode) stage
269
270 // Find clusters and estimate their position/slope
271 FindClusters();
272
273 // Fit "local" tracks through the hit coordinates of each cluster
274 FitTracks();
275
276 // FIXME: The fit may fail, so we might want to call a function here
277 // that deletes points whose clusters have bad fits or are too small
278 // (e.g. 1 hit), etc.
279 // Right now, we keep them, preferring efficiency over precision.
280
281 // Pair U and V clusters
283
284 return 0;
285}
286
287//_____________________________________________________________________________
289{
290 // High precision computation of tracks
291
292 //TODO:
293 //- Redo time-to-distance conversion based on "global" track slope
294 //- Refit clusters of the track (if any), using new drift distances
295 //- Recalculate cluster UV coordinates, so new "global" slopes can be
296 // computed
297
298 return 0;
299}
300
301//_____________________________________________________________________________
int Int_t
#define c(i)
double Double_t
const char Option_t
Option_t Option_t TPoint TPoint const char mode
char name[80]
void Clear(Option_t *option="") override
virtual void SetDebug(Int_t level)
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 Int_t DefineVariables(EMode mode=kDefine)
const TVector3 & GetXax() const
Double_t GetYSize() const
Double_t fSize[3]
Double_t GetXSize() const
Double_t GetZSize() const
const TVector3 & GetYax() const
const TVector3 & GetOrigin() const
const TVector3 & GetZax() const
Double_t fCos_u
THaVDCChamber(const char *name="", const char *description="", THaDetectorBase *parent=nullptr)
virtual Int_t DefineVariables(EMode mode=kDefine)
Double_t fSin_v
virtual void Clear(Option_t *opt="")
TClonesArray * fPoints
THaVDCPlane * fU
Double_t fInv_sin_vu
Double_t UVtoX(Double_t u, Double_t v) const
virtual Int_t CoarseTrack()
Int_t MatchUVClusters()
THaVDCPlane * fV
THaVDCPoint * GetPoint(Int_t i) const
void ApplyTimeCorrection()
Int_t CalcPointCoords() const
Double_t fSin_u
Double_t fSpacing
Double_t GetSpacing() const
virtual void SetDebug(Int_t level)
virtual Int_t FineTrack()
Double_t fCos_v
virtual ~THaVDCChamber()
Int_t GetNPoints() const
virtual Int_t Decode(const THaEvData &evData)
PointCoords_t CalcDetCoords(const THaVDCCluster *u, const THaVDCCluster *v) const
Double_t UVtoY(Double_t u, Double_t v) const
Double_t GetIntercept() const
Double_t GetSlope() const
Double_t GetT0Resolution() const
Definition THaVDCPlane.h:72
virtual Int_t FitTracks()
THaVDCCluster * GetCluster(Int_t i) const
Definition THaVDCPlane.h:42
Int_t GetNWires() const
Definition THaVDCPlane.h:46
Double_t GetSinWAngle() const
Definition THaVDCPlane.h:66
virtual void Clear(Option_t *opt="")
Double_t GetWAngle() const
Definition THaVDCPlane.h:64
Double_t GetWSpac() const
Definition THaVDCPlane.h:63
virtual Int_t FindClusters()
Int_t GetNClusters() const
Definition THaVDCPlane.h:40
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)
Double_t GetZ() const
Definition THaVDCPlane.h:61
Double_t GetWBeg() const
Definition THaVDCPlane.h:62
Double_t GetCosWAngle() const
Definition THaVDCPlane.h:65
void CalcDetCoords()
void Clear(Option_t *option="") override
R__ALWAYS_INLINE Bool_t IsZombie() const
Double_t Abs(Double_t d)
Double_t Sin(Double_t)
Double_t Max(Double_t a, Double_t b)
v0
v
ClassImp(TPyArg)