Hall C ROOT/C++ Analyzer (hcana)
Loading...
Searching...
No Matches
THcDC.cxx
Go to the documentation of this file.
1
15#include "THcDC.h"
16#include "THaEvData.h"
17#include "THaDetMap.h"
18#include "THcDetectorMap.h"
19#include "THcGlobals.h"
20#include "THaCutList.h"
21#include "THcParmList.h"
22#include "THcDCTrack.h"
23#include "VarDef.h"
24#include "VarType.h"
25#include "THaTrack.h"
26#include "TClonesArray.h"
27#include "TMath.h"
28#include "TVectorD.h"
29#include "THaApparatus.h"
31#include "Textvars.h" // vsplit
32
33#include <cstring>
34#include <cstdio>
35#include <cstdlib>
36#include <iostream>
37
38using namespace std;
39
40//_____________________________________________________________________________
42 const char* name, const char* description,
43 THaApparatus* apparatus ) :
44 THaTrackingDetector(name,description,apparatus)
45{
46 // Constructor
47
48 fNPlanes = 0; // No planes until we make them
49
50 fXCenter = nullptr;
51 fYCenter = nullptr;
52 fMinHits = nullptr;
53 fMaxHits = nullptr;
54 fMinCombos = nullptr;
55 fSpace_Point_Criterion = nullptr;
56
57 fTdcWinMin = nullptr;
58 fTdcWinMax = nullptr;
59 fCentralTime = nullptr;
60 fNWires = nullptr;
61 fNChamber = nullptr;
62 fWireOrder = nullptr;
63 fDriftTimeSign = nullptr;
64 fReadoutTB = nullptr;
65 fReadoutLR = nullptr;
66
67 fXPos = nullptr;
68 fYPos = nullptr;
69 fZPos = nullptr;
70 fAlphaAngle = nullptr;
71 fBetaAngle = nullptr;
72 fGammaAngle = nullptr;
73 fPitch = nullptr;
74 fCentralWire = nullptr;
75 fPlaneTimeZero = nullptr;
76 fSigma = nullptr;
77
78 // These should be set to zero (in a parameter file) in order to
79 // replicate historical ENGINE behavior
80 fFixLR = 1;
82 fProjectToChamber = 0; // Use 1 for SOS chambers
83
84 fDCTracks = new TClonesArray( "THcDCTrack", 20 );
85
86 fNChamHits = nullptr;
87 fPlaneEvents = nullptr;
88
89 //The version defaults to 0 (old HMS style). 1 is new HMS style and 2 is SHMS style.
90 fVersion = 0;
91}
92
93//_____________________________________________________________________________
94void THcDC::Setup(const char* name, const char* description)
95{
96
97 Bool_t optional = true;
98
99 // Create the chamber and plane objects using parameters.
100 static const char* const here = "Setup";
101
102 THaApparatus *app = GetApparatus();
103 if(app) {
104 // cout << app->GetName() << endl;
105 fPrefix[0]=tolower(app->GetName()[0]);
106 fPrefix[1]='\0';
107 } else {
108 cout << "No apparatus found" << endl;
109 fPrefix[0]='\0';
110 }
111
112 // For now, decide chamber style from the spectrometer name.
113 // Should override with a paramter
114 //cout<<"HMS Style??\t"<<fHMSStyleChambers<<endl;
115 string planenamelist;
116 DBRequest list[]={
117 {"dc_num_planes",&fNPlanes, kInt},
118 {"dc_num_chambers",&fNChambers, kInt},
119 {"dc_tdc_time_per_channel",&fNSperChan, kDouble},
120 {"dc_wire_velocity",&fWireVelocity,kDouble},
121 {"dc_plane_names",&planenamelist, kString},
122 {"dc_version", &fVersion, kInt, 0, optional},
123 {"dc_tdcrefcut", &fTDC_RefTimeCut, kInt, 0, true},
124 {nullptr}
125 };
126
127 fTDC_RefTimeCut = 0; // Minimum allowed reference times
128 gHcParms->LoadParmValues((DBRequest*)&list,fPrefix);
129
130 if(fVersion==0) {
132 } else {
134 }
135
136
137 cout << "Plane Name List: " << planenamelist << endl;
138 cout << "Drift Chambers: " << fNPlanes << " planes in " << fNChambers << " chambers" << endl;
139
140 vector<string> plane_names = Podd::vsplit(planenamelist);
141
142 if(plane_names.size() != (UInt_t) fNPlanes) {
143 cout << "ERROR: Number of planes " << fNPlanes << " doesn't agree with number of plane names " << plane_names.size() << endl;
144 // Should quit. Is there an official way to quit?
145 }
146 fPlaneNames = new char* [fNPlanes];
147 for(Int_t i=0;i<fNPlanes;i++) {
148 fPlaneNames[i] = new char[plane_names[i].length()+1];
149 strcpy(fPlaneNames[i], plane_names[i].c_str());
150 }
151
152 size_t buflen = strlen(description)+100;
153 char *desc = new char[buflen];
154 fPlanes.clear();
155
156 for(Int_t i=0;i<fNPlanes;i++) {
157 strcpy(desc, description);
158 strcat(desc, " Plane ");
159 strcat(desc, fPlaneNames[i]);
160
161 THcDriftChamberPlane* newplane = new THcDriftChamberPlane(fPlaneNames[i], desc, i+1, this);
162 if( !newplane or newplane->IsZombie() ) {
163 Error( Here(here), "Error creating Drift Chamber plane %s. Call expert.", name);
164 MakeZombie();
165 delete [] desc;
166 return;
167 }
168 fPlanes.push_back(newplane);
169 newplane->SetDebug(fDebug);
170 // cout << "Created Drift Chamber Plane " << fPlaneNames[i] << ", " << desc << endl;
171
172 }
173
174 fChambers.clear();
175 for(UInt_t i=0;i<fNChambers;i++) {
176 snprintf(desc,buflen,"Ch%d",i+1);
177
178 // Should construct a better chamber name
179 THcDriftChamber* newchamber = new THcDriftChamber(desc, desc, i+1, this);
180 fChambers.push_back(newchamber);
181 cout << "Created Drift Chamber " << i+1 << ", " << desc << endl;
182 newchamber->SetHMSStyleFlag(fHMSStyleChambers); // Tell the chamber its style
183 }
184 delete [] desc;
185}
186
187//_____________________________________________________________________________
190{
191 // Constructor
192}
193
194//_____________________________________________________________________________
196{
197 // Register the plane objects with the appropriate chambers.
198 // Trigger ReadDatabase to load the remaining parameters
199 Setup(GetName(), GetTitle()); // Create the subdetectors here
200 EffInit();
201
202 char EngineDID[] = "xDC";
203 EngineDID[0] = toupper(GetApparatus()->GetName()[0]);
204 if( gHcDetectorMap->FillMap(fDetMap, EngineDID) < 0 ) {
205 static const char* const here = "Init()";
206 Error( Here(here), "Error filling detectormap for %s.", EngineDID );
207 return kInitError;
208 }
209
210 // Should probably put this in ReadDatabase as we will know the
211 // maximum number of hits after setting up the detector map
212 cout << " DC tdc ref time cut = " << fTDC_RefTimeCut << endl;
213 InitHitList(fDetMap, "THcRawDCHit", fDetMap->GetTotNumChan()+1,
214 fTDC_RefTimeCut, 0);
215
217
218 EStatus status;
219 // This triggers call of ReadDatabase and DefineVariables
220 if( (status = THaTrackingDetector::Init( date )) )
221 return fStatus=status;
222
223 // Initialize planes and add them to chambers
224 for(Int_t ip=0;ip<fNPlanes;ip++) {
225 if((status = fPlanes[ip]->Init( date ))) {
226 return fStatus=status;
227 } else {
228 Int_t chamber=fNChamber[ip];
229 fChambers[chamber-1]->AddPlane(fPlanes[ip]);
230 }
231 }
232 // Initialize chambers
233 for(UInt_t ic=0;ic<fNChambers;ic++) {
234 if((status = fChambers[ic]->Init ( date ))) {
235 return fStatus=status;
236 }
237 }
238 // Retrieve the fiting coefficients
240 for(Int_t ip=0; ip<fNPlanes;ip++) {
241 fPlaneCoeffs[ip] = fPlanes[ip]->GetPlaneCoef();
242 }
243
246 fLR_best = new Double_t [fNPlanes];
251
252
253 // Replace with what we need for Hall C
254 // const DataDest tmp[NDEST] = {
255 // { &fRTNhit, &fRANhit, fRT, fRT_c, fRA, fRA_p, fRA_c, fROff, fRPed, fRGain },
256 // { &fLTNhit, &fLANhit, fLT, fLT_c, fLA, fLA_p, fLA_c, fLOff, fLPed, fLGain }
257 // };
258 // memcpy( fDataDest, tmp, NDEST*sizeof(DataDest) );
259
260 fPresentP = nullptr;
261 THaVar* vpresent = gHaVars->Find(Form("%s.present",GetApparatus()->GetName()));
262 if(vpresent) {
263 fPresentP = (Bool_t *) vpresent->GetValuePointer();
264 }
265 return fStatus = kOK;
266}
267
268//_____________________________________________________________________________
270{
277 // static const char* const here = "ReadDatabase()";
278
279 delete [] fXCenter; fXCenter = new Double_t [fNChambers];
280 delete [] fYCenter; fYCenter = new Double_t [fNChambers];
281 delete [] fMinHits; fMinHits = new Int_t [fNChambers];
282 delete [] fMaxHits; fMaxHits = new Int_t [fNChambers];
283 delete [] fMinCombos; fMinCombos = new Int_t [fNChambers];
285
286 delete [] fTdcWinMin; fTdcWinMin = new Int_t [fNPlanes];
287 delete [] fTdcWinMax; fTdcWinMax = new Int_t [fNPlanes];
288 delete [] fCentralTime; fCentralTime = new Double_t [fNPlanes];
289 delete [] fNWires; fNWires = new Int_t [fNPlanes];
290 delete [] fNChamber; fNChamber = new Int_t [fNPlanes]; // Which chamber is this plane
291 delete [] fWireOrder; fWireOrder = new Int_t [fNPlanes]; // Wire readout order
292 delete [] fDriftTimeSign; fDriftTimeSign = new Int_t [fNPlanes];
293 delete [] fReadoutLR; fReadoutLR = new Int_t [fNPlanes];
294 delete [] fReadoutTB; fReadoutTB = new Int_t [fNPlanes];
295
296 delete [] fXPos; fXPos = new Double_t [fNPlanes];
297 delete [] fYPos; fYPos = new Double_t [fNPlanes];
298 delete [] fZPos; fZPos = new Double_t [fNPlanes];
299 delete [] fAlphaAngle; fAlphaAngle = new Double_t [fNPlanes];
300 delete [] fBetaAngle; fBetaAngle = new Double_t [fNPlanes];
301 delete [] fGammaAngle; fGammaAngle = new Double_t [fNPlanes];
302 delete [] fPitch; fPitch = new Double_t [fNPlanes];
303 delete [] fCentralWire; fCentralWire = new Double_t [fNPlanes];
305 delete [] fSigma; fSigma = new Double_t [fNPlanes];
306
307 Bool_t optional = true;
308
309 DBRequest list[]={
310 {"dc_tdc_time_per_channel",&fNSperChan, kDouble},
311 {"dc_wire_velocity",&fWireVelocity,kDouble},
312
313 {"dc_xcenter", fXCenter, kDouble, fNChambers},
314 {"dc_ycenter", fYCenter, kDouble, fNChambers},
315 {"min_hit", fMinHits, kInt, fNChambers},
316 {"max_pr_hits", fMaxHits, kInt, fNChambers},
317 {"min_combos", fMinCombos, kInt, fNChambers},
318 {"space_point_criterion", fSpace_Point_Criterion, kDouble, fNChambers},
319
320 {"dc_tdc_min_win", fTdcWinMin, kInt, (UInt_t)fNPlanes},
321 {"dc_tdc_max_win", fTdcWinMax, kInt, (UInt_t)fNPlanes},
322 {"dc_central_time", fCentralTime, kDouble, (UInt_t)fNPlanes},
323 {"dc_nrwire", fNWires, kInt, (UInt_t)fNPlanes},
324 {"dc_chamber_planes", fNChamber, kInt, (UInt_t)fNPlanes},
325 {"dc_wire_counting", fWireOrder, kInt, (UInt_t)fNPlanes},
326 {"dc_drifttime_sign", fDriftTimeSign, kInt, (UInt_t)fNPlanes},
327 {"dc_readoutLR", fReadoutLR, kInt, (UInt_t)fNPlanes, optional},
328 {"dc_readoutTB", fReadoutTB, kInt, (UInt_t)fNPlanes, optional},
329
330 {"dc_zpos", fZPos, kDouble, (UInt_t)fNPlanes},
331 {"dc_alpha_angle", fAlphaAngle, kDouble, (UInt_t)fNPlanes},
332 {"dc_beta_angle", fBetaAngle, kDouble, (UInt_t)fNPlanes},
333 {"dc_gamma_angle", fGammaAngle, kDouble, (UInt_t)fNPlanes},
334 {"dc_pitch", fPitch, kDouble, (UInt_t)fNPlanes},
335 {"dc_central_wire", fCentralWire, kDouble, (UInt_t)fNPlanes},
336 {"dc_plane_time_zero", fPlaneTimeZero, kDouble, (UInt_t)fNPlanes},
337 {"dc_sigma", fSigma, kDouble, (UInt_t)fNPlanes},
338 {"single_stub",&fSingleStub, kInt,0,true},
339 {"ntracks_max_fp", &fNTracksMaxFP, kInt,0,true},
340 {"xt_track_criterion", &fXtTrCriterion, kDouble},
341 {"yt_track_criterion", &fYtTrCriterion, kDouble},
342 {"xpt_track_criterion", &fXptTrCriterion, kDouble},
343 {"ypt_track_criterion", &fYptTrCriterion, kDouble},
344 {"dc_fix_lr", &fFixLR, kInt},
345 {"dc_fix_propcorr", &fFixPropagationCorrection, kInt},
346 {"UseNewFindSpacePoints", &fUseNewFindSpacePoints, kInt,0,true},
347 {"UseNewLinkStubs", &fUseNewLinkStubs, kInt,0,true},
348 {"UseNewTrackFit", &fUseNewTrackFit, kInt,0,true},
349 {"debuglinkstubs", &fdebuglinkstubs, kInt},
350 {"debugprintrawdc", &fdebugprintrawdc, kInt},
351 {"debugprintdecodeddc", &fdebugprintdecodeddc, kInt},
352 {"debugflagpr", &fdebugflagpr, kInt},
353 {"debugflagstubs", &fdebugflagstubs, kInt},
354 {"debugtrackprint", &fdebugtrackprint , kInt},
355 {"TrackLargeResidCut", &fTrackLargeResidCut , kDouble,0,true},
356 {nullptr}
357 };
362 fSingleStub=0;
363 for(Int_t ip=0; ip<fNPlanes;ip++) {
364 fReadoutLR[ip] = 0.0;
365 fReadoutTB[ip] = 0.0;
366 }
367
368
369 gHcParms->LoadParmValues((DBRequest*)&list,fPrefix);
370
371 //Set the default plane x,y positions to those of the chamber
372 for(Int_t ip=0; ip<fNPlanes;ip++) {
373 fXPos[ip] = fXCenter[GetNChamber(ip+1)-1];
374 fYPos[ip] = fYCenter[GetNChamber(ip+1)-1];
375 }
376
377 //Load the x,y positions of the planes if they exist (overwrites defaults)
378 DBRequest listOpt[]={
379 {"dc_xpos", fXPos, kDouble, (UInt_t)fNPlanes, optional},
380 {"dc_ypos", fYPos, kDouble, (UInt_t)fNPlanes, optional},
381 {nullptr}
382 };
383 gHcParms->LoadParmValues((DBRequest*)&listOpt,fPrefix);
384 if(fNTracksMaxFP <= 0) fNTracksMaxFP = 10;
385 // if(fNTracksMaxFP > HNRACKS_MAX) fNTracksMaxFP = NHTRACKS_MAX;
386 cout << "Plane counts:";
387 for(Int_t i=0;i<fNPlanes;i++) {
388 cout << " " << fNWires[i];
389 }
390 cout << endl;
391
392 fIsInit = true;
393
394 return kOK;
395}
396
397//_____________________________________________________________________________
399{
403 if( mode == kDefine && fIsSetup ) return kOK;
404 fIsSetup = ( mode == kDefine );
405
406 // Register variables in global list
407
408 RVarDef vars[] = {
409 { "stubtest", "stub test", "fStubTest" },
410 { "nhit", "Number of DC hits", "fNhits" },
411 { "tnhit", "Number of good DC hits", "fNthits" },
412 { "trawhit", "Number of true raw DC hits", "fN_True_RawHits" },
413 { "ntrack", "Number of Tracks", "fNDCTracks" },
414 { "nsp", "Number of Space Points", "fNSp" },
415 { "track_nsp", "Number of spacepoints in track", "fDCTracks.THcDCTrack.GetNSpacePoints()"},
416 { "track_chisq", "Chisq of track", "fDCTracks.THcDCTrack.GetChisq()"},
417 { "track_nhits", "Number of hits in track", "fDCTracks.THcDCTrack.GetNHits()"},
418 { "track_sp1ID", "CH 1 Spacepoint index", "fDCTracks.THcDCTrack.GetSp1_ID()"},
419 { "track_sp2ID", "CH 2 Spacepoint index", "fDCTracks.THcDCTrack.GetSp2_ID()"},
420 { "x", "X at focal plane", "fDCTracks.THcDCTrack.GetX()"},
421 { "y", "Y at focal plane", "fDCTracks.THcDCTrack.GetY()"},
422 { "xp", "XP at focal plane", "fDCTracks.THcDCTrack.GetXP()"},
423 { "yp", "YP at focal plane", "fDCTracks.THcDCTrack.GetYP()"},
424 { "x_fp", "X at focal plane (golden track)", "fX_fp_best"},
425 { "y_fp", "Y at focal plane( golden track)", "fY_fp_best"},
426 { "xp_fp", "XP at focal plane (golden track)", "fXp_fp_best"},
427 { "yp_fp", "YP at focal plane(golden track) ", "fYp_fp_best"},
428 { "chisq", "chisq/dof (golden track) ", "fChisq_best"},
429 { "sp1_id", " (golden track) ", "fSp1_ID_best"},
430 { "sp2_id", " (golden track) ", "fSp2_ID_best"},
431 { "InsideDipoleExit", " ","fInSideDipoleExit_best"},
432 { "gtrack_nsp", " Number of space points in golden track ", "fNsp_best"},
433 { "pos_best", "Pos (golden track)", "fPos_best"},
434 { "dist_best", "Dist (golden track)", "fDist_best"},
435 { "lr_best", "LR sign (golden track)", "fLR_best"},
436 { "residual", "Residuals", "fResiduals"},
437 { "residualExclPlane", "Residuals", "fResidualsExclPlane"},
438 { "wireHitDid","Wire did have matched track hit", "fWire_hit_did"},
439 { "wireHitShould", "Wire should have matched track hit", "fWire_hit_should"},
440 { nullptr }
441 };
442 return DefineVarsFromList( vars, mode );
443
444}
445
446//_____________________________________________________________________________
448{
449 // Destructor. Remove variables from global list.
450
451 if( fIsSetup )
453 if( fIsInit )
454 DeleteArrays();
455
456 // Delete the plane objects
457 for (vector<THcDriftChamberPlane*>::iterator ip = fPlanes.begin();
458 ip != fPlanes.end(); ++ip) delete *ip;
459 // Delete the chamber objects
460 for (vector<THcDriftChamber*>::iterator ip = fChambers.begin();
461 ip != fChambers.end(); ++ip) delete *ip;
462
463 delete fDCTracks;
464}
465
466//_____________________________________________________________________________
468{
469 // Delete member arrays. Used by destructor.
470
471 delete [] fXCenter; fXCenter = nullptr;
472 delete [] fYCenter; fYCenter = nullptr;
473 delete [] fMinHits; fMinHits = nullptr;
474 delete [] fMaxHits; fMaxHits = nullptr;
475 delete [] fMinCombos; fMinCombos = nullptr;
477
478 delete [] fTdcWinMin; fTdcWinMin = nullptr;
479 delete [] fTdcWinMax; fTdcWinMax = nullptr;
480 delete [] fCentralTime; fCentralTime = nullptr;
481 delete [] fNWires; fNWires = nullptr;
482 delete [] fNChamber; fNChamber = nullptr;
483 delete [] fWireOrder; fWireOrder = nullptr;
484 delete [] fDriftTimeSign; fDriftTimeSign = nullptr;
485 delete [] fReadoutLR; fReadoutLR = nullptr;
486 delete [] fReadoutTB; fReadoutTB = nullptr;
487
488 delete [] fXPos; fXPos = nullptr;
489 delete [] fYPos; fYPos = nullptr;
490 delete [] fZPos; fZPos = nullptr;
491 delete [] fAlphaAngle; fAlphaAngle = nullptr;
492 delete [] fBetaAngle; fBetaAngle = nullptr;
493 delete [] fGammaAngle; fGammaAngle = nullptr;
494 delete [] fPitch; fPitch = nullptr;
495 delete [] fCentralWire; fCentralWire = nullptr;
496 delete [] fPlaneTimeZero; fPlaneTimeZero = nullptr;
497 delete [] fSigma; fSigma = nullptr;
498
499 // Efficiency arrays
500 delete [] fNChamHits; fNChamHits = nullptr;
501 delete [] fPlaneEvents; fPlaneEvents = nullptr;
502
503 for( Int_t i = 0; i<fNPlanes; ++i )
504 delete [] fPlaneNames[i];
505 delete [] fPlaneNames;
506
507 delete [] fPlaneCoeffs; fPlaneCoeffs = 0;
508 delete [] fResiduals; fResiduals = 0;
509 delete [] fDist_best; fDist_best = 0;
510 delete [] fPos_best; fPos_best = 0;
511 delete [] fLR_best; fLR_best = 0;
513 delete [] fWire_hit_did; fWire_hit_did = 0;
514 delete [] fWire_hit_should; fWire_hit_should = 0;
515}
516
517//_____________________________________________________________________________
518inline
520{
521 // Reset per-event data.
522 fNDCTracks=0;
523 fStubTest = 0;
524 fNhits = 0;
525 fNthits = 0;
527 fX_fp_best=-10000.;
528 fY_fp_best=-10000.;
529 fXp_fp_best=-10000.;
530 fYp_fp_best=-10000.;
532 fNsp_best=0;
534 for(UInt_t i=0;i<fNChambers;i++) {
535 fChambers[i]->Clear();
536 }
537
538 for(Int_t i=0;i<fNPlanes;i++) {
539 fResiduals[i] = 1000.0;
540 fPos_best[i] = kBig;
541 fDist_best[i] = kBig;
542 fLR_best[i] = 0;
543 fResidualsExclPlane[i] = 1000.0;
544 fWire_hit_did[i] = 1000.0;
545 fWire_hit_should[i] = 1000.0;
546 }
547
548 // fTrackProj->Clear();
549}
550
551//_____________________________________________________________________________
553{
559 ClearEvent();
560 Int_t num_event = evdata.GetEvNum();
561 if (fdebugprintrawdc ||fdebugprintdecodeddc || fdebuglinkstubs || fdebugtrackprint) cout << " event num = " << num_event << endl;
562 // Get the Hall C style hitlist (fRawHitList) for this event
563
564 Bool_t present = kTRUE; // Suppress reference time warnings
565 if(fPresentP) { // if this spectrometer not part of trigger
566 present = *fPresentP;
567 }
568 fNhits = DecodeToHitList(evdata, !present);
569
570
571
572 if(!gHaCuts->Result("Pedestal_event")) {
573 // Let each plane get its hits
574 Int_t nexthit = 0;
575 for(Int_t ip=0;ip<fNPlanes;ip++) {
576 nexthit = fPlanes[ip]->ProcessHits(fRawHitList, nexthit);
577 fN_True_RawHits += fPlanes[ip]->GetNRawhits();
578
579 }
580
581 // fRawHitList is TClones array of THcRawDCHit objects
582 Int_t counter=0;
583 if (fdebugprintrawdc) {
584 cout << " RAW_TOT_HITS = " << fNRawHits << endl;
585 cout << " Hit # " << "Plane " << " Wire " << " Raw TDC " << endl;
586 for(UInt_t ihit = 0; ihit < fNRawHits ; ihit++) {
587 auto* hit = (THcRawDCHit *) fRawHitList->At(ihit);
588 for(UInt_t imhit = 0; imhit < hit->GetRawTdcHit().GetNHits(); imhit++) {
589 counter++;
590 cout << counter << " " << hit->fPlane << " " << hit->fCounter << " " << hit->GetRawTdcHit().GetTimeRaw(imhit) << endl;
591 }
592 }
593 cout << endl;
594 }
595 Eff(); // Accumlate statistics
596 }
597 return fNhits;
598}
599
600//_____________________________________________________________________________
602{
603 return(0);
604}
605
606//_____________________________________________________________________________
608{
615 // Subtract starttimes from each plane hit
616 for(Int_t ip=0;ip<fNPlanes;ip++) {
617 fPlanes[ip]->SubtractStartTime();
618 }
619 //
620 // Let each chamber get its hits
621 for(UInt_t ic=0;ic<fNChambers;ic++) {
622 fChambers[ic]->ProcessHits();
623 fNthits += fChambers[ic]->GetNHits();
624 if (fdebugprintdecodeddc)fChambers[ic]->PrintDecode();
625 }
626 //
627 for(UInt_t i=0;i<fNChambers;i++) {
628 if (fUseNewFindSpacePoints)fChambers[i]->NewFindSpacePoints();
629 if (!fUseNewFindSpacePoints) fChambers[i]->FindSpacePoints();
630 fChambers[i]->CorrectHitTimes();
631 fChambers[i]->LeftRight();
632 }
635 // Now link the stubs between chambers
636 if (fUseNewLinkStubs) {
637 NewLinkStubs();
638 } else {
639 LinkStubs();
640 }
641 if(fNDCTracks > 0) {
643 // Copy tracks into podd tracks list
644 for(UInt_t itrack=0;itrack<fNDCTracks;itrack++) {
645 if (fUseNewTrackFit) NewTrackFit(itrack);
646 THaTrack* theTrack = nullptr;
647 theTrack = AddTrack(tracks, 0.0, 0.0, 0.0, 0.0); // Leaving off trackID
648 // Should we add stubs with AddCluster? Could we do this
649 // by having stubs inherit from cluster
650
651 THcDCTrack *tr = static_cast<THcDCTrack*>( fDCTracks->At(itrack));
652 theTrack->Set(tr->GetX(), tr->GetY(), tr->GetXP(), tr->GetYP());
653 theTrack->SetFlag((UInt_t) 0);
654 // Need to look at how engine does chi2 and track selection. Reduced?
655 theTrack->SetChi2(tr->GetChisq(),tr->GetNFree());
656 // CalcFocalPlaneCoords. Aren't our tracks already in focal plane coords
657 // We should have some kind of track ID so that the THaTrack can be
658 // associate back with the DC track
659 // Assign the track number
660 theTrack->SetTrkNum(itrack+1);
661 }
663 }
664 //
665 if(fNDCTracks > 0) {
666 THcHallCSpectrometer *spectro = dynamic_cast<THcHallCSpectrometer*>(GetApparatus());
667 for (Int_t it=0;it<tracks.GetLast()+1;it++) {
668 THaTrack* track = static_cast<THaTrack*>( tracks[it] );
669 Double_t xptar=kBig,yptar=kBig,ytar=kBig,delta=kBig;
670 Double_t xtar=0;
671 spectro->CalculateTargetQuantities(track,xtar,xptar,ytar,yptar,delta);
672 // Transfer results to track
673 // No beam raster yet
674 //; In transport coordinates phi = hyptar = dy/dz and theta = hxptar = dx/dz
675 //; but for unknown reasons the yp offset is named htheta_offset
676 //; and the xp offset is named hphi_offset
677
678 track->SetTarget(0.0, ytar*100.0, xptar, yptar);
679 track->SetDp(delta*100.0); // Percent.
680 Double_t ptemp = spectro->GetPcentral()*(1+track->GetDp()/100.0);
681 track->SetMomentum(ptemp);
682 TVector3 pvect_temp;
683 spectro->TransportToLab(track->GetP(),track->GetTTheta(),track->GetTPhi(),pvect_temp);
684 track->SetPvect(pvect_temp);
685 }
686 }
687 //
688
690
691 return 0;
692}
693
694//_____________________________________________________________________________
696{
697 return 0;
698}
699//
700void THcDC::SetFocalPlaneBestTrack(Int_t golden_track_index)
701{
702 THcDCTrack *tr1 = static_cast<THcDCTrack*>( fDCTracks->At(golden_track_index));
703 fX_fp_best=tr1->GetX();
704 fY_fp_best=tr1->GetY();
705 fXp_fp_best=tr1->GetXP();
706 fYp_fp_best=tr1->GetYP();
708 fInSideDipoleExit_best = app->InsideDipoleExitWindow(fX_fp_best, fXp_fp_best ,fY_fp_best,fYp_fp_best);
709 fSp1_ID_best=tr1->GetSp1_ID();
710 fSp2_ID_best=tr1->GetSp2_ID();
711 fChisq_best=tr1->GetChisq();
713 for (UInt_t ihit = 0; ihit < UInt_t (tr1->GetNHits()); ihit++) {
714 THcDCHit *hit = tr1->GetHit(ihit);
715 Int_t plane = hit->GetPlaneNum() - 1;
716 fResiduals[plane] = tr1->GetResidual(plane);
717 fResidualsExclPlane[plane] = tr1->GetResidualExclPlane(plane);
718 fDist_best[plane] = tr1->GetHitDist(ihit);
719 fLR_best[plane] = tr1->GetHitLR(ihit);
720 fPos_best[plane] = hit->GetPos();
721 }
722 EfficiencyPerWire(golden_track_index);
723}
724//
725void THcDC::EfficiencyPerWire(Int_t golden_track_index)
726{
727 THcDCTrack *tr1 = static_cast<THcDCTrack*>( fDCTracks->At(golden_track_index));
728 Double_t track_pos;
729 for (UInt_t ihit = 0; ihit < UInt_t (tr1->GetNHits()); ihit++) {
730 THcDCHit *hit = tr1->GetHit(ihit);
731 Int_t plane = hit->GetPlaneNum() - 1;
732 track_pos=tr1->GetCoord(plane);
733 Int_t wire_num = hit->GetWireNum();
734 Int_t wire_track_num=round(fPlanes[plane]->CalcWireFromPos(track_pos));
735 if ( (wire_num-wire_track_num) ==0) fWire_hit_did[plane]=wire_num;
736 }
737 for(Int_t ip=0; ip<fNPlanes;ip++) {
738 track_pos=tr1->GetCoord(ip);
739 Int_t wire_should = round(fPlanes[ip]->CalcWireFromPos(track_pos));
740 fWire_hit_should[ip]=wire_should;
741 }
742}
743//
745{
746 for(UInt_t ich=0;ich<fNChambers;ich++) {
747 printf("%s %2d %s %3d %s %3d \n"," chamber = ",fChambers[ich]->GetChamberNum()," number of hits = ",fChambers[ich]->GetNHits()," number of spacepoints = ",fChambers[ich]->GetNSpacePoints());
748 printf("%6s %-8s %-8s %6s %6s \n"," "," "," ","Number","Number");
749 printf("%6s %-8s %-8s %6s %6s \n","Point","x","y"," hits ","combos");
750 TClonesArray* spacepointarray = fChambers[ich]->GetSpacePointsP();
751 for(Int_t isp=0;isp<fChambers[ich]->GetNSpacePoints();isp++) {
752 THcSpacePoint* sp = (THcSpacePoint*)(spacepointarray->At(isp));
753 printf("%5d %8.5f %8.5f %5d %5d \n",isp+1,sp->GetX(),sp->GetY(),sp->GetNHits(),sp->GetCombos()) ;
754 for (Int_t ii=0;ii<sp->GetNHits();ii++) {
755 THcDCHit* hittemp = (THcDCHit*)(sp->GetHit(ii));
756 printf("%3d %3d %3d \n",hittemp->GetPlaneNum(),hittemp->GetWireNum(),hittemp->GetLR());
757 }
758 }
759 }
760}
761//
762//
764{
765 for(UInt_t ich=0;ich<fNChambers;ich++) {
766 printf("%s %3d \n"," Stub fit results Chamber = ",ich+1);
767 printf("%-5s %-18s %-18s %-18s %-18s\n","point","x_t","y_t","xp_t","yp_t");
768 printf("%-5s %-18s %-18s %-18s %-18s\n"," ","[cm]","[cm]","[cm]","[cm]");
769 TClonesArray* spacepointarray = fChambers[ich]->GetSpacePointsP();
770 for(Int_t isp=0;isp<fChambers[ich]->GetNSpacePoints();isp++) {
771 THcSpacePoint* sp = (THcSpacePoint*)(spacepointarray->At(isp));
772 Double_t *spstubt=sp->GetStubP();
773 printf("%-5d % 15.10e % 15.10e % 15.10e % 15.10e \n",isp+1,spstubt[0],spstubt[1],spstubt[2],spstubt[3]);
774 }
775 }
776}
777//
779{
780 fNDCTracks=0; // Number of Focal Plane tracks found
781 fDCTracks->Delete();
782 if ( fChambers[0]->GetNSpacePoints() >0 && fChambers[1]->GetNSpacePoints() >0 ) {
783 std::vector<THcSpacePoint*> fSp; // vector of spacepoints from ch1 and ch2
784 fNSp=0;
785 fSp.clear();
786 // create array of spacepoints
787 for(UInt_t ich=0;ich<fNChambers;ich++) {
788 Int_t nchamber=fChambers[ich]->GetChamberNum();
789 TClonesArray* spacepointarray = fChambers[ich]->GetSpacePointsP();
790 for(Int_t isp=0;isp<fChambers[ich]->GetNSpacePoints();isp++) {
791 fSp.push_back(static_cast<THcSpacePoint*>(spacepointarray->At(isp)));
792 fSp[fNSp]->fNChamber = nchamber;
793 fSp[fNSp]->fNChamber_spnum = isp;
794 fNSp++;
795 }
796 }
797 //
798 for(Int_t isp1=0;isp1<fChambers[0]->GetNSpacePoints();isp1++) {
799 THcSpacePoint* sp1 = fSp[isp1];
800 Int_t ch1_nsp = fChambers[0]->GetNSpacePoints();
801 Int_t ch2_nsp = fChambers[1]->GetNSpacePoints();
802 for(Int_t isp2=ch1_nsp;isp2<ch1_nsp+ch2_nsp;isp2++) {
803 THcSpacePoint* sp2=fSp[isp2];
804 if (sp1->GetSetStubFlag()&&sp2->GetSetStubFlag()) {
805 Double_t *spstub1=sp1->GetStubP();
806 Double_t *spstub2=sp2->GetStubP();
807 Double_t dposx = spstub1[0] - spstub2[0];
808 Double_t dposy;
810 Double_t y1=spstub1[1]+fChambers[sp1->fNChamber]->GetZPos()*spstub1[3];
811 Double_t y2=spstub2[1]+fChambers[sp2->fNChamber]->GetZPos()*spstub2[3];
812 dposy = y1-y2;
813 } else {
814 dposy = spstub1[1] - spstub2[1];
815 }
816 Double_t dposxp = spstub1[2] - spstub2[2];
817 Double_t dposyp = spstub1[3] - spstub2[3];
818 if((TMath::Abs(dposx) < fXtTrCriterion)
819 && (TMath::Abs(dposy) < fYtTrCriterion)
820 && (TMath::Abs(dposxp) < fXptTrCriterion)
821 && (TMath::Abs(dposyp) < fYptTrCriterion)) {
822 fStubTest = 1;
823 if(fNDCTracks < MAXTRACKS) {
824 THcDCTrack *theDCTrack = new( (*fDCTracks)[fNDCTracks++]) THcDCTrack(fNPlanes);
825 theDCTrack->AddSpacePoint(sp1);
826 theDCTrack->AddSpacePoint(sp2);
827 theDCTrack->SetSp1_ID(sp1->fNChamber_spnum);
828 theDCTrack->SetSp2_ID(sp2->fNChamber_spnum);
829 }
830 }
831 }
832 }// isp2
833 } // isp1
834 //
835 if (fNDCTracks == 0) { // make tracks from all combinations of spacepoints
836 for(Int_t isp1=0;isp1<fChambers[0]->GetNSpacePoints();isp1++) {
837 THcSpacePoint* sp1 = fSp[isp1];
838 Int_t ch1_nsp = fChambers[0]->GetNSpacePoints();
839 Int_t ch2_nsp = fChambers[1]->GetNSpacePoints();
840 for(Int_t isp2=ch1_nsp;isp2<ch1_nsp+ch2_nsp;isp2++) {
841 THcSpacePoint* sp2=fSp[isp2];
842 fStubTest = 1;
843 if(fNDCTracks < MAXTRACKS) {
844 THcDCTrack *theDCTrack = new( (*fDCTracks)[fNDCTracks++]) THcDCTrack(fNPlanes);
845 theDCTrack->AddSpacePoint(sp1);
846 theDCTrack->AddSpacePoint(sp2);
847 theDCTrack->SetSp1_ID(sp1->fNChamber_spnum);
848 theDCTrack->SetSp2_ID(sp2->fNChamber_spnum);
849 }
850 }
851 }
852 }
853 //
854 } // both chambers have spacepoints
855//
856
857}
858//_____________________________________________________________________________
860{
876 std::vector<THcSpacePoint*> fSp;
877 fNSp=0;
878 fSp.clear();
879 fSp.reserve(100);
880 fNDCTracks=0; // Number of Focal Plane tracks found
881 fDCTracks->Delete();
882 // Make a vector of pointers to the SpacePoints
883 //if (fChambers[0]->GetNSpacePoints()+fChambers[1]->GetNSpacePoints()>10) return;
884
885 for(UInt_t ich=0;ich<fNChambers;ich++) {
886 Int_t nchamber=fChambers[ich]->GetChamberNum();
887 TClonesArray* spacepointarray = fChambers[ich]->GetSpacePointsP();
888 for(Int_t isp=0;isp<fChambers[ich]->GetNSpacePoints();isp++) {
889 fSp.push_back(static_cast<THcSpacePoint*>(spacepointarray->At(isp)));
890 fSp[fNSp]->fNChamber = nchamber;
891 fSp[fNSp]->fNChamber_spnum = isp;
892 fNSp++;
893 if (ich==0 && fNSp>50) break;
894 if (fNSp>100) break;
895 }
896 }
897 Double_t stubminx = 999999;
898 Double_t stubminy = 999999;
899 Double_t stubminxp = 999999;
900 Double_t stubminyp = 999999;
901 Int_t stub_tracks[MAXTRACKS];
902 if(fSingleStub==0) {
903 for(Int_t isp1=0;isp1<fNSp-1;isp1++) { // isp1 is index/id in total list of space points
904 THcSpacePoint* sp1 = fSp[isp1];
905 Int_t sptracks=0;
906 // Now make sure this sp is not already used in a sp.
907 // Could this be done by having a sp point to the track it is in?
908 Int_t tryflag=1;
909 for(UInt_t itrack=0;itrack<fNDCTracks;itrack++) {
910 THcDCTrack *theDCTrack = static_cast<THcDCTrack*>( fDCTracks->At(itrack));
911 for(Int_t isp=0;isp<theDCTrack->GetNSpacePoints();isp++) {
912 // isp is index into list of space points attached to theDCTrack
913 if(theDCTrack->GetSpacePoint(isp) == sp1) {
914 tryflag=0;
915 }
916 }
917 }
918 if(tryflag) { // SP not already part of a track
919 Int_t newtrack=1;
920 for(Int_t isp2=isp1+1;isp2<fNSp;isp2++) {
921 THcSpacePoint* sp2=fSp[isp2];
922 if(sp1->fNChamber!=sp2->fNChamber&&sp1->GetSetStubFlag()&&sp2->GetSetStubFlag()) {
923 Double_t *spstub1=sp1->GetStubP();
924 Double_t *spstub2=sp2->GetStubP();
925 Double_t dposx = spstub1[0] - spstub2[0];
926 Double_t dposy;
927 if(fProjectToChamber) { // From SOS s_link_stubs
928 // Since single chamber resolution is ~50mr, and the maximum y`
929 // angle is about 30mr, use differenece between y AT CHAMBERS, rather
930 // than at focal plane. (Project back to chamber, to take out y' uncertainty)
931 // (Should this be done for SHMS and HMS too?)
932 Double_t y1=spstub1[1]+fChambers[sp1->fNChamber]->GetZPos()*spstub1[3];
933 Double_t y2=spstub2[1]+fChambers[sp2->fNChamber]->GetZPos()*spstub2[3];
934 dposy = y1-y2;
935 } else {
936 dposy = spstub1[1] - spstub2[1];
937 }
938 Double_t dposxp = spstub1[2] - spstub2[2];
939 Double_t dposyp = spstub1[3] - spstub2[3];
940
941 // What is the point of saving these stubmin values. They
942 // Don't seem to be used anywhere except that they can be
943 // printed out if hbypass_track_eff_files is zero.
944 if(TMath::Abs(dposx)<TMath::Abs(stubminx)) stubminx = dposx;
945 if(TMath::Abs(dposy)<TMath::Abs(stubminy)) stubminy = dposy;
946 if(TMath::Abs(dposxp)<TMath::Abs(stubminxp)) stubminxp = dposxp;
947 if(TMath::Abs(dposyp)<TMath::Abs(stubminyp)) stubminyp = dposyp;
948
949 // if hbypass_track_eff_files == 0 then
950 // Print out each stubminX that is less that its criterion
951
952 if((TMath::Abs(dposx) < fXtTrCriterion)
953 && (TMath::Abs(dposy) < fYtTrCriterion)
954 && (TMath::Abs(dposxp) < fXptTrCriterion)
955 && (TMath::Abs(dposyp) < fYptTrCriterion)) {
956 if(newtrack) {
957 assert(sptracks==0);
958 fStubTest = 1;
959 //stubtest=1; Used in h_track_tests.f
960 // Make a new track if there are not to many
961 if(fNDCTracks < MAXTRACKS) {
962 sptracks=0; // Number of tracks with this seed
963 stub_tracks[sptracks++] = fNDCTracks;
964 THcDCTrack *theDCTrack = new( (*fDCTracks)[fNDCTracks++]) THcDCTrack(fNPlanes);
965 theDCTrack->AddSpacePoint(sp1);
966 theDCTrack->AddSpacePoint(sp2);
967 if (sp1->fNChamber==1) theDCTrack->SetSp1_ID(sp1->fNChamber_spnum);
968 if (sp1->fNChamber==2) theDCTrack->SetSp2_ID(sp1->fNChamber_spnum);
969 if (sp2->fNChamber==1) theDCTrack->SetSp1_ID(sp2->fNChamber_spnum);
970 if (sp2->fNChamber==2) theDCTrack->SetSp2_ID(sp2->fNChamber_spnum);
971 newtrack = 0; // Make no more tracks in this loop
972 // (But could replace a SP?)
973 } else {
974 if (fHMSStyleChambers) {
975 fNDCTracks=0;
976 return;
977 }
978 }
979 } else {
980 // Check if there is another space point in the same chamber
981 for(Int_t itrack=0;itrack<sptracks;itrack++) {
982 Int_t track=stub_tracks[itrack];
983 THcDCTrack *theDCTrack = static_cast<THcDCTrack*>( fDCTracks->At(track));
984 Int_t spoint=-1;
985 Int_t duppoint=0;
986 for(Int_t isp=0;isp<theDCTrack->GetNSpacePoints();isp++) {
987 // isp is index of space points in theDCTrack
988 if(sp2->fNChamber ==
989 theDCTrack->GetSpacePoint(isp)->fNChamber) {
990 spoint=isp;
991 }
992 if(sp2==theDCTrack->GetSpacePoint(isp)) {
993 duppoint=1;
994 }
995 } // End loop over sp in tracks with isp1
996 // If there is no other space point in this chamber
997 // add this space point to current track(2)
998 if(!duppoint) {
999 if(spoint<0) {
1000 theDCTrack->AddSpacePoint(sp2);
1001 if (sp2->fNChamber==1) theDCTrack->SetSp1_ID(sp2->fNChamber_spnum);
1002 if (sp2->fNChamber==2) theDCTrack->SetSp2_ID(sp2->fNChamber_spnum);
1003 } else {
1004 // If there is another point in the same chamber
1005 // in this track create a new track with all the
1006 // same space points except spoint
1007 if(fNDCTracks < MAXTRACKS) {
1008 stub_tracks[sptracks++] = fNDCTracks;
1009 THcDCTrack *newDCTrack = new( (*fDCTracks)[fNDCTracks++]) THcDCTrack(fNPlanes);
1010 for(Int_t isp=0;isp<theDCTrack->GetNSpacePoints();isp++) {
1011 if(isp!=spoint) {
1012 newDCTrack->AddSpacePoint(theDCTrack->GetSpacePoint(isp));
1013 if (theDCTrack->GetSpacePoint(isp)->fNChamber==1) newDCTrack->SetSp1_ID(theDCTrack->GetSpacePoint(isp)->fNChamber_spnum);
1014 if (theDCTrack->GetSpacePoint(isp)->fNChamber==2) newDCTrack->SetSp2_ID(theDCTrack->GetSpacePoint(isp)->fNChamber_spnum);
1015 } else {
1016 newDCTrack->AddSpacePoint(sp2);
1017 if (sp2->fNChamber==1) newDCTrack->SetSp1_ID(sp2->fNChamber_spnum);
1018 if (sp2->fNChamber==2) newDCTrack->SetSp2_ID(sp2->fNChamber_spnum);
1019 } // End check for dup on copy
1020 } // End copy of track
1021 } else {
1022 if (fHMSStyleChambers) {
1023 if (fdebuglinkstubs) cout << "EPIC FAIL 2: Too many tracks found in THcDC::LinkStubs maxtracks = " << MAXTRACKS << endl;
1024 fNDCTracks=0;
1025 return; // Max # of allowed tracks
1026 }
1027 }
1028 } // end if on same chamber
1029 } // end if on duplicate point
1030 } // end for over tracks with isp1
1031 } // else newtrack
1032 } // criterion
1033 } // end test on same chamber
1034 } // end isp2 loop over new space points
1035 } // end test on tryflag
1036 } // end isp1 outer loop over space points
1037 //
1038 //
1039 } else { // Make track out of each single space point
1040 for(Int_t isp=0;isp<fNSp;isp++) {
1041 if(fNDCTracks<MAXTRACKS) {
1042 // Need some constructed t thingy
1043 if (fSp[isp]->GetSetStubFlag()) {
1044 THcDCTrack *newDCTrack = new( (*fDCTracks)[fNDCTracks++]) THcDCTrack(fNPlanes);
1045 newDCTrack->AddSpacePoint(fSp[isp]);
1046 }
1047 } else {
1048 if (fdebuglinkstubs) cout << "EPIC FAIL 3: Too many tracks found in THcDC::LinkStubs" << endl;
1049 fNDCTracks=0;
1050 // Do something here to fail this event
1051 return; // Max # of allowed tracks
1052 }
1053 }
1054 }
1056 if (fdebuglinkstubs) {
1057 cout << " Number of tracks from link stubs = " << fNDCTracks << endl;
1058 printf("%s %s \n","Track","Plane Wire ");
1059 for (UInt_t itrack=0;itrack<fNDCTracks;itrack++) {
1060 THcDCTrack *tempTrack = (THcDCTrack*)( fDCTracks->At(itrack));
1061 printf("%-5d ",itrack+1);
1062 for (Int_t ihit=0;ihit<tempTrack->GetNHits();ihit++) {
1063 THcDCHit* hit=(THcDCHit*)(tempTrack->GetHit(ihit));
1064 printf("%3d %3d",hit->GetPlaneNum(),hit->GetWireNum());
1065 }
1066 printf("\n");
1067 }
1068 }
1069}
1070//_____________________________________________________________________________
1071void THcDC::FitLineToTrack(Int_t TrackHits,Double_t coords[],Int_t planes[],Double_t wiresigma[], Double_t TrackCoord[], Double_t save_ray[])
1072{
1073 const Int_t raycoeffmap[]={4,5,2,3};
1074 TVectorD TT(4);
1075 TMatrixD AA(4,4);
1076 for(Int_t irayp=0;irayp<NUM_FPRAY;irayp++) {
1077 TT[irayp] = 0.0;
1078 for(Int_t ihit=0;ihit < TrackHits;ihit++) {
1079 TT[irayp] += (coords[ihit]*fPlaneCoeffs[planes[ihit]][raycoeffmap[irayp]])/pow(wiresigma[ihit],2);
1080 }
1081 }
1082 //
1083 for(Int_t irayp=0;irayp<NUM_FPRAY;irayp++) {
1084 for(Int_t jrayp=0;jrayp<NUM_FPRAY;jrayp++) {
1085 AA[irayp][jrayp] = 0.0;
1086 if(jrayp<irayp) { // Symmetric
1087 AA[irayp][jrayp] = AA[jrayp][irayp];
1088 } else {
1089 for(Int_t ihit=0;ihit <TrackHits ;ihit++) {
1090 AA[irayp][jrayp] += fPlaneCoeffs[planes[ihit]][raycoeffmap[irayp]]*
1091 fPlaneCoeffs[planes[ihit]][raycoeffmap[jrayp]]/pow(wiresigma[ihit],2);
1092 } //end ihit loop
1093 }
1094 }
1095 }
1096
1097 // Solve 4x4 equations
1098 TVectorD dray(NUM_FPRAY);
1099 AA.Invert();
1100 dray = AA*TT;
1101 for(Int_t ir=0;ir<NUM_FPRAY;ir++) save_ray[ir]=dray[ir];
1102 for(Int_t iplane=0;iplane < fNPlanes; iplane++) {
1103 Double_t sumcoord=0.0;
1104 for(Int_t ir=0;ir<NUM_FPRAY;ir++) {
1105 sumcoord += fPlaneCoeffs[iplane][raycoeffmap[ir]]*dray[ir];
1106 }
1107 TrackCoord[iplane]=sumcoord;
1108 }
1109}
1110//_____________________________________________________________________________
1112{
1113 THcDCTrack *tr = static_cast<THcDCTrack*>( fDCTracks->At(TrackIndex));
1114 Double_t minchi2=100000.;
1115 Int_t MAXHITS=12;
1116 Int_t plusminus[MAXHITS];
1117 Int_t plusminusknown[MAXHITS];
1118 Int_t TrackHits = tr->GetNHits();
1119 Double_t coords[TrackHits];
1120 Double_t TrackCoord[TrackHits];
1121 Int_t planes[TrackHits];
1122 Double_t wiresigma[TrackHits];
1123 Double_t save_ray[NUM_FPRAY];
1124 for(Int_t ihit=0;ihit<MAXHITS;ihit++) {
1125 plusminusknown[ihit]=0;
1126 }
1127 Int_t nplusminus=1<<TrackHits;
1128 for(Int_t pmloop=0;pmloop<nplusminus;pmloop++) {
1129 Int_t iswhit = 1;
1130 for(Int_t ihit=0;ihit<TrackHits;ihit++) {
1131 if(plusminusknown[ihit]!=0) {
1132 plusminus[ihit] = plusminusknown[ihit];
1133 } else {
1134 if(pmloop & iswhit) {
1135 plusminus[ihit] = 1;
1136 } else {
1137 plusminus[ihit] = -1;
1138 }
1139 iswhit <<= 1;
1140 }
1141 }
1142 for(Int_t ihit=0;ihit<TrackHits;ihit++) {
1143 THcDCHit* hit=tr->GetHit(ihit);
1144 coords[ihit] =hit->GetPos()+ plusminus[ihit]*tr->GetHitDist(ihit);
1145 planes[ihit] =hit->GetPlaneNum()-1;
1146 wiresigma[ihit] = hit->GetWireSigma();
1147 }
1148 FitLineToTrack(TrackHits,coords,planes,wiresigma,TrackCoord,save_ray);
1149 Double_t chi2 = 0.0;
1150 for(Int_t ihit=0;ihit < TrackHits;ihit++) {
1151 Double_t residual = coords[ihit] - TrackCoord[planes[ihit]];
1152 chi2 += pow(residual/wiresigma[ihit],2);
1153 }
1154 if (chi2 < minchi2) {
1155 minchi2 = chi2;
1156 tr->SetVector(save_ray[0], save_ray[1], 0.0, save_ray[2], save_ray[3]);
1157 for(Int_t ihit=0;ihit < TrackHits;ihit++) {
1158 tr->SetHitLR(ihit, plusminus[ihit]);
1159 tr->SetCoord(planes[ihit], TrackCoord[planes[ihit]]);
1160 Double_t residual = coords[ihit] -TrackCoord[planes[ihit]] ;
1161 tr->SetResidual(planes[ihit], residual);
1162 tr->SetChisq(chi2);
1163 tr->SetNFree(TrackHits - 4);
1164 }
1165 //cout << " fit = " << chi2 << " " << dray[0] << " " << dray[1] << " "<< dray[2] << " "<< dray[3] << " " << endl ;
1166 }
1167 } // pmloop
1168 //
1169 // calculate ray without a plane in track
1170 const Int_t raycoeffmap[]={4,5,2,3};
1171 for(Int_t ihit=0;ihit<TrackHits;ihit++) {
1172 THcDCHit* hit=tr->GetHit(ihit);
1173 coords[ihit] =hit->GetPos()+ tr->GetHitLR(ihit)*tr->GetHitDist(ihit);
1174 planes[ihit] =hit->GetPlaneNum()-1;
1175 wiresigma[ihit] = hit->GetWireSigma();
1176 }
1177 for(Int_t ipl_hit=0;ipl_hit < tr->GetNHits();ipl_hit++) {
1178
1179
1180 if(tr->GetNFree() > 0) {
1181 TVectorD TT(NUM_FPRAY);
1183 for(Int_t irayp=0;irayp<NUM_FPRAY;irayp++) {
1184 TT[irayp] = 0.0;
1185 for(Int_t ihit=0;ihit < tr->GetNHits();ihit++) {
1186
1187
1188 THcDCHit* hit=tr->GetHit(ihit);
1189
1190 if (ihit != ipl_hit) {
1191 TT[irayp] += (coords[ihit]*
1192 fPlaneCoeffs[planes[ihit]][raycoeffmap[irayp]])
1193 /pow(hit->GetWireSigma(),2);
1194
1195 }
1196 }
1197 }
1198 for(Int_t irayp=0;irayp<NUM_FPRAY;irayp++) {
1199 for(Int_t jrayp=0;jrayp<NUM_FPRAY;jrayp++) {
1200 AA[irayp][jrayp] = 0.0;
1201 if(jrayp<irayp) { // Symmetric
1202 AA[irayp][jrayp] = AA[jrayp][irayp];
1203 } else {
1204
1205 for(Int_t ihit=0;ihit < tr->GetNHits();ihit++) {
1206
1207 THcDCHit* hit=tr->GetHit(ihit);
1208
1209
1210 if (ihit != ipl_hit) {
1211 AA[irayp][jrayp] += fPlaneCoeffs[planes[ihit]][raycoeffmap[irayp]]*
1212 fPlaneCoeffs[planes[ihit]][raycoeffmap[jrayp]]/
1213 pow(hit->GetWireSigma(),2);
1214
1215 }
1216 }
1217 }
1218 }
1219 }
1220 //
1221 // Solve 4x4 equations
1222 // Should check that it is invertable
1223 TVectorD dray(NUM_FPRAY);
1224 AA.Invert();
1225 dray = AA*TT;
1226 Double_t coord=0.0;
1227 for(Int_t ir=0;ir<NUM_FPRAY;ir++) {
1228 coord += fPlaneCoeffs[planes[ipl_hit]][raycoeffmap[ir]]*dray[ir];
1229 // cout << "ir = " << ir << ", dray[ir] = " << dray[ir] << endl;
1230 }
1231 Double_t residual = coords[ipl_hit] - coord;
1232 tr->SetResidualExclPlane(planes[ipl_hit], residual);
1233 }//Calculate residual without plane
1234 //
1235 }
1236 //
1237
1238}
1239
1240//_____________________________________________________________________________
1242{
1247 // Number of ray parameters in focal plane.
1248 const Int_t raycoeffmap[]={4,5,2,3};
1249
1250 Double_t dummychi2 = 1.0E4;
1251 for(UInt_t itrack=0;itrack<fNDCTracks;itrack++) {
1252 // Double_t chi2 = dummychi2;
1253 // Int_t htrack_fit_num = itrack;
1254 Bool_t CheckResid = kTRUE;
1255 Int_t NLoops=0;
1256 Int_t NLoopsMax=1;
1257 while (CheckResid) {
1258 Int_t MaxResidHitIndex = -1;
1259 Double_t MaxResid = -1000;
1260 THcDCTrack *theDCTrack = static_cast<THcDCTrack*>( fDCTracks->At(itrack));
1261
1262 Double_t coords[theDCTrack->GetNHits()];
1263 Int_t planes[theDCTrack->GetNHits()];
1264 for(Int_t ihit=0;ihit < theDCTrack->GetNHits();ihit++) {
1265 THcDCHit* hit=theDCTrack->GetHit(ihit);
1266 planes[ihit]=hit->GetPlaneNum()-1;
1267 if(fFixLR) {
1269 coords[ihit] = hit->GetPos()
1270 + theDCTrack->GetHitLR(ihit)*theDCTrack->GetHitDist(ihit);
1271 } else {
1272 coords[ihit] = hit->GetPos()
1273 + theDCTrack->GetHitLR(ihit)*hit->GetDist();
1274 }
1275
1276 } else {
1278 coords[ihit] = hit->GetPos()
1279 + hit->GetLR()*theDCTrack->GetHitDist(ihit);
1280 } else {
1281 coords[ihit] = hit->GetCoord();
1282 }
1283 }
1284
1285
1286 } //end loop over hits
1287
1288 theDCTrack->SetNFree(theDCTrack->GetNHits() - NUM_FPRAY);
1289 Double_t chi2 = dummychi2;
1290 if(theDCTrack->GetNFree() > 0) {
1291 TVectorD TT(NUM_FPRAY);
1293 for(Int_t irayp=0;irayp<NUM_FPRAY;irayp++) {
1294 TT[irayp] = 0.0;
1295 for(Int_t ihit=0;ihit < theDCTrack->GetNHits();ihit++) {
1296
1297 THcDCHit* hit=theDCTrack->GetHit(ihit);
1298 TT[irayp] += (coords[ihit]*fPlaneCoeffs[planes[ihit]][raycoeffmap[irayp]])/pow(hit->GetWireSigma(),2);
1299 // if (hit->GetPlaneNum()==5)
1300 // {
1301 // // cout << "Plane: " << hit->GetPlaneNum() << endl;
1302 // //cout << "Wire: " <<hit->GetWireNum() << endl;
1303 // //cout << "Sigma: " << hit->GetWireSigma() << endl;
1304 // }
1305
1306 } //end hit loop
1307 }
1308 /*
1309 for(Int_t ihit=0;ihit < theDCTrack->GetNHits();ihit++) {
1310 cout << " hit = " << ihit << " " << coords[ihit] << " " << planes[ihit] ;
1311 }
1312 cout << " ****" << endl;
1313 */
1314 for(Int_t irayp=0;irayp<NUM_FPRAY;irayp++) {
1315 for(Int_t jrayp=0;jrayp<NUM_FPRAY;jrayp++) {
1316 AA[irayp][jrayp] = 0.0;
1317 if(jrayp<irayp) { // Symmetric
1318 AA[irayp][jrayp] = AA[jrayp][irayp];
1319 } else {
1320 for(Int_t ihit=0;ihit < theDCTrack->GetNHits();ihit++) {
1321
1322 THcDCHit* hit=theDCTrack->GetHit(ihit);
1323
1324
1325 AA[irayp][jrayp] += fPlaneCoeffs[planes[ihit]][raycoeffmap[irayp]]*
1326 fPlaneCoeffs[planes[ihit]][raycoeffmap[jrayp]]/
1327 pow(hit->GetWireSigma(),2);
1328
1329 } //end ihit loop
1330 }
1331 }
1332 }
1333
1334 // Solve 4x4 equations
1335 TVectorD dray(NUM_FPRAY);
1336 // Should check that it is invertable
1337 AA.Invert();
1338 dray = AA*TT;
1339
1340 // Make sure fCoords, fResiduals, and fDoubleResiduals are clear
1341 for(Int_t iplane=0;iplane < fNPlanes; iplane++) {
1342 Double_t coord=0.0;
1343 for(Int_t ir=0;ir<NUM_FPRAY;ir++) {
1344 coord += fPlaneCoeffs[iplane][raycoeffmap[ir]]*dray[ir];
1345 // cout << "ir = " << ir << ", dray[ir] = " << dray[ir] << endl;
1346 }
1347 theDCTrack->SetCoord(iplane,coord);
1348 }
1349 // Compute Chi2 and residuals
1350 chi2 = 0.0;
1351 for(Int_t ihit=0;ihit < theDCTrack->GetNHits();ihit++) {
1352
1353 THcDCHit* hit=theDCTrack->GetHit(ihit);
1354
1355
1356 Double_t residual = coords[ihit] - theDCTrack->GetCoord(planes[ihit]);
1357 theDCTrack->SetResidual(planes[ihit], residual);
1358 if (abs(residual) > MaxResid) {
1359 MaxResidHitIndex =ihit;
1360 MaxResid = abs(residual);
1361 }
1362
1363 // double track_coord = theDCTrack->GetCoord(planes[ihit]);
1364//cout<<planes[ihit]<<"\t"<<track_coord<<"\t"<<coords[ihit]<<"\t"<<residual<<endl;
1365 chi2 += pow(residual/hit->GetWireSigma(),2);
1366
1367 }
1368
1369 theDCTrack->SetVector(dray[0], dray[1], 0.0, dray[2], dray[3]);
1370 }
1371 theDCTrack->SetChisq(chi2);
1372 // calculate ray without a plane in track
1373 for(Int_t ipl_hit=0;ipl_hit < theDCTrack->GetNHits();ipl_hit++) {
1374
1375
1376 if(theDCTrack->GetNFree() > 0) {
1377 TVectorD TT(NUM_FPRAY);
1379 for(Int_t irayp=0;irayp<NUM_FPRAY;irayp++) {
1380 TT[irayp] = 0.0;
1381 for(Int_t ihit=0;ihit < theDCTrack->GetNHits();ihit++) {
1382
1383
1384 THcDCHit* hit=theDCTrack->GetHit(ihit);
1385
1386 if (ihit != ipl_hit) {
1387 TT[irayp] += (coords[ihit]*
1388 fPlaneCoeffs[planes[ihit]][raycoeffmap[irayp]])
1389 /pow(hit->GetWireSigma(),2);
1390
1391 }
1392 }
1393 }
1394 for(Int_t irayp=0;irayp<NUM_FPRAY;irayp++) {
1395 for(Int_t jrayp=0;jrayp<NUM_FPRAY;jrayp++) {
1396 AA[irayp][jrayp] = 0.0;
1397 if(jrayp<irayp) { // Symmetric
1398 AA[irayp][jrayp] = AA[jrayp][irayp];
1399 } else {
1400
1401 for(Int_t ihit=0;ihit < theDCTrack->GetNHits();ihit++) {
1402
1403 THcDCHit* hit=theDCTrack->GetHit(ihit);
1404
1405
1406 if (ihit != ipl_hit) {
1407 AA[irayp][jrayp] += fPlaneCoeffs[planes[ihit]][raycoeffmap[irayp]]*
1408 fPlaneCoeffs[planes[ihit]][raycoeffmap[jrayp]]/
1409 pow(hit->GetWireSigma(),2);
1410
1411 }
1412 }
1413 }
1414 }
1415 }
1416 //
1417 // Solve 4x4 equations
1418 // Should check that it is invertable
1419 TVectorD dray(NUM_FPRAY);
1420 AA.Invert();
1421 dray = AA*TT;
1422 Double_t coord=0.0;
1423 for(Int_t ir=0;ir<NUM_FPRAY;ir++) {
1424 coord += fPlaneCoeffs[planes[ipl_hit]][raycoeffmap[ir]]*dray[ir];
1425 // cout << "ir = " << ir << ", dray[ir] = " << dray[ir] << endl;
1426 }
1427 Double_t residual = coords[ipl_hit] - coord;
1428 theDCTrack->SetResidualExclPlane(planes[ipl_hit], residual);
1429 }
1430 }
1431 //
1432 if (fTrackLargeResidCut == -1) {
1433 CheckResid = kFALSE;
1434 } else {
1435 CheckResid = kFALSE;
1436 if ( MaxResid > fTrackLargeResidCut) {
1437 CheckResid = kTRUE;
1438 theDCTrack->RemoveHit(MaxResidHitIndex);
1439 }
1440 }
1441 //
1442 if (NLoops == NLoopsMax) CheckResid = kFALSE;
1443 NLoops++;
1444 } // while
1445 } // loop over tracks
1446 //Calculate residual without plane
1447
1448 // Calculate residuals for each chamber if in single stub mode
1449 // and there was a track found in each chamber
1450 // Specific for two chambers. Can/should it be generalized?
1451
1452 if(fSingleStub != 0) {
1453 if(fNDCTracks == 2) {
1454 THcDCTrack *theDCTrack1 = static_cast<THcDCTrack*>( fDCTracks->At(0));
1455 THcDCTrack *theDCTrack2 = static_cast<THcDCTrack*>( fDCTracks->At(1));
1456 // Int_t itrack=0;
1457 Int_t ihit=0;
1458 THcDCHit* hit=theDCTrack1->GetHit(ihit);
1459 Int_t plane=hit->GetPlaneNum()-1;
1460 Int_t chamber=fNChamber[plane];
1461 if(chamber==1) {
1462 // itrack=1;
1463 hit=theDCTrack2->GetHit(ihit);
1464 plane=hit->GetPlaneNum()-1;
1465 chamber=fNChamber[plane];
1466 if(chamber==2) {
1467 Double_t ray1[4];
1468 Double_t ray2[4];
1469 theDCTrack1->GetRay(ray1);
1470 theDCTrack2->GetRay(ray2);
1471 // itrack = 1;
1472 // Loop over hits in second chamber
1473 for(Int_t ihit=0;ihit < theDCTrack2->GetNHits();ihit++) {
1474 // Calculate residual in second chamber from first chamber track
1475 THcDCHit* hit=theDCTrack2->GetHit(ihit);
1476 Int_t plane=hit->GetPlaneNum()-1;
1477 Double_t pos = DpsiFun(ray1,plane);
1478 Double_t coord;
1479 if(fFixLR) {
1481 coord = hit->GetPos()
1482 + theDCTrack2->GetHitLR(ihit)*theDCTrack2->GetHitDist(ihit);
1483 } else {
1484 coord = hit->GetPos()
1485 + theDCTrack2->GetHitLR(ihit)*hit->GetDist();
1486 }
1487 } else {
1489 coord = hit->GetPos()
1490 + hit->GetLR()*theDCTrack2->GetHitDist(ihit);
1491 } else {
1492 coord = hit->GetCoord();
1493 }
1494 }
1495 theDCTrack1->SetDoubleResidual(plane,coord - pos);
1496 // hdc_dbl_res(pln) = hdc_double_residual(1,pln) for hists
1497 }
1498 // itrack=0;
1499 // Loop over hits in first chamber
1500 for(Int_t ihit=0;ihit < theDCTrack1->GetNHits();ihit++) {
1501 // Calculate residual in first chamber from second chamber track
1502 THcDCHit* hit=theDCTrack1->GetHit(ihit);
1503 Int_t plane=hit->GetPlaneNum()-1;
1504 Double_t pos = DpsiFun(ray1,plane);
1505 Double_t coord;
1506 if(fFixLR) {
1508 coord = hit->GetPos()
1509 + theDCTrack1->GetHitLR(ihit)*theDCTrack1->GetHitDist(ihit);
1510 } else {
1511 coord = hit->GetPos()
1512 + theDCTrack1->GetHitLR(ihit)*hit->GetDist();
1513 }
1514 } else {
1516 coord = hit->GetPos()
1517 + hit->GetLR()*theDCTrack1->GetHitDist(ihit);
1518 } else {
1519 coord = hit->GetCoord();
1520 }
1521 }
1522 theDCTrack2->SetDoubleResidual(plane,coord - pos);
1523 // hdc_dbl_res(pln) = hdc_double_residual(1,pln) for hists
1524 }
1525 }
1526 }
1527 }
1528 }
1529 //
1530
1531 //
1532}
1533//
1535 printf("%5s %-14s %-14s %-14s %-14s %-10s %-10s \n","Track","x_t","y_t","xp_t","yp_t","chi2","DOF");
1536 printf("%5s %-14s %-14s %-14s %-14s %-10s %-10s \n"," ","[cm]","[cm]","[rad]","[rad]"," "," ");
1537 for(UInt_t itr=0;itr < fNDCTracks;itr++) {
1538 THcDCTrack *theDCTrack = static_cast<THcDCTrack*>( fDCTracks->At(itr));
1539 printf("%-5d %14.6e %14.6e %14.6e %14.6e %10.3e %3d \n", itr+1,theDCTrack->GetX(),theDCTrack->GetY(),theDCTrack->GetXP(),theDCTrack->GetYP(),theDCTrack->GetChisq(),theDCTrack->GetNFree());
1540 }
1541 for(UInt_t itr=0;itr < fNDCTracks;itr++) {
1542 printf("%s %5d \n","Hit info for track number = ",itr+1);
1543 printf("%5s %-15s %-15s %-15s \n","Plane","WIRE_COORD","Fit postiion","Residual");
1544 THcDCTrack *theDCTrack = static_cast<THcDCTrack*>( fDCTracks->At(itr));
1545 for(Int_t ihit=0;ihit < theDCTrack->GetNHits();ihit++) {
1546 THcDCHit* hit=theDCTrack->GetHit(ihit);
1547 Int_t plane=hit->GetPlaneNum()-1;
1548 Double_t coords_temp;
1549 if(fFixLR) {
1551 coords_temp = hit->GetPos()
1552 + theDCTrack->GetHitLR(ihit)*theDCTrack->GetHitDist(ihit);
1553 } else {
1554 coords_temp = hit->GetPos()
1555 + theDCTrack->GetHitLR(ihit)*hit->GetDist();
1556 }
1557 } else {
1559 coords_temp = hit->GetPos()
1560 + hit->GetLR()*theDCTrack->GetHitDist(ihit);
1561 } else {
1562 coords_temp = hit->GetCoord();
1563 }
1564 }
1565 printf("%-5d %15.7e %15.7e %15.7e \n",plane+1,coords_temp,theDCTrack->GetCoord(plane),theDCTrack->GetResidual(plane));
1566 }
1567 }
1568}
1569//
1571{
1590 Double_t infinity = 1.0E+20;
1591 Double_t cinfinity = 1/infinity;
1593 ray[2]*ray[1]*fPlaneCoeffs[plane][0] +
1594 ray[3]*ray[0]*fPlaneCoeffs[plane][1] +
1595 ray[2]*fPlaneCoeffs[plane][2] +
1596 ray[3]*fPlaneCoeffs[plane][3] +
1597 ray[0]*fPlaneCoeffs[plane][4] +
1598 ray[1]*fPlaneCoeffs[plane][5];
1599 Double_t denom = ray[2]*fPlaneCoeffs[plane][6]
1600 + ray[3]*fPlaneCoeffs[plane][7]
1601 + fPlaneCoeffs[plane][8];
1602 if(TMath::Abs(denom) < cinfinity) {
1603 DpsiFun = infinity;
1604 } else {
1605 DpsiFun = DpsiFun/denom;
1606 }
1607 return(DpsiFun);
1608}
1609
1610//_____________________________________________________________________________
1612{
1613 // EffCalc();
1614 MissReport(Form("%s.%s", GetApparatus()->GetName(), GetName()));
1615 return 0;
1616}
1617
1618//_____________________________________________________________________________
1620{
1627 delete [] fNChamHits; fNChamHits = new Int_t [fNChambers];
1628 delete [] fPlaneEvents; fPlaneEvents = new Int_t [fNPlanes];
1629
1630 fTotEvents = 0;
1631 for(UInt_t i=0;i<fNChambers;i++) {
1632 fNChamHits[i] = 0;
1633 }
1634 for(Int_t i=0;i<fNPlanes;i++) {
1635 fPlaneEvents[i] = 0;
1636 }
1637 gHcParms->Define(Form("%sdc_tot_events",fPrefix),"Total DC Events",fTotEvents);
1638 gHcParms->Define(Form("%sdc_cham_hits[%d]",fPrefix,fNChambers),"N events with hits per chamber",*fNChamHits);
1639 gHcParms->Define(Form("%sdc_events[%d]",fPrefix,fNPlanes),"N events with hits per plane",*fPlaneEvents);
1640}
1641
1642//_____________________________________________________________________________
1644{
1649 fTotEvents++;
1650 for(UInt_t i=0;i<fNChambers;i++) {
1651 if(fChambers[i]->GetNHits()>0) fNChamHits[i]++;
1652 }
1653 for(Int_t i=0;i<fNPlanes;i++) {
1654 if(fPlanes[i]->GetNHits() > 0) fPlaneEvents[i]++;
1655 }
1656 return;
1657}
1658
int Int_t
unsigned int UInt_t
const Data_t kBig
bool Bool_t
const Bool_t kFALSE
double Double_t
const Bool_t kTRUE
Option_t Option_t TPoint TPoint const char mode
Option_t Option_t TPoint TPoint const char y2
Option_t Option_t TPoint TPoint const char y1
R__EXTERN class THaVarList * gHaVars
R__EXTERN class THaCutList * gHaCuts
ClassImp(VDC::AnalyticTTDConv) using namespace std
#define NUM_FPRAY
Definition THcDC.h:18
string::size_type pos
R__EXTERN class THcParmList * gHcParms
Definition THcGlobals.h:11
R__EXTERN class THcDetectorMap * gHcDetectorMap
Definition THcGlobals.h:12
Int_t fDebug
char * Form(const char *fmt,...)
#define snprintf
void Delete(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 const char * Here(const char *) const
Int_t RemoveVariables()
virtual Int_t Result(const char *cutname="", EWarnMode mode=kWarn)
UInt_t GetTotNumChan() const
THaDetMap * fDetMap
THaApparatus * GetApparatus() const
UInt_t GetEvNum() const
Double_t GetPcentral() const
virtual void TransportToLab(Double_t p, Double_t th, Double_t ph, TVector3 &pvect) const
void SetPvect(const TVector3 &pvect)
void Set(Double_t x, Double_t y, Double_t theta, Double_t phi)
Double_t GetTPhi() const
Double_t GetDp() const
void SetMomentum(Double_t p)
void SetDp(Double_t dp)
void SetTrkNum(Int_t n)
void SetFlag(UInt_t flag)
void SetTarget(Double_t x, Double_t y, Double_t theta, Double_t phi)
Double_t GetP() const
Double_t GetTTheta() const
void SetChi2(Double_t chi2, Int_t ndof)
virtual THaTrack * AddTrack(TClonesArray &tracks, Double_t x, Double_t y, Double_t theta, Double_t phi, THaTrackID *ID=nullptr)
THaVar * Define(const char *name, const Byte_t &var, const Int_t *count=nullptr)
virtual THaVar * Find(const char *name) const
const void * GetValuePointer() const
Drift chamber wire hit info.
Definition THcDCHit.h:16
Int_t GetPlaneNum() const
Definition THcDCHit.h:58
Double_t GetCoord() const
Definition THcDCHit.h:42
Int_t GetWireNum() const
Definition THcDCHit.h:36
Int_t GetLR()
Definition THcDCHit.h:55
Double_t GetPos() const
Definition THcDCHit.h:41
Double_t GetDist() const
Definition THcDCHit.h:40
Double_t GetWireSigma() const
Definition THcDCHit.h:35
Class representing a track found from linking DC Space points.
Definition THcDCTrack.h:19
Double_t GetXP() const
Definition THcDCTrack.h:44
Int_t GetNFree() const
Definition THcDCTrack.h:36
void SetResidual(Int_t ip, Double_t coord)
Definition THcDCTrack.h:51
Double_t GetYP() const
Definition THcDCTrack.h:45
Double_t GetResidualExclPlane(Int_t ip) const
Definition THcDCTrack.h:40
Double_t GetSp2_ID() const
Definition THcDCTrack.h:47
void SetChisq(Double_t chi2)
Definition THcDCTrack.h:56
void SetSp2_ID(Int_t isp2)
Definition THcDCTrack.h:58
virtual void RemoveHit(Int_t RemoveHitIndex)
virtual void AddSpacePoint(THcSpacePoint *sp)
THcDCHit * GetHit(Int_t i) const
Definition THcDCTrack.h:32
Double_t GetX() const
Definition THcDCTrack.h:42
Double_t GetHitDist(Int_t ihit)
Definition THcDCTrack.h:33
void GetRay(Double_t *ray) const
Definition THcDCTrack.h:41
Int_t GetNSpacePoints() const
Definition THcDCTrack.h:29
void SetNFree(Int_t nfree)
Definition THcDCTrack.h:49
Int_t GetHitLR(Int_t ihit)
Definition THcDCTrack.h:34
void SetCoord(Int_t ip, Double_t coord)
Definition THcDCTrack.h:50
Double_t GetCoord(Int_t ip) const
Definition THcDCTrack.h:37
void SetSp1_ID(Int_t isp1)
Definition THcDCTrack.h:57
void SetResidualExclPlane(Int_t ip, Double_t coord)
Definition THcDCTrack.h:52
Int_t GetNHits() const
Definition THcDCTrack.h:35
void SetVector(Double_t x, Double_t y, Double_t z, Double_t xp, Double_t yp)
Definition THcDCTrack.h:54
Double_t GetChisq() const
Definition THcDCTrack.h:48
void SetDoubleResidual(Int_t ip, Double_t coord)
Definition THcDCTrack.h:53
Double_t GetSp1_ID() const
Definition THcDCTrack.h:46
THcSpacePoint * GetSpacePoint(Int_t i) const
Definition THcDCTrack.h:31
Double_t GetResidual(Int_t ip) const
Definition THcDCTrack.h:38
void SetHitLR(Int_t ihit, Int_t lrtemp)
Definition THcDCTrack.h:59
Double_t GetY() const
Definition THcDCTrack.h:43
Analyze a package of horizontal drift chambers.
Definition THcDC.h:23
virtual Int_t End(THaRunBase *run=0)
Definition THcDC.cxx:1611
Bool_t fInSideDipoleExit_best
Definition THcDC.h:186
Bool_t * fPresentP
Definition THcDC.h:194
Int_t GetNChamber(Int_t plane) const
Definition THcDC.h:45
virtual Int_t ReadDatabase(const TDatime &date)
Definition THcDC.cxx:269
Int_t fNSp
Definition THcDC.h:123
Int_t fTDC_RefTimeCut
Definition THcDC.h:96
Int_t fdebugprintrawdc
Definition THcDC.h:90
Int_t fUseNewTrackFit
Definition THcDC.h:142
Int_t fdebugflagstubs
Definition THcDC.h:92
Int_t fTotEvents
Definition THcDC.h:188
Double_t * fLR_best
Definition THcDC.h:126
Int_t fN_True_RawHits
Definition THcDC.h:122
virtual EStatus Init(const TDatime &run_time)
Definition THcDC.cxx:195
Int_t fSp2_ID_best
Definition THcDC.h:185
Int_t * fNWires
Definition THcDC.h:159
Double_t * fZPos
Definition THcDC.h:169
Int_t fFixPropagationCorrection
Definition THcDC.h:110
void NewLinkStubs()
Definition THcDC.cxx:778
Double_t * fBetaAngle
Definition THcDC.h:171
void PrintStubs()
Definition THcDC.cxx:763
Int_t fUseNewLinkStubs
Definition THcDC.h:141
Double_t fY_fp_best
Definition THcDC.h:180
Double_t * fGammaAngle
Definition THcDC.h:172
void Setup(const char *name, const char *description)
Definition THcDC.cxx:94
Int_t * fReadoutTB
Definition THcDC.h:163
Int_t * fPlaneEvents
Definition THcDC.h:190
THcDC()
Definition THcDC.cxx:188
void ClearEvent()
Definition THcDC.cxx:519
Int_t fNPlanes
Definition THcDC.h:105
Int_t fVersion
Definition THcDC.h:143
Int_t fNTracksMaxFP
Definition THcDC.h:136
void NewTrackFit(UInt_t TrackIndex)
Definition THcDC.cxx:1111
Double_t ** fPlaneCoeffs
Definition THcDC.h:177
Double_t fNSperChan
Definition THcDC.h:133
Double_t * fYPos
Definition THcDC.h:168
Int_t fProjectToChamber
Definition THcDC.h:114
Double_t * fCentralWire
Definition THcDC.h:174
Double_t fXtTrCriterion
Definition THcDC.h:137
Double_t fXptTrCriterion
Definition THcDC.h:139
Int_t * fMinHits
Definition THcDC.h:148
virtual ~THcDC()
Definition THcDC.cxx:447
std::vector< THcDriftChamberPlane * > fPlanes
Definition THcDC.h:202
Double_t * fSpace_Point_Criterion
Definition THcDC.h:151
Double_t fWireVelocity
Definition THcDC.h:134
Int_t fStubTest
Definition THcDC.h:119
TClonesArray * fDCTracks
Definition THcDC.h:100
std::vector< THcDriftChamber * > fChambers
Definition THcDC.h:203
Int_t fFixLR
Definition THcDC.h:108
char fPrefix[2]
Definition THcDC.h:104
void EffInit()
Definition THcDC.cxx:1619
virtual Int_t ApplyCorrections(void)
Definition THcDC.cxx:601
Int_t * fWireOrder
Definition THcDC.h:161
Double_t * fPlaneTimeZero
Definition THcDC.h:175
Int_t * fMinCombos
Definition THcDC.h:150
Double_t fYptTrCriterion
Definition THcDC.h:140
Int_t fSingleStub
Definition THcDC.h:135
void Eff()
Definition THcDC.cxx:1643
Double_t fX_fp_best
Definition THcDC.h:179
void TrackFit()
Definition THcDC.cxx:1241
Double_t * fResiduals
Definition THcDC.h:128
Int_t fSp1_ID_best
Definition THcDC.h:184
static const UInt_t MAXTRACKS
Definition THcDC.h:200
Double_t * fResidualsExclPlane
Definition THcDC.h:129
Double_t DpsiFun(Double_t ray[4], Int_t plane)
Definition THcDC.cxx:1570
Double_t fChisq_best
Definition THcDC.h:183
void PrintSpacePoints()
Definition THcDC.cxx:744
Int_t * fTdcWinMin
Definition THcDC.h:156
Double_t * fAlphaAngle
Definition THcDC.h:170
Double_t * fWire_hit_should
Definition THcDC.h:131
virtual Int_t Decode(const THaEvData &)
Definition THcDC.cxx:552
Double_t * fPitch
Definition THcDC.h:173
Int_t * fNChamber
Definition THcDC.h:160
Int_t * fReadoutLR
Definition THcDC.h:164
Int_t fdebugprintdecodeddc
Definition THcDC.h:94
Double_t * fSigma
Definition THcDC.h:176
void EfficiencyPerWire(Int_t golden_track_index)
Definition THcDC.cxx:725
virtual Int_t CoarseTrack(TClonesArray &tracks)
Definition THcDC.cxx:607
Int_t fdebugflagpr
Definition THcDC.h:91
Double_t * fXCenter
Definition THcDC.h:146
Int_t * fMaxHits
Definition THcDC.h:149
Double_t * fCentralTime
Definition THcDC.h:158
void FitLineToTrack(Int_t TrackHits, Double_t coords[], Int_t planes[], Double_t wiresigma[], Double_t TrackCoord[], Double_t save_ray[])
Definition THcDC.cxx:1071
Double_t * fXPos
Definition THcDC.h:167
Int_t fHMSStyleChambers
Definition THcDC.h:95
Double_t fXp_fp_best
Definition THcDC.h:181
UInt_t fNChambers
Definition THcDC.h:107
Double_t fTrackLargeResidCut
Definition THcDC.h:97
Int_t * fDriftTimeSign
Definition THcDC.h:162
Double_t * fPos_best
Definition THcDC.h:127
Int_t * fTdcWinMax
Definition THcDC.h:157
Double_t * fWire_hit_did
Definition THcDC.h:130
virtual Int_t FineTrack(TClonesArray &tracks)
Definition THcDC.cxx:695
Int_t fUseNewFindSpacePoints
Definition THcDC.h:88
Int_t fNhits
Definition THcDC.h:120
UInt_t fNDCTracks
Definition THcDC.h:99
void PrintTrack()
Definition THcDC.cxx:1534
Double_t * fDist_best
Definition THcDC.h:125
void LinkStubs()
Definition THcDC.cxx:859
Int_t fNsp_best
Definition THcDC.h:124
Int_t * fNChamHits
Definition THcDC.h:189
Double_t fYp_fp_best
Definition THcDC.h:182
virtual Int_t DefineVariables(EMode mode=kDefine)
Definition THcDC.cxx:398
Double_t fYtTrCriterion
Definition THcDC.h:138
Int_t fdebugtrackprint
Definition THcDC.h:93
Int_t fdebuglinkstubs
Definition THcDC.h:89
Double_t * fYCenter
Definition THcDC.h:147
Int_t fNthits
Definition THcDC.h:121
char ** fPlaneNames
Definition THcDC.h:106
void DeleteArrays()
Definition THcDC.cxx:467
void SetFocalPlaneBestTrack(Int_t golden_track_index)
Definition THcDC.cxx:700
virtual Int_t FillMap(THaDetMap *detmap, const char *detectorname)
Class for a a single Hall C horizontal drift chamber plane.
Subdetector class for a single drift chamber with several planes.
void SetHMSStyleFlag(Int_t flag)
A standard Hall C spectrometer apparatus.
virtual void CalculateTargetQuantities(THaTrack *track, Double_t &gbeam_y, Double_t &xptar, Double_t &ytar, Double_t &yptar, Double_t &delta)
Transport focal plane track to target.
virtual Int_t DecodeToHitList(const THaEvData &evdata, Bool_t suppress=kFALSE)
Populate the hitlist from the raw event data.
TClonesArray * fRawHitList
Definition THcHitList.h:51
void CreateMissReportParms(const char *prefix) const
void InitHitList(THaDetMap *detmap, const char *hitclass, Int_t maxhits, Int_t tdcref_cut=0, Int_t adcref_cut=0)
Save the electronics module to detector mapping and initialize a hit array of hits of class hitclass.
UInt_t fNRawHits
Definition THcHitList.h:45
void MissReport(const char *name) const
Int_t LoadParmValues(const DBRequest *list, const char *prefix="")
Retrieve parameter values from the parameter cache.
Class representing for drift chamber wire (or other device with a single multihit TDC channel per det...
Definition THcRawDCHit.h:8
Class representing a single hit DC.
Double_t GetY()
Int_t fNChamber_spnum
THcDCHit * GetHit(Int_t ihit)
Double_t * GetStubP()
Bool_t GetSetStubFlag()
Double_t GetX()
Int_t GetCombos()
TMatrixT< Element > & Invert(Double_t *det=nullptr)
const char * GetName() const override
const char * GetTitle() const override
TObject * At(Int_t idx) const override
R__ALWAYS_INLINE Bool_t IsZombie() const
virtual void Error(const char *method, const char *msgfmt,...) const
void MakeZombie()
RVec< PromoteType< T > > abs(const RVec< T > &v)
RVec< PromoteType< T > > round(const RVec< T > &v)
RVec< PromoteTypes< T0, T1 > > pow(const RVec< T0 > &v, const T1 &y)
Double_t Abs(Double_t d)
STL namespace.
void tracks()