Hall C ROOT/C++ Analyzer (hcana)
Loading...
Searching...
No Matches
THcDriftChamber.cxx
Go to the documentation of this file.
1
13#include "THcDriftChamber.h"
14#include "THcDCPlaneCluster.h"
15#include "THcDC.h"
16#include "THcDCHit.h"
17#include "THcGlobals.h"
18#include "THcParmList.h"
19#include "VarDef.h"
20#include "VarType.h"
21#include "THaTrack.h"
22#include "TClonesArray.h"
23#include "TMath.h"
24#include "TVectorD.h"
25#include "THcSpacePoint.h"
26#include "THaApparatus.h"
27
28#include "THaTrackProj.h"
29
30#include <cstring>
31#include <cstdio>
32#include <cstdlib>
33#include <iostream>
34#include <iomanip>
35
36using namespace std;
37
38//_____________________________________________________________________________
40 const char* name, const char* description,
41 const Int_t chambernum, THaDetectorBase* parent ) :
42 THaSubDetector(name,description,parent)
43{
44 // Constructor
45
46 // fTrackProj = new TClonesArray( "THaTrackProj", 5 );
47 fTrackProj = NULL;
48 fNPlanes = 0; // No planes until we make them
49
50 fChamberNum = chambernum;
51
52 fSpacePoints = new TClonesArray("THcSpacePoint",10);
53 fXPlaneClusters = new TClonesArray("THcDCPlaneCluster",10);
54 fUPlaneClusters = new TClonesArray("THcDCPlaneCluster",10);
55 fVPlaneClusters = new TClonesArray("THcDCPlaneCluster",10);
56 fUXPlaneClusters = new TClonesArray("THcDCPlaneCluster",10);
57 fVXPlaneClusters = new TClonesArray("THcDCPlaneCluster",10);
58
59
60 fHMSStyleChambers = 0; // Default
61}
62
63//_____________________________________________________________________________
66{
67 // Constructor
68 fTrackProj = NULL;
69 fSpacePoints = NULL;
70 fXPlaneClusters = NULL;
71 fUPlaneClusters = NULL;
72 fVPlaneClusters = NULL;
73 fUXPlaneClusters = NULL;
74 fVXPlaneClusters = NULL;
75 fIsInit = 0;
76
77}
78//_____________________________________________________________________________
79void THcDriftChamber::Setup(const char* name, const char* description)
80{
81
82}
83
84//_____________________________________________________________________________
86{
87 return 0;
88}
89
90//_____________________________________________________________________________
92{
93 // static const char* const here = "Init()";
94
96
97 EStatus status;
98 // This triggers call of ReadDatabase and DefineVariables
99 if( (status = THaSubDetector::Init( date )) )
100 return fStatus=status;
101
102 return fStatus = kOK;
103}
104
106{
107 plane->SetPlaneIndex(fNPlanes);
108 fPlanes.push_back(plane);
109 // HMS Specific
110 // Hard code Y plane numbers. Eventually need to get from database
111 if (fHMSStyleChambers) {
112 if(fChamberNum == 1) {
113 YPlaneNum=2;
114 YPlanePNum=5;
115 if(plane->GetPlaneNum() == 2) YPlaneInd = fNPlanes;
116 if(plane->GetPlaneNum() == 5) YPlanePInd = fNPlanes;
117 } else {
118 YPlaneNum=8;
119 YPlanePNum=11;
120 if(plane->GetPlaneNum() == 8) YPlaneInd = fNPlanes;
121 if(plane->GetPlaneNum() == 11) YPlanePInd = fNPlanes;
122 }
123 } else {
124 // SOS Specific
125 // Hard code X plane numbers. Eventually need to get from database
126 if(fChamberNum == 1) {
127 XPlaneNum=3;
128 XPlanePNum=4;
129 if(plane->GetPlaneNum() == 3) XPlaneInd = fNPlanes;
130 if(plane->GetPlaneNum() == 4) XPlanePInd = fNPlanes;
131 } else {
132 XPlaneNum=9;
133 XPlanePNum=10;
134 if(plane->GetPlaneNum() == 9) XPlaneInd = fNPlanes;
135 if(plane->GetPlaneNum() == 10) XPlanePInd = fNPlanes;
136 }
137 }
138 fNPlanes++;
139 return;
140}
141
142//_____________________________________________________________________________
144{
145
146 // cout << "THcDriftChamber::ReadDatabase()" << endl;
147 char prefix[2];
148 prefix[0]=tolower(GetApparatus()->GetName()[0]);
149 prefix[1]='\0';
150 DBRequest list[]={
151 {"_remove_sppt_if_one_y_plane",&fRemove_Sppt_If_One_YPlane, kInt,0,1},
152 {"dc_wire_velocity", &fWireVelocity, kDouble},
153 {"SmallAngleApprox", &fSmallAngleApprox, kInt,0,1},
154 {"stub_max_xpdiff", &fStubMaxXPDiff, kDouble,0,1},
155 {"debugflagpr", &fhdebugflagpr, kInt},
156 {"debugstubchisq", &fdebugstubchisq, kInt},
157 {Form("dc_%d_zpos",fChamberNum), &fZPos, kDouble},
158 {0}
159 };
162 TString SHMS="p";
163 TString HMS="h";
164 TString test=prefix[0];
165 if (test==SHMS ) fRatio_xpfp_to_xfp=0.0018; // SHMS
166 if (test == HMS ) fRatio_xpfp_to_xfp=0.0011; // HMS
167 fRemove_Sppt_If_One_YPlane = 0; // Default
168 fStubMaxXPDiff = 999.; //
169 gHcParms->LoadParmValues((DBRequest*)&list,prefix);
170 // Get parameters parent knows about
171 fParent = GetParent();
172 fMinHits = static_cast<THcDC*>(fParent)->GetMinHits(fChamberNum);
173 fMaxHits = static_cast<THcDC*>(fParent)->GetMaxHits(fChamberNum);
174 fMinCombos = static_cast<THcDC*>(fParent)->GetMinCombos(fChamberNum);
175 fFixPropagationCorrection = static_cast<THcDC*>(fParent)->GetFixPropagationCorrectionFlag();
176
177 fSpacePointCriterion = static_cast<THcDC*>(fParent)->GetSpacePointCriterion(fChamberNum);
178 fMaxDist = TMath::Sqrt(fSpacePointCriterion/2.0); // For easy space points
179
180 if (fhdebugflagpr) cout << " cham = " << fChamberNum << " Set yplane num " << YPlaneNum << " "<< YPlanePNum << endl;
181 // Generate the HAA3INV matrix for all the acceptable combinations
182 // of hit planes. Try to make it as generic as possible
183 // pindex=0 -> Plane 1 missing, pindex5 -> plane 6 missing. Won't
184 // replicate the exact values used in the ENGINE, because the engine
185 // had one big list of matrices for both chambers, while here we will
186 // have a list just for one chamber. Also, call pindex, pmindex as
187 // we tend to use pindex as a plane index.
188 fCosBeta = new Double_t [fNPlanes];
189 fSinBeta = new Double_t [fNPlanes];
190 fTanBeta = new Double_t [fNPlanes];
191 fSigma = new Double_t [fNPlanes];
192 fPsi0 = new Double_t [fNPlanes];
194 Int_t allplanes=0;
195 for(Int_t ip=0;ip<fNPlanes;ip++) {
196 fCosBeta[ip] = TMath::Cos(fPlanes[ip]->GetBeta());
197 fSinBeta[ip] = TMath::Sin(fPlanes[ip]->GetBeta());
198 fTanBeta[ip] = fSinBeta[ip]/fCosBeta[ip];
199 fSigma[ip] = fPlanes[ip]->GetSigma();
200 fPsi0[ip] = fPlanes[ip]->GetPsi0();
201 fStubCoefs[ip] = fPlanes[ip]->GetStubCoef();
202 allplanes |= 1<<ip;
203 }
204 // Unordered map introduced in C++-11
205 // Can use unordered_map if using C++-11
206 // May not want to use map a all for performance, but using it now
207 // for code clarity
208 for(Int_t ipm1=0;ipm1<fNPlanes+1;ipm1++) { // Loop over missing plane1
209 for(Int_t ipm2=ipm1;ipm2<fNPlanes+1;ipm2++) {
210 if(ipm1==ipm2 && ipm1<fNPlanes) continue;
211 TMatrixD AA3(3,3);
212 for(Int_t i=0;i<3;i++) {
213 for(Int_t j=i;j<3;j++) {
214 AA3[i][j] = 0.0;
215 for(Int_t ip=0;ip<fNPlanes;ip++) {
216 if(ipm1 != ip && ipm2 != ip) {
217 AA3[i][j] += fStubCoefs[ip][i]*fStubCoefs[ip][j];
218 }
219 }
220 AA3[j][i] = AA3[i][j];
221 }
222 }
223 Int_t bitpat = allplanes & ~(1<<ipm1) & ~(1<<ipm2);
224 // Should check that it is invertable
225 // if (fhdebugflagpr) cout << bitpat << " Determinant: " << AA3->Determinant() << endl;
226 AA3.Invert();
227 fAA3Inv[bitpat].ResizeTo(AA3);
228 fAA3Inv[bitpat] = AA3;
229 }
230 }
231
232 fIsInit = true;
233 return kOK;
234}
235
236//_____________________________________________________________________________
238{
239 // Initialize global variables and lookup table for decoder
240
241 if( mode == kDefine && fIsSetup ) return kOK;
242 fIsSetup = ( mode == kDefine );
243 // Register variables in global list
244
245 RVarDef vars[] = {
246 { "maxhits", "Maximum hits allowed", "fMaxHits" },
247 { "spacepoints", "Space points of DC", "fNSpacePoints" },
248 { "nhit", "Number of DC hits", "fNhits" },
249 { "trawhit", "Number of True Raw hits", "fN_True_RawHits" },
250 { "sp_nhits", "", "fSpacePoints.THcSpacePoint.GetNHits()" },
251 { "stub_x", "", "fSpacePoints.THcSpacePoint.GetStubX()" },
252 { "stub_xp", "", "fSpacePoints.THcSpacePoint.GetStubXP()" },
253 { "stub_y", "", "fSpacePoints.THcSpacePoint.GetStubY()" },
254 { "stub_yp", "", "fSpacePoints.THcSpacePoint.GetStubYP()" },
255 { "ncombos", "", "fSpacePoints.THcSpacePoint.GetCombos()" },
256 { "U_pos", "", "fUPlaneClusters.THcDCPlaneCluster.GetX()" },
257 { "X_pos", "", "fXPlaneClusters.THcDCPlaneCluster.GetX()" },
258 { "V_pos", "", "fVPlaneClusters.THcDCPlaneCluster.GetX()" },
259 { "UX_posx", "", "fUXPlaneClusters.THcDCPlaneCluster.GetX()" },
260 { "UX_posy", "", "fUXPlaneClusters.THcDCPlaneCluster.GetY()" },
261 { "VX_posx", "", "fVXPlaneClusters.THcDCPlaneCluster.GetX()" },
262 { "VX_posy", "", "fVXPlaneClusters.THcDCPlaneCluster.GetY()" },
263 { 0 }
264 };
265 DefineVarsFromList( vars, mode );
266 //
267 std::vector<RVarDef> ve;
268 ve.push_back( { "sphit", "", "fSpHit.SpNHits" });
269 ve.push_back( { "sphit_index", "", "fSpHit.SpHitIndex" });
270 ve.push_back({0}); // Needed to specify the end of list
271 return DefineVarsFromList( ve.data(), mode );
272
273 //return kOK;
274
275}
277{
278 // Make a list of hits for whole chamber
279 fNhits = 0;
280 fHits.clear();
281 fHits.reserve(40);
282 fSpHit.clear();
283 for(Int_t ip=0;ip<fNPlanes;ip++) {
284 TClonesArray* hitsarray = fPlanes[ip]->GetHits();
285 for(Int_t ihit=0;ihit<fPlanes[ip]->GetNHits();ihit++) {
286 fHits.push_back(static_cast<THcDCHit*>(hitsarray->At(ihit)));
287 fNhits++;
288 }
289 }
290 // if (fhdebugflagpr) cout << "ThcDriftChamber::ProcessHits() " << fNhits << " hits" << endl;
291}
292
293
295{
296 cout << " Num of nits = " << fNhits << endl;
297 cout << " Num " << " Plane " << " Wire " << " Wire-Center " << " RAW TDC " << " Drift time" << endl;
298 for(Int_t ihit=0;ihit<fNhits;ihit++) {
299 THcDCHit* thishit = fHits[ihit];
300 cout << ihit << " " <<thishit->GetPlaneNum() << " " << thishit->GetWireNum() << " " << thishit->GetPos() << " " << thishit->GetRawTime() << " " << thishit->GetTime() << endl;
301 }
302}
303
305{
306 //
315 if(fNhits >= fMinHits && fNhits < fMaxHits) {
316 if (fhdebugflagpr) cout << " in newfindspacepoints " << endl;
317 // first look for clusters pairs in X,U,V
318 for(Int_t ihit1=0;ihit1<fNhits;ihit1++) {
319 THcDCHit* hit1=fHits[ihit1];
320 THcDriftChamberPlane* plane1 = hit1->GetWirePlane();
321 for(Int_t ihit2=ihit1+1;ihit2<fNhits;ihit2++) {
322 THcDCHit* hit2=fHits[ihit2];
323 THcDriftChamberPlane* plane2 = hit2->GetWirePlane();
324 Double_t determinate = plane1->GetXsp()*plane2->GetYsp() - plane1->GetYsp()*plane2->GetXsp();
325 if (hit1->GetPlaneIndex() != hit2->GetPlaneIndex() &&
326 TMath::Abs(determinate) < 0.1 && abs(hit1->GetPos() - hit2->GetPos() ) < 0.6) {
327 auto idx = hit1->GetPlaneIndex();
328 if (idx >= 0 && idx < 6) {
329 THcDCPlaneCluster* clus = nullptr;;
333 clus->AddHit(hit1);
334 clus->AddHit(hit2);
335 if (fhdebugflagpr)
336 cout << ihit1 << " " << hit1->GetPos()<< " " << ihit2 << " " << hit2->GetPos() << " "
337 << plane1->GetYsp()<< " " << plane1->GetXsp()<< " "
338 << plane2->GetYsp()<< " " << plane2->GetXsp() << endl;
339 Double_t x = (hit1->GetPos() + hit2->GetPos())/2.;
340 Double_t y = 0.;
341 clus->SetXY(x,y);
342 hit1->IncrNPlaneClust();
343 hit2->IncrNPlaneClust();
344 }
345 }
346 }
347 if (hit1->GetNPlaneClust() == 0) {
348 if (fhdebugflagpr)
349 cout << " No match found" << endl;
350 auto idx = hit1->GetPlaneIndex();
351 if (idx >= 0 && idx < 6) {
352 THcDCPlaneCluster* clus = nullptr;;
356 clus->AddHit(hit1);
357 if (fhdebugflagpr)
358 cout << ihit1 << " " << hit1->GetPos()<< " " << plane1->GetYsp()<< " " << plane1->GetXsp()<< endl;
359 Double_t x = hit1->GetPos();
360 Double_t y = 0.;
361 clus->SetXY(x,y);
362 hit1->IncrNPlaneClust();
363 }
364 }
365 }
366 if (fhdebugflagpr) {
367 cout << " NUmber of U clusters = " << fNUPlaneClusters << endl;
368 for (Int_t nc=0;nc<fNUPlaneClusters;nc++) {
370 cout << " nc = " << nc << " nhits = " << clus->GetNHits() << endl;
371 }
372 cout << " NUmber of V clusters = " << fNVPlaneClusters<< endl;
373 for (Int_t nc=0;nc<fNVPlaneClusters;nc++) {
375 cout << " nc = " << nc << " nhits = " << clus->GetNHits()<< endl;
376 }
377 cout << " NUmber of X clusters = " << fNXPlaneClusters << endl;
378 for (Int_t nc=0;nc<fNXPlaneClusters;nc++) {
380 cout << " nc = " << nc << " nhits = " << clus->GetNHits()<< endl;
381 }
382 }
383 //
388 if (fNUPlaneClusters>0) {
389 for (Int_t ncu=0;ncu<fNUPlaneClusters;ncu++) {
390 THcDCPlaneCluster *uclus = static_cast<THcDCPlaneCluster*>(fUPlaneClusters->At(ncu));
391 THcDCHit* uhit = uclus->GetHit(0);
392 Double_t upos=uclus->GetX();
393 THcDriftChamberPlane* uplane = uhit->GetWirePlane();
394 for (Int_t ncx=0;ncx<fNXPlaneClusters;ncx++) {
395 THcDCPlaneCluster *xclus = static_cast<THcDCPlaneCluster*>(fXPlaneClusters->At(ncx));
396 THcDCHit* xhit = xclus->GetHit(0);
397 Double_t xpos=xclus->GetX();
398 if ((uclus->GetNHits()+xclus->GetNHits())>2) {
399 THcDriftChamberPlane* xplane = xhit->GetWirePlane();
400 Double_t determinate = uplane->GetXsp()*xplane->GetYsp()-uplane->GetYsp()*xplane->GetXsp();
401 Double_t x = (upos*xplane->GetYsp()- xpos*uplane->GetYsp())/determinate;
402 Double_t y = (xpos*uplane->GetXsp()- upos*xplane->GetXsp())/determinate;
404 for (Int_t ih=0;ih<uclus->GetNHits();ih++) {
405 THcDCHit* temphit = uclus->GetHit(ih);
406 clus->AddHit(temphit);
407 }
408 for (Int_t ih=0;ih<xclus->GetNHits();ih++) {
409 THcDCHit* temphit = xclus->GetHit(ih);
410 clus->AddHit(temphit);
411 }
412 clus->SetXY(x,y);
413 }
414 }
415 }
416 }
417 //
418 if (fNVPlaneClusters>0) {
419 for (Int_t ncv=0;ncv<fNVPlaneClusters;ncv++) {
420 THcDCPlaneCluster *vclus = static_cast<THcDCPlaneCluster*>(fVPlaneClusters->At(ncv));
421 THcDCHit* vhit = vclus->GetHit(0);
422 Double_t vpos=vclus->GetX();
423 THcDriftChamberPlane* vplane = vhit->GetWirePlane();
424 for (Int_t ncx=0;ncx<fNXPlaneClusters;ncx++) {
425 THcDCPlaneCluster *xclus = static_cast<THcDCPlaneCluster*>(fXPlaneClusters->At(ncx));
426 THcDCHit* xhit = xclus->GetHit(0);
427 Double_t xpos=xclus->GetX();
428 if ((vclus->GetNHits()+xclus->GetNHits())>2) {
429 THcDriftChamberPlane* xplane = xhit->GetWirePlane();
430 Double_t determinate = vplane->GetXsp()*xplane->GetYsp()-vplane->GetYsp()*xplane->GetXsp();
431 Double_t x = (vpos*xplane->GetYsp()- xpos*vplane->GetYsp())/determinate;
432 Double_t y = (xpos*vplane->GetXsp()- vpos*xplane->GetXsp())/determinate;
434 for (Int_t ih=0;ih<vclus->GetNHits();ih++) {
435 THcDCHit* temphit = vclus->GetHit(ih);
436 clus->AddHit(temphit);
437 }
438 for (Int_t ih=0;ih<xclus->GetNHits();ih++) {
439 THcDCHit* temphit = xclus->GetHit(ih);
440 clus->AddHit(temphit);
441 }
442 clus->SetXY(x,y);
443 }
444 }
445 }
446 }
447 //
448 if (fhdebugflagpr) {
449 cout << " NUmber of UX clusters = " << fNUXPlaneClusters << endl;
450 for (Int_t nc=0;nc<fNUXPlaneClusters;nc++) {
452 cout << " nc = " << nc << " nhits = " << clus->GetNHits() << endl;
453 }
454 cout << " NUmber of VX clusters = " << fNVXPlaneClusters<< endl;
455 for (Int_t nc=0;nc<fNVXPlaneClusters;nc++) {
457 cout << " nc = " << nc << " nhits = " << clus->GetNHits()<< endl;
458 }
459 }
460 //
461 vector <THcDCHit*> UX_uHits;
462 vector <THcDCHit*> UX_xHits;
463 vector <THcDCHit*> VX_vHits;
464 vector <THcDCHit*> VX_xHits;
465 for (Int_t nc=0;nc<fNUXPlaneClusters;nc++) {
466 UX_uHits.clear();
467 UX_xHits.clear();
468 THcDCPlaneCluster *uxclus = static_cast<THcDCPlaneCluster*>(fUXPlaneClusters->At(nc));
469 Double_t ux_posx = uxclus->GetX();
470 Double_t ux_posy = uxclus->GetY();
471 for (Int_t ih=0;ih<uxclus->GetNHits();ih++) {
472 THcDCHit* hit1 = uxclus->GetHit(ih);
473 if (hit1->GetPlaneIndex() ==0 || hit1->GetPlaneIndex() ==1) UX_uHits.push_back(hit1);
474 if (hit1->GetPlaneIndex() ==2 || hit1->GetPlaneIndex() ==3) UX_xHits.push_back(hit1);
475 }
476 for (Int_t nc2=0;nc2<fNVXPlaneClusters;nc2++) {
477 VX_vHits.clear();
478 VX_xHits.clear();
479 THcDCPlaneCluster *vxclus = static_cast<THcDCPlaneCluster*>(fVXPlaneClusters->At(nc2));
480 Double_t vx_posx = vxclus->GetX();
481 Double_t vx_posy = vxclus->GetY();
482 Double_t dist2= pow(ux_posx-vx_posx,2) + pow(ux_posy-vx_posy,2) ;
483 for (Int_t ih=0;ih<vxclus->GetNHits();ih++) {
484 THcDCHit* hit1 = vxclus->GetHit(ih);
485 if (hit1->GetPlaneIndex() ==4 || hit1->GetPlaneIndex() ==5) VX_vHits.push_back(hit1);
486 if (hit1->GetPlaneIndex() ==2 || hit1->GetPlaneIndex() ==3) VX_xHits.push_back(hit1);
487 }
488 Bool_t Xhits_match = kFALSE ;
489 if (fhdebugflagpr) cout << " vx_xhits = " << VX_xHits.size() << " ux_xhits = "<< UX_xHits.size() << endl;
490 if (VX_xHits.size() == UX_xHits.size() && UX_xHits.size()==1) {
491 if (VX_xHits[0] == UX_xHits[0]) Xhits_match = kTRUE;
492 } else if (VX_xHits.size() == UX_xHits.size() && UX_xHits.size()==2) {
493 if (VX_xHits[0] == UX_xHits[0] && VX_xHits[1] == UX_xHits[1]) Xhits_match = kTRUE;
494 if (VX_xHits[0] == UX_xHits[1] && VX_xHits[1] == UX_xHits[0]) Xhits_match = kTRUE;
495 }
496 if (fhdebugflagpr) cout << " Xhits_match = " << Xhits_match << " dist2 = " << dist2 << " spcrit = " << fSpacePointCriterion << " VX_x " << vx_posx<< " UX_x " << ux_posx<< " VX_y " << vx_posy<< " UX_y " << ux_posy << endl;
497 Int_t TotHits = UX_uHits.size()+UX_xHits.size()+VX_vHits.size();
498 if (dist2 <= fSpacePointCriterion && Xhits_match && TotHits>=fMinHits) {
499 if (fhdebugflagpr) cout << " Make spacepoint " << endl;
501 sp->Clear();
502 Double_t xt = (ux_posx + vx_posx)/2.;
503 Double_t yt = (ux_posy + vx_posy)/2.;
504 sp->SetXY(xt, yt);
505 sp->SetCombos(1);
506 for (UInt_t ih=0;ih<UX_uHits.size();ih++) { sp->AddHit(UX_uHits[ih]);}
507 for (UInt_t ih=0;ih<UX_xHits.size();ih++) { sp->AddHit(UX_xHits[ih]);}
508 for (UInt_t ih=0;ih<VX_vHits.size();ih++) { sp->AddHit(VX_vHits[ih]);}
509 fSpHit.SpNHits.push_back(sp->GetNHits());
510 for(Int_t ihit1=0;ihit1<fNhits;ihit1++) {
511 THcDCHit* hit1=fHits[ihit1];
512 for(Int_t isph=0;isph<sp->GetNHits();isph++) {
513 THcDCHit* hitsp=sp->GetHit(isph);
514 if (hitsp==hit1) fSpHit.SpHitIndex.push_back(ihit1);
515 }
516 }
517 }
518 }
519 }
520 //
521
522 //
523 } //fNhits >= fMinHits && fNhits < fMaxH
524 if (fhdebugflagpr) cout << " leave newfindspacepoints nsp = " << fNSpacePoints<< endl;
525
526 //
527 return(fNSpacePoints);
528}
529
531{
571 // fSpacePoints->Clear();
572 if (fhdebugflagpr) cout << " In old findspacepoint " << endl;
574
575 Int_t plane_hitind=0;
576 Int_t planep_hitind=0;
577
578
580 fEasySpacePoint = 0;
581 if(fNhits >= fMinHits && fNhits < fMaxHits) {
582 for(Int_t ihit=0;ihit<fNhits;ihit++) {
583 THcDCHit* thishit = fHits[ihit];
584 Int_t ip=thishit->GetPlaneNum(); // This is the absolute plane mumber
585 if(ip==YPlaneNum && fHMSStyleChambers) plane_hitind = ihit;
586 if(ip==YPlanePNum && fHMSStyleChambers) planep_hitind = ihit;
587 if(ip==XPlaneNum && !fHMSStyleChambers) plane_hitind = ihit;
588 if(ip==XPlanePNum && !fHMSStyleChambers) planep_hitind = ihit;
589 }
590 Int_t PlaneInd=0,PlanePInd=0;
591 if (fHMSStyleChambers) {
592 PlaneInd=YPlaneInd;
593 PlanePInd=YPlanePInd;
594 } else {
595 PlaneInd=XPlaneInd;
596 PlanePInd=XPlanePInd;
597 }
598 if(fPlanes[PlaneInd]->GetNHits() == 1 && fPlanes[PlanePInd]->GetNHits() == 1
599 && pow( (fHits[plane_hitind]->GetPos() - fHits[planep_hitind]->GetPos()),2)
601 && fNhits <= 6) { // An easy case, probably one hit per plane
602 if(fHMSStyleChambers) fEasySpacePoint = FindEasySpacePoint_HMS(plane_hitind, planep_hitind);
603 if(!fHMSStyleChambers) fEasySpacePoint = FindEasySpacePoint_SOS(plane_hitind, planep_hitind);
604 }
605 if(!fEasySpacePoint) { // It's not easy
607 }
608 // We have our space points for this chamber
609 if(fNSpacePoints > 0) {
612 // The routine is specific to HMS
613 //Int_t ndest=
614 DestroyPoorSpacePoints(); // Only for HMS?
615 // Loop over space points and remove those with less than 4 planes
616 // hit and missing hits in Y,Y' planes
617 }
618 Int_t nadded=SpacePointMultiWire(); // Only for HMS?
619 if (nadded) if (fhdebugflagpr) cout << nadded << " Space Points added with SpacePointMultiWire()" << endl;
620 }
623 // if(fNSpacePoints == 0) if (fhdebugflagpr) cout << "SelectSpacePoints() killed SP" << endl;
624 }
625 // if (fhdebugflagpr) cout << fNSpacePoints << " Space Points remain" << endl;
626 }
627 //
628 for (Int_t nsp=0;nsp<fNSpacePoints;nsp++) {
630 fSpHit.SpNHits.push_back(sp->GetNHits());
631 for(Int_t ihit1=0;ihit1<fNhits;ihit1++) {
632 THcDCHit* hit1=fHits[ihit1];
633 for(Int_t isph=0;isph<sp->GetNHits();isph++) {
634 THcDCHit* hitsp=sp->GetHit(isph);
635 if (hitsp==hit1) fSpHit.SpHitIndex.push_back(ihit1);
636 }
637 }
638 }
639
640 //
641 return(fNSpacePoints);
642}
643
644//_____________________________________________________________________________
645// HMS Specific
647{
657 Int_t easy_space_point=0;
658 Double_t yt = (fHits[yplane_hitind]->GetPos() + fHits[yplanep_hitind]->GetPos())/2.0;
659 Double_t xt = 0.0;
660 Int_t num_xhits = 0;
662
663 for(Int_t ihit=0;ihit<fNhits;ihit++) {
664 THcDCHit* thishit = fHits[ihit];
665 if(ihit!=yplane_hitind && ihit!=yplanep_hitind) { // x-like hit
666 // ysp and xsp are from h_generate_geometry
667 x_pos[ihit] = (thishit->GetPos()
668 -yt*thishit->GetWirePlane()->GetYsp())
669 /thishit->GetWirePlane()->GetXsp();
670 xt += x_pos[ihit];
671 num_xhits++;
672 } else {
673 x_pos[ihit] = 0.0;
674 }
675 }
676 xt = (num_xhits>0?xt/num_xhits:0.0);
677 easy_space_point = 1; // Assume we have an easy space point
678 // Rule it out if x points don't cluster well enough
679 for(Int_t ihit=0;ihit<fNhits;ihit++) {
680 if(ihit!=yplane_hitind && ihit!=yplanep_hitind) { // select x-like hit
681 if(TMath::Abs(xt-x_pos[ihit]) >= fMaxDist)
682 { easy_space_point=0; break;}
683 }
684 }
685 if(easy_space_point) { // Register the space point
687 sp->Clear();
688 sp->SetXY(xt, yt);
689 sp->SetCombos(0);
690 for(Int_t ihit=0;ihit<fNhits;ihit++) {
691 sp->AddHit(fHits[ihit]);
692 }
693 }
694 return(easy_space_point);
695}
696// SOS Specific
698{
708 Int_t easy_space_point=0;
709 Double_t xt = (fHits[xplane_hitind]->GetPos() + fHits[xplanep_hitind]->GetPos())/2.0;
710 Double_t yt = 0.0;
711 Int_t num_yhits = 0;
713
714 for(Int_t ihit=0;ihit<fNhits;ihit++) {
715 THcDCHit* thishit = fHits[ihit];
716 if(ihit!=xplane_hitind && ihit!=xplanep_hitind) { // y-like hit
717 // ysp and xsp are from h_generate_geometry
718 y_pos[ihit] = (thishit->GetPos()
719 -xt*thishit->GetWirePlane()->GetXsp())
720 /thishit->GetWirePlane()->GetYsp();
721 yt += y_pos[ihit];
722 num_yhits++;
723 } else {
724 y_pos[ihit] = 0.0;
725 }
726 }
727 yt = (num_yhits>0?yt/num_yhits:0.0);
728 easy_space_point = 1; // Assume we have an easy space point
729 // Rule it out if x points don't cluster well enough
730 for(Int_t ihit=0;ihit<fNhits;ihit++) {
731 if(ihit!=xplane_hitind && ihit!=xplanep_hitind) { // select y-like hit
732 if(TMath::Abs(yt-y_pos[ihit]) >= fMaxDist)
733 { easy_space_point=0; break;}
734 }
735 }
736 if(easy_space_point) { // Register the space point
738 sp->Clear();
739 sp->SetXY(xt, yt);
740 sp->SetCombos(0);
741 for(Int_t ihit=0;ihit<fNhits;ihit++) {
742 sp->AddHit(fHits[ihit]);
743 }
744 }
745 return(easy_space_point);
746}
747
748//_____________________________________________________________________________
749// Generic
751{
752 Int_t MAX_NUMBER_PAIRS=1000; // Where does this get set?
753 struct Pair {
754 THcDCHit* hit1;
755 THcDCHit* hit2;
756 Double_t x, y;
757 };
758 Pair pairs[MAX_NUMBER_PAIRS];
759 //
760 Int_t ntest_points=0;
761 for(Int_t ihit1=0;ihit1<fNhits-1;ihit1++) {
762 THcDCHit* hit1=fHits[ihit1];
763 THcDriftChamberPlane* plane1 = hit1->GetWirePlane();
764 for(Int_t ihit2=ihit1+1;ihit2<fNhits;ihit2++) {
765 if(ntest_points < MAX_NUMBER_PAIRS) {
766 THcDCHit* hit2=fHits[ihit2];
767 THcDriftChamberPlane* plane2 = hit2->GetWirePlane();
768 Double_t determinate = plane1->GetXsp()*plane2->GetYsp()
769 -plane1->GetYsp()*plane2->GetXsp();
770 if(TMath::Abs(determinate) > 0.3) { // 0.3 is sin(alpha1-alpha2)=sin(17.5)
771 pairs[ntest_points].hit1 = hit1;
772 pairs[ntest_points].hit2 = hit2;
773 pairs[ntest_points].x = (hit1->GetPos()*plane2->GetYsp()
774 - hit2->GetPos()*plane1->GetYsp())
775 /determinate;
776 pairs[ntest_points].y = (hit2->GetPos()*plane1->GetXsp()
777 - hit1->GetPos()*plane2->GetXsp())
778 /determinate;
779 ntest_points++;
780 }
781 }
782 }
783 }
784 Int_t ncombos=0;
785 struct Combo {
786 Pair* pair1;
787 Pair* pair2;
788 };
789 Combo combos[10*MAX_NUMBER_PAIRS];
790 for(Int_t ipair1=0;ipair1<ntest_points-1;ipair1++) {
791 for(Int_t ipair2=ipair1+1;ipair2<ntest_points;ipair2++) {
792 if(ncombos < 10*MAX_NUMBER_PAIRS) {
793 Double_t dist2 = pow(pairs[ipair1].x - pairs[ipair2].x,2)
794 + pow(pairs[ipair1].y - pairs[ipair2].y,2);
795 if(dist2 <= fSpacePointCriterion) {
796 combos[ncombos].pair1 = &pairs[ipair1];
797 combos[ncombos].pair2 = &pairs[ipair2];
798 ncombos++;
799 }
800 }
801 }
802 }
803 // Loop over all valid combinations and build space points
804 //if (fhdebugflagpr) cout << "looking for hard Space Point combos = " << ncombos << endl;
805 for(Int_t icombo=0;icombo<ncombos;icombo++) {
806 THcDCHit* hits[4];
807 hits[0]=combos[icombo].pair1->hit1;
808 hits[1]=combos[icombo].pair1->hit2;
809 hits[2]=combos[icombo].pair2->hit1;
810 hits[3]=combos[icombo].pair2->hit2;
811 // Get Average Space point xt, yt
812 Double_t xt = (combos[icombo].pair1->x + combos[icombo].pair2->x)/2.0;
813 Double_t yt = (combos[icombo].pair1->y + combos[icombo].pair2->y)/2.0;
814 // Loop over space points
815
816 if(fNSpacePoints > 0) {
817 Int_t add_flag=1;
818 for(Int_t ispace=0;ispace<fNSpacePoints;ispace++) {
819 THcSpacePoint* sp = (THcSpacePoint*)(*fSpacePoints)[ispace];
820 if(sp->GetNHits() > 0) {
821 Double_t sqdist_test = pow(xt - sp->GetX(),2) + pow(yt - sp->GetY(),2);
822 // I (who is I) want to be careful if sqdist_test is bvetween 1 and
823 // 3 fSpacePointCriterion. Let me ignore not add a new point the
824 if(sqdist_test < 3*fSpacePointCriterion) {
825 add_flag = 0; // do not add a new space point
826 }
827 if(sqdist_test < fSpacePointCriterion) {
828 // This is a real match
829 // Add the new hits to the existing space point
830 Int_t iflag[4];
831 iflag[0]=0;iflag[1]=0;iflag[2]=0;iflag[3]=0;
832 // Find out which of the four hits in the combo are already
833 // in the space point under consideration so that we don't
834 // add duplicate hits to the space point
835 for(Int_t isp_hit=0;isp_hit<sp->GetNHits();isp_hit++) {
836 for(Int_t icm_hit=0;icm_hit<4;icm_hit++) { // Loop over combo hits
837 if(sp->GetHit(isp_hit)==hits[icm_hit]) {
838 iflag[icm_hit] = 1;
839 }
840 }
841 }
842 // Remove duplicated pionts in the combo so we don't add
843 // duplicate hits to the space point
844 for(Int_t icm1=0;icm1<3;icm1++) {
845 for(Int_t icm2=icm1+1;icm2<4;icm2++) {
846 if(hits[icm1]==hits[icm2]) {
847 iflag[icm2] = 1;
848 }
849 }
850 }
851 // Add the unique combo hits to the space point
852 for(Int_t icm=0;icm<4;icm++) {
853 if(iflag[icm]==0) {
854 sp->AddHit(hits[icm]);
855 }
856 }
857 sp->IncCombos();
858 // cout << " number of combos = " << sp->GetCombos() << endl;
859 // Terminate loop since this combo can only belong to one space point
860 break;
861 }
862 }
863 }// End of loop over existing space points
864 // Create a new space point if more than 2*space_point_criteria
866 if(add_flag) {
867 //if (fhdebugflagpr) cout << " add glag = " << add_flag << " space pts = " << fNSpacePoints << endl ;
869 sp->Clear();
870 sp->SetXY(xt, yt);
871 sp->SetCombos(1);
872 sp->AddHit(hits[0]);
873 sp->AddHit(hits[1]);
874 if(hits[0] != hits[2] && hits[1] != hits[2]) {
875 sp->AddHit(hits[2]);
876 }
877 if(hits[0] != hits[3] && hits[1] != hits[3]) {
878 sp->AddHit(hits[3]);
879 }
880 }
881 }
882 } else {// Create first space point
883 // This duplicates code above. Need to see if we can restructure
884 // to avoid
886 sp->Clear();
887 sp->SetXY(xt, yt);
888 sp->SetCombos(1);
889 sp->AddHit(hits[0]);
890 sp->AddHit(hits[1]);
891 if(hits[0] != hits[2] && hits[1] != hits[2]) {
892 sp->AddHit(hits[2]);
893 }
894 if(hits[0] != hits[3] && hits[1] != hits[3]) {
895 sp->AddHit(hits[3]);
896 }
897 //if (fhdebugflagpr) cout << "1st hard Space Point " << xt << " " << yt << " Space point # =" << fNSpacePoints << " combos = " << sp->GetCombos() << endl;
898 }//End check on 0 space points
899 }//End loop over combos
900 //if (fhdebugflagpr) cout << " finished findspacept # of sp pts = " << fNSpacePoints << endl;
901 return(fNSpacePoints);
902}
903
904//_____________________________________________________________________________
905// HMS Specific?
907{
908 Int_t nhitsperplane[fNPlanes];
909
910 Int_t spacepointsgood[fNSpacePoints];
911 Int_t ngood=0;
912
913 for(Int_t i=0;i<fNSpacePoints;i++) {
914 spacepointsgood[i] = 0;
915 }
916 for(Int_t isp=0;isp<fNSpacePoints;isp++) {
917 Int_t nplanes_hit = 0;
918 for(Int_t ip=0;ip<fNPlanes;ip++) {
919 nhitsperplane[ip] = 0;
920 }
921 // Count # hits in each plane for this space point
923 for(Int_t ihit=0;ihit<sp->GetNHits();ihit++) {
924 THcDCHit* hit=sp->GetHit(ihit);
925 // hit_order(hit) = ihit;
926 Int_t ip = hit->GetPlaneIndex();
927 nhitsperplane[ip]++;
928 }
929 // Count # planes that have hits
930 for(Int_t ip=0;ip<fNPlanes;ip++) {
931 if(nhitsperplane[ip] > 0) {
932 nplanes_hit++;
933 }
934 }
935 if(nplanes_hit >= fMinHits && nhitsperplane[YPlaneInd]>0
936 && nhitsperplane[YPlanePInd] > 0) {
937 spacepointsgood[ngood++] = isp; // Build list of good points
938 } else {
939 // if (fhdebugflagpr) cout << "Missing Y-hit!!";
940 }
941 }
942
943 // Remove the bad space points
944 Int_t nremoved=fNSpacePoints-ngood;
945 fNSpacePoints = ngood;
946 for(Int_t isp=0;isp<fNSpacePoints;isp++) { // New index num ber
947 Int_t osp=spacepointsgood[isp]; // Original index number
948 if(osp > isp) {
949 // Does this work, or do we have to copy each member?
950 // If it doesn't we should overload the = operator
951 //(*fSpacePoints)[isp] = (*fSpacePoints)[osp];
954 spi->Clear();
955 Double_t xt,yt;
956 xt=spo->GetX();
957 yt=spo->GetY();
958 spi->SetXY(xt, yt);
959 for(Int_t ihit=0;ihit<spo->GetNHits();ihit++) {
960 THcDCHit* hit = spo->GetHit(ihit);
961 spi->AddHit(hit);
962 }
963 }
964 }
965 return nremoved;
966}
967
968//_____________________________________________________________________________
969// HMS Specific?
971{
984 Int_t nhitsperplane[fNPlanes];
986
987 Int_t nsp_check;
988 //Int_t nplanes_single;
989
990 Int_t nsp_tot=fNSpacePoints;
991 Int_t nsp_totl=fNSpacePoints;
992 //if (fhdebugflagpr) cout << "Start Multiwire # of sp pts = " << nsp_totl << endl;
993
994 for(Int_t isp=0;isp<nsp_totl;isp++) {
995 Int_t nplanes_hit = 0; // Number of planes with hits
996 Int_t nplanes_mult = 0; // Number of planes with multiple hits
997 Int_t nsp_new = 1;
998 Int_t newsp_num=0;
999 //if (fhdebugflagpr) cout << "Looping thru space pts at # = " << isp << " total = " << fNSpacePoints << endl;
1000
1001 for(Int_t ip=0;ip<fNPlanes;ip++) {
1002 nhitsperplane[ip] = 0;
1003 for(Int_t ih=0;ih<MAX_HITS_PER_POINT;ih++) {
1004 hits_plane[ip][ih] = 0;
1005 }
1006 }
1007 // Sort Space Points hits by plane
1009 for(Int_t ihit=0;ihit<sp->GetNHits();ihit++) { // All hits in SP
1010 THcDCHit* hit=sp->GetHit(ihit);
1011 // hit_order Make a hash
1012 // hash(hit) = ihit;
1013 Int_t ip = hit->GetPlaneIndex();
1014 hits_plane[ip][nhitsperplane[ip]++] = hit;
1015 //if (fhdebugflagpr) cout << " hit = " << ihit+1 << " plane index = " << ip << " nhitsperplane = " << nhitsperplane[ip] << endl;
1016 }
1017 for(Int_t ip=0;ip<fNPlanes;ip++) {
1018 if(nhitsperplane[ip] > 0) {
1019 nplanes_hit++;
1020 nsp_new *= nhitsperplane[ip];
1021 if(nhitsperplane[ip] > 1) nplanes_mult++;
1022 //if (fhdebugflagpr) cout << "Found plane with multi hits plane =" << ip+1 << " nplane_hit = "<< nplanes_hit << " nsp_new = " <<nsp_new << " nplane_mult = "<< nplanes_mult << endl;
1023 }
1024 }
1025 --nsp_new;
1026 nsp_check=nsp_tot + nsp_new;
1027 //nplanes_single = nplanes_hit - nplanes_mult;
1028 //if (fhdebugflagpr) cout << " # of new space points = " << nsp_new << " total now = " << nsp_tot<< endl;
1029 // Check if cloning conditions are met
1030 Int_t ntot = 0;
1031 if(nplanes_hit >= 4 && nplanes_mult < 4 && nplanes_mult >0
1032 && nsp_check < 20) {
1033 //if (fhdebugflagpr) cout << " Cloning space point " << endl;
1034 // Order planes by decreasing # of hits
1035
1036 Int_t maxplane[fNPlanes];
1037 for(Int_t ip=0;ip<fNPlanes;ip++) {
1038 maxplane[ip] = ip;
1039 }
1040 // Sort by decreasing # of hits
1041 for(Int_t ip1=0;ip1<fNPlanes-1;ip1++) {
1042 for(Int_t ip2=ip1+1;ip2<fNPlanes;ip2++) {
1043 if(nhitsperplane[maxplane[ip2]] > nhitsperplane[maxplane[ip1]]) {
1044 Int_t temp = maxplane[ip1];
1045 maxplane[ip1] = maxplane[ip2];
1046 maxplane[ip2] = temp;
1047 }
1048 }
1049 }
1050 // First fill clones with 1 hit each from the 3 planes with the most hits
1051 for(Int_t n1=0;n1<nhitsperplane[maxplane[0]];n1++) {
1052 for(Int_t n2=0;n2<nhitsperplane[maxplane[1]];n2++) {
1053 for(Int_t n3=0;n3<nhitsperplane[maxplane[2]];n3++) {
1054 ntot++;
1055 newsp_num = fNSpacePoints; //
1056 //if (fhdebugflagpr) cout << " new space pt num = " << newsp_num << " " << fNSpacePoints << endl;
1057 //THcSpacePoint* newsp;
1058 if(n1==0 && n2==0 && n3==0) {
1059 newsp_num = isp; // Copy over the original SP
1060 THcSpacePoint* newsp = (THcSpacePoint*)fSpacePoints->ConstructedAt(newsp_num);//= (THcSpacePoint*)(*fSpacePoints)[newsp_num];
1061 //if (fhdebugflagpr) cout << " Copy over original SP " << endl;
1062 // newsp = sp;
1063 Int_t combos_save=sp->GetCombos();
1064 newsp->Clear(); // Clear doesn't clear X, Y
1065 // if (fhdebugflagpr) cout << " original sp #hits combos X y " << sp->GetCombos() << sp->GetNHits() << sp->GetX() << " " << sp->GetY() << endl;
1066 newsp->SetXY(sp->GetX(), sp->GetY());
1067 newsp->SetCombos(combos_save);
1068 newsp->AddHit(hits_plane[maxplane[0]][n1]);
1069 newsp->AddHit(hits_plane[maxplane[1]][n2]);
1070 newsp->AddHit(hits_plane[maxplane[2]][n3]);
1071 newsp->AddHit(hits_plane[maxplane[3]][0]);
1072 if(nhitsperplane[maxplane[4]] == 1) {
1073 newsp->AddHit(hits_plane[maxplane[4]][0]);
1074 if(nhitsperplane[maxplane[5]] == 1)
1075 newsp->AddHit(hits_plane[maxplane[5]][0]);
1076 }
1077 } else {
1078 // if (fhdebugflagpr) cout << " setting other sp " << "# space pts now = " << fNSpacePoints << endl;
1080 fNSpacePoints++;
1081 Int_t combos_save=sp->GetCombos();
1082 newsp->Clear();
1083 newsp->SetXY(sp->GetX(), sp->GetY());
1084 newsp->SetCombos(combos_save);
1085 newsp->AddHit(hits_plane[maxplane[0]][n1]);
1086 newsp->AddHit(hits_plane[maxplane[1]][n2]);
1087 newsp->AddHit(hits_plane[maxplane[2]][n3]);
1088 newsp->AddHit(hits_plane[maxplane[3]][0]);
1089 if(nhitsperplane[maxplane[4]] == 1) {
1090 newsp->AddHit(hits_plane[maxplane[4]][0]);
1091 if(nhitsperplane[maxplane[5]] == 1)
1092 newsp->AddHit(hits_plane[maxplane[5]][0]);
1093 }
1094 }
1095 }
1096 }
1097 }
1098 // In the ENGINE, we loop over the clones and order the hits in the
1099 // same order as the parent SP. It is not done here because it is a little
1100 // tricky. Is it necessary?
1101 nsp_tot += (ntot-1);
1102 } else {
1103 ntot=1; // space point not to be cloned
1104 }
1105 }
1106 assert (nsp_tot > 0); // program terminates if nsp_tot <=0
1107 Int_t nadded=0;
1108 if(nsp_tot <= 20) {
1109 nadded = nsp_tot - fNSpacePoints;
1110 // fNSpacePoints = nsp_tot;
1111 }
1112 //if (fhdebugflagpr) cout << " Added space pts " << nadded << " total space pts = " << fNSpacePoints << endl;
1113
1114 // In Fortran, fill in zeros.
1115 return(nadded);
1116}
1117
1118//_____________________________________________________________________________
1119// Generic
1121{
1126 for(Int_t isp=0;isp<fNSpacePoints;isp++) {
1128 Int_t startnum = sp->GetNHits();
1129 Int_t goodhit[startnum];
1130 for(Int_t ihit=0;ihit<startnum;ihit++) {
1131 goodhit[ihit] = 1;
1132 }
1133 // For each plane, mark all hits longer than the shortest drift time
1134 for(Int_t ihit1=0;ihit1<startnum-1;ihit1++) {
1135 THcDCHit* hit1 = sp->GetHit(ihit1);
1136 Int_t plane1=hit1->GetPlaneIndex();
1137 Double_t tdrift1 = hit1->GetTime();
1138 for(Int_t ihit2=ihit1+1;ihit2<startnum;ihit2++) {
1139 THcDCHit* hit2 = sp->GetHit(ihit2);
1140 Int_t plane2=hit2->GetPlaneIndex();
1141 Double_t tdrift2 = hit2->GetTime();
1142 if(plane1 == plane2) {
1143 if(tdrift1 > tdrift2) {
1144 goodhit[ihit1] = 0;
1145 } else {
1146 goodhit[ihit2] = 0;
1147 }
1148 }
1149 }
1150 }
1151 // Gather the remaining hits
1152 Int_t finalnum = 0;
1153 for(Int_t ihit=0;ihit<startnum;ihit++) {
1154 if(goodhit[ihit] > 0) { // Keep this hit
1155 if (ihit > finalnum) { // Move hit
1156 sp->ReplaceHit(finalnum++, sp->GetHit(ihit));
1157 } else {
1158 finalnum++ ;
1159 }
1160 }
1161 }
1162 sp->SetNHits(finalnum);
1163 // if (fhdebugflagpr) cout << " choose single hit start # of hits = " << startnum << " final # = " <<finalnum << endl;
1164 }
1165}
1166//_____________________________________________________________________________
1167// Generic
1169{
1176 Int_t sp_count=0;
1177 //
1178 //
1179 for(Int_t isp=0;isp<fNSpacePoints;isp++) {
1180 // Include fEasySpacePoint because ncombos not filled in
1182 if(sp->GetCombos() >= fMinCombos || fEasySpacePoint) {
1183 if(sp->GetNHits() >= fMinHits) {
1184 if(isp > sp_count) {
1185 // (*fSpacePoints)[sp_count] = (*fSpacePoints)[isp];
1186 THcSpacePoint* sp1 = (THcSpacePoint*)(*fSpacePoints)[sp_count];
1187 sp1->Clear();
1188 Double_t xt,yt;
1189 xt=sp->GetX();
1190 yt=sp->GetY();
1191 sp1->SetXY(xt, yt);
1192 sp1->SetCombos(sp->GetCombos());
1193 for(Int_t ihit=0;ihit<sp->GetNHits();ihit++) {
1194 THcDCHit* hit = sp->GetHit(ihit);
1195 sp1->AddHit(hit);
1196 }
1197 }
1198 sp_count++;
1199 }
1200 }
1201 }
1202 fNSpacePoints = sp_count;
1203 //for(Int_t isp=0;isp<fNSpacePoints;isp++) {
1204 // THcSpacePoint* sp = (THcSpacePoint*)(*fSpacePoints)[isp];
1205 //if (fhdebugflagpr) cout << " sp pt = " << isp+1 << " # of hits = " << sp->GetNHits() << endl;
1206 //for(Int_t ihit=0;ihit<sp->GetNHits();ihit++) {
1207 //THcDCHit* hit = sp->GetHit(ihit);
1208 //THcDriftChamberPlane* plane=hit->GetWirePlane();
1209 // if (fhdebugflagpr) cout << ihit+1 << "selecting " << plane->GetPlaneNum() << " " << plane->GetChamberNum() << " " << hit->GetTime() << " " << hit->GetDist() << " " << plane->GetCentralTime() << " " << plane->GetDriftTimeSign() << endl;
1210 // }
1211 // }
1212}
1213
1215{
1228 //if (fhdebugflagpr) cout << "In correcthittimes fNSpacePoints = " << fNSpacePoints << endl;
1229
1230 for(Int_t isp=0;isp<fNSpacePoints;isp++) {
1232 Double_t x = sp->GetX();
1233 Double_t y = sp->GetY();
1234 for(Int_t ihit=0;ihit<sp->GetNHits();ihit++) {
1235 THcDCHit* hit = sp->GetHit(ihit);
1236 THcDriftChamberPlane* plane=hit->GetWirePlane();
1237
1238 // How do we know this correction only gets applied once? Is
1239 // it determined that a given hit can only belong to one space point?
1240 Double_t time_corr = plane->GetReadoutX() ?
1241 y*plane->GetReadoutCorr()/fWireVelocity :
1243
1244 // This applies the wire velocity correction for new SHMS chambers --hszumila, SEP17
1245 if (!fHMSStyleChambers){
1246 Int_t pln = hit->GetPlaneNum();
1247 Int_t readoutSide = hit->GetReadoutSide();
1248
1249 Double_t posn = hit->GetPos();
1250 //The following values are determined from param file as permutations on planes 5 and 10
1251 Int_t readhoriz = plane->GetReadoutLR();
1252 Int_t readvert = plane->GetReadoutTB();
1253
1254 //+x is up and +y is beam right!
1255 double alpha = static_cast<THcDC*>(fParent)->GetAlphaAngle(pln);
1256 double xc = posn*TMath::Sin(alpha);
1257 double yc = posn*TMath::Cos(alpha);
1258
1259 Double_t wireDistance = plane->GetReadoutX() ?
1260 (abs(y-yc))*abs(plane->GetReadoutCorr()) :
1261 (abs(x-xc))*abs(plane->GetReadoutCorr());
1262
1263 //Readout side is based off wiring diagrams
1264 switch (readoutSide){
1265 case 1: //readout from top of chamber
1266 if (x>xc){wireDistance = -readvert*wireDistance;}
1267 else{wireDistance = readvert*wireDistance;}
1268
1269 break;
1270 case 2://readout from right of chamber
1271 if (y>yc){wireDistance = -readhoriz*wireDistance;}
1272 else{wireDistance = readhoriz*wireDistance;}
1273
1274 break;
1275 case 3: //readout from bottom of chamber
1276 if (xc>x){wireDistance= -readvert*wireDistance;}
1277 else{wireDistance = readvert*wireDistance;}
1278
1279 break;
1280 case 4: //readout from left of chamber
1281 if(yc>y){wireDistance = -readhoriz*wireDistance;}
1282 else{wireDistance = readhoriz*wireDistance;}
1283
1284 break;
1285 default:
1286 wireDistance = 0.0;
1287 }
1288
1289 Double_t timeCorrection = wireDistance/fWireVelocity;
1290
1291 if(fFixPropagationCorrection==0) { // ENGINE behavior
1292 Double_t time=hit->GetTime();
1293 hit->SetTime(time - timeCorrection);
1294 hit->ConvertTimeToDist();
1295 } else {
1296 Double_t time=hit->GetTime();
1297 Double_t dist=hit->GetDist();
1298
1299 hit->SetTime(time - timeCorrection);
1300 //double usingOldTime = time-time_corr;
1301 //hit->SetTime(time- time_corr);
1302
1303 hit->ConvertTimeToDist();
1304 sp->SetHitDist(ihit, hit->GetDist());
1305
1306 hit->SetTime(time); // Restore time
1307 hit->SetDist(dist); // Restore distance
1308 }
1309
1310 } else {
1312
1313 // if (fhdebugflagpr) cout << "Correcting hit " << hit << " " << plane->GetPlaneNum() << " " << isp << "/" << ihit << " " << x << "," << y << endl;
1314 // Fortran ENGINE does not do this check, so hits can get "corrected"
1315 // multiple times if they belong to multiple space points.
1316 // To reproduce the precise ENGINE behavior, remove this if condition.
1317 if(fFixPropagationCorrection==0) { // ENGINE behavior
1318 hit->SetTime(hit->GetTime() - plane->GetCentralTime()
1319 + plane->GetDriftTimeSign()*time_corr);
1320 hit->ConvertTimeToDist();
1321 // hit->SetCorrectedStatus(1);
1322 } else {
1323 // New behavior: Save corrected distance with the hit in the space point
1324 // so that the same hit can have a different correction depending on
1325 // which space point it is in.
1326 //
1327 // This is a hack now because the converttimetodist method is connected to the hit
1328 // so I compute the corrected time and distance, and then restore the original
1329 // time and distance. Can probably add a method to hit that does a conversion on a time
1330 // but does not modify the hit data.
1331 Double_t time=hit->GetTime();
1332 Double_t dist=hit->GetDist();
1333 hit->SetTime(time - plane->GetCentralTime()
1334 + plane->GetDriftTimeSign()*time_corr);
1335 hit->ConvertTimeToDist();
1336 sp->SetHitDist(ihit, hit->GetDist());
1337 hit->SetTime(time); // Restore time
1338 hit->SetDist(dist); // Restore distance
1339 }
1340 }
1341 }
1342 }
1343}
1345// From http://graphics.stanford.edu/~seander/bithacks.html
1346{
1347 x = x - ((x >> 1) & 0x55555555);
1348 x = (x & 0x33333333) + ((x >> 2) & 0x33333333);
1349 return (((x + (x >> 4)) & 0x0F0F0F0F) * 0x01010101) >> 24;
1350}
1351
1352//_____________________________________________________________________________
1354{
1361 for(Int_t isp=0; isp<fNSpacePoints; isp++) {
1362 // Build a bit pattern of which planes are hit
1364 Int_t nhits = sp->GetNHits();
1365 UInt_t bitpat = 0; // Bit pattern of which planes are hit
1366 Double_t maxchi2= 1.0e10;
1367 Double_t minchi2 = maxchi2;
1368 Double_t tmp_minchi2=maxchi2;
1369 Double_t minxp = 0.25;
1370 Int_t hasy1 = -1;
1371 Int_t hasy2 = -1;
1372 Int_t plusminusknown[nhits];
1373 Int_t plusminusbest[nhits];
1374 Int_t plusminus[nhits]; // ENGINE makes this array float. Why?
1375 Int_t tmp_plusminus[nhits];
1376 Int_t plane_list[nhits];
1377 Double_t stub[4];
1378 Double_t tmp_stub[4];
1379 Int_t nplusminus;
1380
1381 if(nhits < 0) {
1382 if (fhdebugflagpr) cout << "THcDriftChamber::LeftRight() nhits < 0" << endl;
1383 } else if (nhits==0) {
1384 if (fhdebugflagpr) cout << "THcDriftChamber::LeftRight() nhits = 0" << endl;
1385 }
1386 for(Int_t ihit=0;ihit < nhits;ihit++) {
1387 THcDCHit* hit = sp->GetHit(ihit);
1388 Int_t pindex = hit->GetPlaneIndex();
1389 plane_list[ihit] = pindex;
1390
1391 bitpat |= 1<<pindex;
1392
1393 plusminusknown[ihit] = 0;
1394
1395 if(pindex == YPlaneInd) hasy1 = ihit;
1396 if(pindex == YPlanePInd) hasy2 = ihit;
1397 }
1398 nplusminus = 1<<nhits;
1399 if(fHMSStyleChambers) {
1400 Int_t smallAngOK = (hasy1>=0) && (hasy2>=0);
1401 if(fSmallAngleApprox !=0 && smallAngOK) { // to small Angle L/R for Y,Y' planes
1402
1403 if(sp->GetHit(hasy2)->GetPos() <=
1404 sp->GetHit(hasy1)->GetPos()) {
1405 plusminusknown[hasy1] = -1;
1406 plusminusknown[hasy2] = 1;
1407 } else {
1408 plusminusknown[hasy1] = 1;
1409 plusminusknown[hasy2] = -1;
1410 }
1411 nplusminus = 1<<(nhits-2);
1412 // if (fhdebugflagpr) cout << " Small angle approx = " << smallAngOK << " " << plusminusknown[hasy1] << endl;
1413 //if (fhdebugflagpr) cout << "pm = " << plusminusknown[hasy2] << " " << hasy1 << " " << hasy2 << endl;
1414 //if (fhdebugflagpr) cout << " Plane index " << YPlaneInd << " " << YPlanePInd << endl;
1415 }
1416 } else { // SOS Style
1417 if(fSmallAngleApprox !=0) {
1418 // Brookhaven style chamber L/R code
1419 Int_t npaired=0;
1420 for(Int_t ihit1=0;ihit1 < nhits;ihit1++) {
1421 THcDCHit* hit1 = sp->GetHit(ihit1);
1422 Int_t pindex1=hit1->GetPlaneIndex();
1423 if((pindex1%2)==0) { // Odd plane (or even index)
1424 for(Int_t ihit2=0;ihit2<nhits;ihit2++) {
1425 THcDCHit* hit2 = sp->GetHit(ihit2);
1426 if(hit2->GetPlaneIndex()-pindex1 == 1 && TMath::Abs(hit2->GetPos()-hit1->GetPos())<0.51) { // Adjacent plane
1427 if(hit2->GetPos() <= hit1->GetPos() ) {
1428 plusminusknown[ihit1] = -1;
1429 plusminusknown[ihit2] = 1;
1430 } else {
1431 plusminusknown[ihit1] = 1;
1432 plusminusknown[ihit2] = -1;
1433 }
1434 npaired+=2;
1435 }
1436 }
1437 }
1438 }
1439 nplusminus = 1 << (nhits-npaired);
1440 }//end if fSmallAngleApprox!=0
1441 }//end else SOS style
1442 if(nhits < 2) {
1443 if (fdebugstubchisq) cout << "THcDriftChamber::LeftRight: numhits-2 < 0" << endl;
1444 } else if (nhits == 2) {
1445 if (fdebugstubchisq) cout << "THcDriftChamber::LeftRight: numhits-2 = 0" << endl;
1446 }
1447 Int_t nplaneshit = Count1Bits(bitpat);
1448 if (fhdebugflagpr) cout << " num of pm = " << nplusminus << " num of hits =" << nhits << endl;
1449 // Use bit value of integer word to set + or -
1450 // Loop over all combinations of left right.
1451 for(Int_t pmloop=0;pmloop<nplusminus;pmloop++) {
1452 Int_t iswhit = 1;
1453 for(Int_t ihit=0;ihit<nhits;ihit++) {
1454 if(plusminusknown[ihit]!=0) {
1455 plusminus[ihit] = plusminusknown[ihit];
1456 } else {
1457 // Max hits per point has to be less than 32.
1458 if(pmloop & iswhit) {
1459 plusminus[ihit] = 1;
1460 } else {
1461 plusminus[ihit] = -1;
1462 }
1463 iswhit <<= 1;
1464 }
1465 }
1466 if ( (nplaneshit >= fNPlanes-1) || (nplaneshit >= fNPlanes-2 && !fHMSStyleChambers)) {
1467 Double_t chi2;
1468 chi2 = FindStub(nhits, sp,plane_list, bitpat, plusminus, stub);
1469 if (fdebugstubchisq) cout << " pmloop = " << pmloop << " chi2 = " << chi2 << endl;
1470 if(chi2 < minchi2) {
1471 if (fStubMaxXPDiff<100. ) {
1472 // Take best chi2 IF x' of the stub agrees with x' as expected from x.
1473 // Sometimes an incorrect x' gives a good chi2 for the stub, even though it is
1474 // not the correct left/right combination for the real track.
1475 // Rotate x'(=stub(3)) to hut coordinates and compare to x' expected from x.
1476 // THIS ASSUMES STANDARD HMS TUNE!!!!, for which x' is approx. x/875.
1477 if(stub[2]*fTanBeta[plane_list[0]]==-1.0) {
1478 if (fhdebugflagpr) cout << "THcDriftChamber::LeftRight() Error 3" << endl;
1479 }
1480 Double_t xp_fit=stub[2]-fTanBeta[plane_list[0]]
1481 /(1+stub[2]*fTanBeta[plane_list[0]]);
1482 Double_t xp_expect = sp->GetX()*fRatio_xpfp_to_xfp;
1483 if(TMath::Abs(xp_fit-xp_expect)<fStubMaxXPDiff) {
1484 minchi2 = chi2;
1485 for(Int_t ihit=0;ihit<nhits;ihit++) {
1486 plusminusbest[ihit] = plusminus[ihit];
1487 }
1488 sp->SetStub(stub);
1489 } else { // Record best stub failing angle cut
1490 if (chi2 < tmp_minchi2) {
1491 tmp_minchi2 = chi2;
1492 for(Int_t ihit=0;ihit<nhits;ihit++) {
1493 tmp_plusminus[ihit] = plusminus[ihit];
1494 }
1495 for(Int_t i=0;i<4;i++) {
1496 tmp_stub[i] = stub[i];
1497 }
1498 }
1499 }
1500 } else { // Not HMS specific
1501 minchi2 = chi2;
1502 for(Int_t ihit=0;ihit<nhits;ihit++) {
1503 plusminusbest[ihit] = plusminus[ihit];
1504 }
1505 sp->SetStub(stub);
1506 }
1507 }
1509 } else if (nplaneshit >= fNPlanes-2 && fHMSStyleChambers) { // Two planes missing
1510 Double_t chi2 = FindStub(nhits, sp,plane_list, bitpat, plusminus, stub);
1511 //if(debugging)
1512 if (fhdebugflagpr) cout << "pmloop=" << pmloop << " Chi2=" << chi2 << endl;
1513 // Isn't this a bad idea, doing == with reals
1514 if(stub[2]*fTanBeta[plane_list[0]] == -1.0) {
1515 if (fhdebugflagpr) cout << "THcDriftChamber::LeftRight() Error 3" << endl;
1516 }
1517 Double_t xp_fit=stub[2]-fTanBeta[plane_list[0]]
1518 /(1+stub[2]*fTanBeta[plane_list[0]]);
1519 if(TMath::Abs(xp_fit) <= minxp) {
1520 minxp = TMath::Abs(xp_fit);
1521 minchi2 = chi2;
1522 for(Int_t ihit=0;ihit<nhits;ihit++) {
1523 plusminusbest[ihit] = plusminus[ihit];
1524 }
1525 sp->SetStub(stub);
1526 }
1527 } else {
1528 if (fhdebugflagpr) cout << "Insufficient planes hit in THcDriftChamber::LeftRight()" << bitpat <<endl;
1529 }
1530 } // End loop of pm combinations
1531
1532 if (minchi2 == maxchi2 && tmp_minchi2 == maxchi2) {
1533
1534 } else {
1535 if(minchi2 == maxchi2 ) { // No track passed angle cut
1536 minchi2 = tmp_minchi2;
1537 for(Int_t ihit=0;ihit<nhits;ihit++) {
1538 plusminusbest[ihit] = tmp_plusminus[ihit];
1539 }
1540 sp->SetStub(tmp_stub);
1541 }
1542 Double_t *spstub = sp->GetStubP();
1543
1544 // Calculate final coordinate based on plusminusbest
1545 // Update the hit positions in the space points
1546 for(Int_t ihit=0; ihit<nhits; ihit++) {
1547 // Save left/right status with the hit and in the spaleftce point
1548 // In THcDC will decide which to used based on fix_lr flag
1549 sp->GetHit(ihit)->SetLeftRight(plusminusbest[ihit]);
1550 sp->SetHitLR(ihit, plusminusbest[ihit]);
1551 }
1552
1553 // Stubs are calculated in rotated coordinate system
1554 // (I think this rotates in case chambers not perpendicular to central ray)
1555 Int_t pindex=plane_list[0];
1556 if(spstub[2] - fTanBeta[pindex] == -1.0) {
1557 if (fhdebugflagpr) cout << "THcDriftChamber::LeftRight(): stub3 error" << endl;
1558 }
1559 stub[2] = (spstub[2] - fTanBeta[pindex])
1560 / (1.0 + spstub[2]*fTanBeta[pindex]);
1561 if(spstub[2]*fSinBeta[pindex] == -fCosBeta[pindex]) {
1562 if (fhdebugflagpr) cout << "THcDriftChamber::LeftRight(): stub4 error" << endl;
1563 }
1564 stub[3] = spstub[3]
1565 / (spstub[2]*fSinBeta[pindex]+fCosBeta[pindex]);
1566 stub[0] = spstub[0]*fCosBeta[pindex]
1567 - spstub[0]*stub[2]*fSinBeta[pindex];
1568 stub[1] = spstub[1]
1569 - spstub[1]*stub[3]*fSinBeta[pindex];
1570 sp->SetStub(stub);
1571 if (fhdebugflagpr) cout << " Left/Right space pt " << isp+1 << " " << stub[0]<< " " << stub[1] << " " << stub[2]<< " " << stub[3] << endl;
1572 }
1573 }
1574 // Option to print stubs
1575}
1576//_____________________________________________________________________________
1578 Int_t* plane_list, UInt_t bitpat,
1579 Int_t* plusminus, Double_t* stub)
1580{
1581 // For a given combination of L/R, fit a stub to the space point
1582 // This method does a linear least squares fit of a line to the
1583 // hits in an individual chamber. It assumes that the y slope is 0
1584 // The wire coordinate is calculated by
1585 // wire center + plusminus*(drift distance).
1586 // Method is called in a loop over all combinations of plusminus
1587 Double_t zeros[] = {0.0,0.0,0.0};
1588 TVectorD TT; TT.Use(3, zeros); // X, X', Y
1589 Double_t dpos[nhits];
1590 for(Int_t ihit=0;ihit<nhits; ihit++) {
1591 dpos[ihit] = sp->GetHit(ihit)->GetPos()
1592 + plusminus[ihit]*sp->GetHitDist(ihit)
1593 - fPsi0[plane_list[ihit]];
1594 for(Int_t index=0;index<3;index++) {
1595 TT[index]+= dpos[ihit]*fStubCoefs[plane_list[ihit]][index]
1596 /fSigma[plane_list[ihit]];
1597 }
1598 }
1599 // fAA3Inv[bitpat].Print();
1600 // if (fhdebugflagpr) cout << TT[0] << " " << TT[1] << " " << TT[2] << endl;
1601 // TT.Print();
1602
1603 TT *= fAA3Inv[bitpat];
1604 // if (fhdebugflagpr) cout << TT[0] << " " << TT[1] << " " << TT[2] << endl;
1605
1606 // Calculate Chi2. Remember one power of sigma is in fStubCoefs
1607 stub[0] = TT[0];
1608 stub[1] = TT[1];
1609 stub[2] = TT[2];
1610 stub[3] = 0.0;
1611 Double_t chi2=0.0;
1612 for(Int_t ihit=0;ihit<nhits; ihit++) {
1613 chi2 += pow( dpos[ihit]/fSigma[plane_list[ihit]]
1614 - fStubCoefs[plane_list[ihit]][0]*stub[0]
1615 - fStubCoefs[plane_list[ihit]][1]*stub[1]
1616 - fStubCoefs[plane_list[ihit]][2]*stub[2]
1617 , 2);
1618 }
1619 return(chi2);
1620}
1621
1622//_____________________________________________________________________________
1624{
1625 // Destructor. Remove variables from global list.
1626
1627 if (fhdebugflagpr) cout << fChamberNum << ": delete fSpacePoints: " << hex << fSpacePoints << endl;
1628 delete fSpacePoints;
1629 delete fXPlaneClusters;
1630 delete fUPlaneClusters;
1631 delete fVPlaneClusters;
1632 delete fUXPlaneClusters;
1633 delete fVXPlaneClusters;
1634 // Should delete the matrices
1635
1636 if( fIsSetup )
1638 if( fIsInit )
1639 DeleteArrays();
1640 if (fTrackProj) {
1641 fTrackProj->Clear();
1642 delete fTrackProj; fTrackProj = 0;
1643 }
1644}
1645
1646
1647//_____________________________________________________________________________
1649{
1650 // Delete member arrays. Used by destructor.
1651 delete [] fCosBeta; fCosBeta = NULL;
1652 delete [] fSinBeta; fSinBeta = NULL;
1653 delete [] fTanBeta; fTanBeta = NULL;
1654 delete [] fSigma; fSigma = NULL;
1655 delete [] fPsi0; fPsi0 = NULL;
1656 delete [] fStubCoefs; fStubCoefs = NULL;
1657
1658}
1659
1660//_____________________________________________________________________________
1661inline
1663{
1664 // Reset per-event data.
1665 // fNhits = 0;
1666
1667 // fTrackProj->Clear();
1668 fNhits = 0;
1669
1670}
1671
1672//_____________________________________________________________________________
1674{
1675 return(0);
1676}
1677
int Int_t
unsigned int UInt_t
uint32_t time
bool Bool_t
const Bool_t kFALSE
double Double_t
const Bool_t kTRUE
const char Option_t
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void xpos
Option_t Option_t TPoint TPoint const char mode
ClassImp(VDC::AnalyticTTDConv) using namespace std
int ncx
#define MAX_HITS_PER_POINT
#define MAX_SPACE_POINTS
R__EXTERN class THcParmList * gHcParms
Definition THcGlobals.h:11
char * Form(const char *fmt,...)
void Clear(Option_t *option="") override
void Delete(Option_t *option="") override
TObject * ConstructedAt(Int_t idx)
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="")
Int_t RemoveVariables()
THaDetectorBase * GetParent() const
THaApparatus * GetApparatus() const
Drift chamber wire hit info.
Definition THcDCHit.h:16
Int_t GetPlaneNum() const
Definition THcDCHit.h:58
Int_t GetWireNum() const
Definition THcDCHit.h:36
Int_t GetReadoutSide()
Definition THcDCHit.h:63
Int_t GetRawTime() const
Definition THcDCHit.h:37
virtual Double_t ConvertTimeToDist()
Definition THcDCHit.cxx:36
void SetDist(Double_t dist)
Definition THcDCHit.h:51
void SetLeftRight(Int_t lr)
Definition THcDCHit.h:52
void SetTime(Double_t time)
Definition THcDCHit.h:50
Int_t GetNPlaneClust() const
Definition THcDCHit.h:54
Double_t GetPos() const
Definition THcDCHit.h:41
Double_t GetDist() const
Definition THcDCHit.h:40
THcDriftChamberPlane * GetWirePlane() const
Definition THcDCHit.h:46
Int_t GetPlaneIndex() const
Definition THcDCHit.h:59
Double_t GetTime() const
Definition THcDCHit.h:39
void IncrNPlaneClust()
Definition THcDCHit.h:53
Class for clusters in the same orientation planes (V,X or U)
void AddHit(THcDCHit *hit)
THcDCHit * GetHit(Int_t ihit)
void SetXY(Double_t x, Double_t y)
Analyze a package of horizontal drift chambers.
Definition THcDC.h:23
Class for a a single Hall C horizontal drift chamber plane.
void SetPlaneIndex(Int_t index)
Subdetector class for a single drift chamber with several planes.
Int_t fRemove_Sppt_If_One_YPlane
std::map< int, TMatrixD > fAA3Inv
SpacePointHitOutputData fSpHit
virtual Int_t DefineVariables(EMode mode=kDefine)
virtual void ProcessHits(void)
Int_t SpacePointMultiWire(void)
Int_t DestroyPoorSpacePoints(void)
TClonesArray * fTrackProj
TClonesArray * fVPlaneClusters
Int_t FindEasySpacePoint_SOS(Int_t xplane_hitind, Int_t xplanep_hitind)
Int_t FindEasySpacePoint_HMS(Int_t yplane_hitind, Int_t yplanep_hitind)
TClonesArray * fUXPlaneClusters
virtual void PrintDecode(void)
virtual void CorrectHitTimes(void)
std::vector< THcDCHit * > fHits
Double_t ** fStubCoefs
virtual Int_t ApplyCorrections(void)
TClonesArray * fSpacePoints
Int_t fFixPropagationCorrection
std::vector< THcDriftChamberPlane * > fPlanes
Int_t FindHardSpacePoints(void)
virtual Int_t Decode(const THaEvData &)
virtual Int_t NewFindSpacePoints(void)
UInt_t Count1Bits(UInt_t x)
virtual Int_t FindSpacePoints(void)
void Setup(const char *name, const char *description)
virtual void Clear(Option_t *opt="")
Double_t fRatio_xpfp_to_xfp
THaDetectorBase * fParent
TClonesArray * fUPlaneClusters
TClonesArray * fXPlaneClusters
TClonesArray * fVXPlaneClusters
virtual void AddPlane(THcDriftChamberPlane *plane)
virtual void LeftRight(void)
void SelectSpacePoints(void)
virtual Int_t ReadDatabase(const TDatime &date)
Int_t GetNHits() const
Double_t fSpacePointCriterion
Double_t FindStub(Int_t nhits, THcSpacePoint *sp, Int_t *plane_list, UInt_t bitpat, Int_t *plusminus, Double_t *stub)
Int_t LoadParmValues(const DBRequest *list, const char *prefix="")
Retrieve parameter values from the parameter cache.
Class representing a single hit DC.
void SetHitLR(Int_t ihit, Int_t lr)
void SetStub(Double_t stub[4])
void SetCombos(Int_t ncombos)
Double_t GetY()
void SetNHits(Int_t nhits)
void SetXY(Double_t x, Double_t y)
void AddHit(THcDCHit *hit)
THcDCHit * GetHit(Int_t ihit)
Double_t * GetStubP()
void ReplaceHit(Int_t ihit, THcDCHit *hit)
Double_t GetHitDist(Int_t ihit)
void Clear(Option_t *opt="")
Double_t GetX()
Int_t GetCombos()
void SetHitDist(Int_t ihit, Double_t dist)
TMatrixT< Element > & Invert(Double_t *det=nullptr)
const char * GetName() const override
const char * GetTitle() const override
TObject * At(Int_t idx) const override
const TVectorT< Element > & Use(const TVectorT< Element > &v) const
RVec< PromoteType< T > > abs(const RVec< T > &v)
RVec< PromoteTypes< T0, T1 > > pow(const RVec< T0 > &v, const T1 &y)
Double_t y[n]
Double_t x[n]
double dist(AxisAngle const &r1, AxisAngle const &r2)
Double_t Sqrt(Double_t x)
Double_t Cos(Double_t)
Double_t Abs(Double_t d)
Double_t Sin(Double_t)
STL namespace.
std::vector< Int_t > SpHitIndex
std::vector< Int_t > SpNHits