Hall C ROOT/C++ Analyzer (hcana)
Loading...
Searching...
No Matches
THcDriftChamberPlane.cxx
Go to the documentation of this file.
1
8#include "THcDC.h"
10#include "THcDCWire.h"
11#include "THcDCHit.h"
12#include "THcDCLookupTTDConv.h"
13#include "THcSignalHit.h"
14#include "THcGlobals.h"
15#include "THcParmList.h"
16#include "THcHitList.h"
17#include "THaApparatus.h"
18#include "THcHodoscope.h"
19#include "TClass.h"
20
21
22#include <cstring>
23#include <cstdio>
24#include <cstdlib>
25#include <iostream>
26
27using namespace std;
28
30
31//______________________________________________________________________________
33 const char* description,
34 const Int_t planenum,
35 THaDetectorBase* parent )
36: THaSubDetector(name,description,parent), fTTDConv(0)
37{
38 // Normal constructor with name and description
39 fHits = new TClonesArray("THcDCHit",100);
40 fFirstPassHits = new TClonesArray("THcDCHit",100);
41 fRawHits = new TClonesArray("THcDCHit",100);
42 fWires = new TClonesArray("THcDCWire", 100);
43
44 fPlaneNum = planenum;
45}
46
47//_____________________________________________________________________________
50{
51 // Constructor
52 fHits = NULL;
53 fFirstPassHits = NULL;
54 fRawHits = NULL;
55 fWires = NULL;
56 fTTDConv = NULL;
57}
58//______________________________________________________________________________
60{
61 // Destructor
62 delete [] fTzeroWire;
63 delete [] fSigmaWire;
64 delete fWires;
65 delete fHits;
66 delete fFirstPassHits;
67 delete fRawHits;
68 delete fTTDConv;
69
70}
72{
73 // Extra initialization for scintillator plane: set up DataDest map
74
75 // cout << "THcDriftChamberPlane::Init called " << GetName() << endl;
76
77 if( IsZombie())
78 return fStatus = kInitError;
79
80 // How to get information for parent
81 // if( GetParent() )
82 // fOrigin = GetParent()->GetOrigin();
83
84 EStatus status;
85 if( (status=THaSubDetector::Init( date )) )
86 return fStatus = status;
87
88 return fStatus = kOK;
89
90}
91
92//_____________________________________________________________________________
94{
95
102 char prefix[2];
103 UInt_t NumDriftMapBins;
104 Double_t DriftMapFirstBin;
105 Double_t DriftMapBinSize;
108 prefix[0]=tolower(GetParent()->GetPrefix()[0]);
109 prefix[1]='\0';
110 DBRequest list[]={
111 {"driftbins", &NumDriftMapBins, kInt},
112 {"drift1stbin", &DriftMapFirstBin, kDouble},
113 {"driftbinsz", &DriftMapBinSize, kDouble},
114 {"_using_tzero_per_wire", &fUsingTzeroPerWire, kInt,0,1},
115 {"_using_sigma_per_wire", &fUsingSigmaPerWire, kInt,0,1},
116 {"min_drifttime", &fMin_DriftTime, kDouble,0,1},
117 {"max_drifttime", &fMax_DriftTime, kDouble,0,1},
118 {0}
119 };
120
121
122 fMin_DriftTime = -50.; //ns
123 fMax_DriftTime = 200.; //ns
124
125 gHcParms->LoadParmValues((DBRequest*)&list,prefix);
126
127 Double_t *DriftMap = new Double_t[NumDriftMapBins];
128 DBRequest list2[]={
129 {Form("wc%sfract",GetName()),DriftMap,kDouble,NumDriftMapBins},
130 {0}
131 };
132 gHcParms->LoadParmValues((DBRequest*)&list2,prefix);
133
134
135 // Retrieve parameters we need from parent class
136 THcDC* fParent;
137
138 fParent = (THcDC*) GetParent();
139 // These are single variables here, but arrays in THcDriftChamber.
140 fSigma = fParent->GetSigma(fPlaneNum);
141 fChamberNum = fParent->GetNChamber(fPlaneNum);
142 fNWires = fParent->GetNWires(fPlaneNum);
143 fWireOrder = fParent->GetWireOrder(fPlaneNum);
144 fPitch = fParent->GetPitch(fPlaneNum);
145 fCentralWire = fParent->GetCentralWire(fPlaneNum);
146 fTdcWinMin = fParent->GetTdcWinMin(fPlaneNum);
147 fTdcWinMax = fParent->GetTdcWinMax(fPlaneNum);
148 fPlaneTimeZero = fParent->GetPlaneTimeZero(fPlaneNum);
149 fCenter = fParent->GetCenter(fPlaneNum);
150 fCentralTime = fParent->GetCentralTime(fPlaneNum);
151 fDriftTimeSign = fParent->GetDriftTimeSign(fPlaneNum);
152 fReadoutLR = fParent->GetReadoutLR(fPlaneNum);
153 fReadoutTB = fParent->GetReadoutTB(fPlaneNum);
154 fVersion = fParent->GetVersion();
155
156 fNSperChan = fParent->GetNSperChan();
157
158
161
162
163 if (fUsingTzeroPerWire==1) {
164 DBRequest list3[]={
165 {Form("tzero%s",GetName()),fTzeroWire,kDouble,(UInt_t) fNWires},
166 {0}
167 };
168 gHcParms->LoadParmValues((DBRequest*)&list3,prefix);
169
170 } else {
171 for (Int_t iw=0;iw < fNWires;iw++) {
172 fTzeroWire[iw]=0.0;
173 }
174 }
175
176 if (fUsingSigmaPerWire==1) {
177 DBRequest list4[]={
178 {Form("wire_sigma%s",GetName()),fSigmaWire,kDouble,(UInt_t) fNWires},
179 {0}
180 };
181 gHcParms->LoadParmValues((DBRequest*)&list4,prefix);
182 } else {
183 for (Int_t iw=0;iw < fNWires;iw++) {
184 fSigmaWire[iw]=fSigma;
185 }
186 }
187 // Calculate Geometry Constants
188 // Do we want to move all this to the Chamber of DC Package leve
189 // as that is where these things will be needed?
190 Double_t z0 = fParent->GetZPos(fPlaneNum);
191 Double_t alpha = fParent->GetAlphaAngle(fPlaneNum);
192 Double_t beta = fParent->GetBetaAngle(fPlaneNum);
193 fBeta = beta;
194 Double_t gamma = fParent->GetGammaAngle(fPlaneNum);
195 Double_t cosalpha = TMath::Cos(alpha);
196 Double_t sinalpha = TMath::Sin(alpha);
197 Double_t cosbeta = TMath::Cos(beta);
198 Double_t sinbeta = TMath::Sin(beta);
199 Double_t cosgamma = TMath::Cos(gamma);
200 Double_t singamma = TMath::Sin(gamma);
201
202 Double_t hzchi = -cosalpha*sinbeta + sinalpha*cosbeta*singamma;
203 Double_t hzpsi = sinalpha*sinbeta + cosalpha*cosbeta*singamma;
204 Double_t hxchi = -cosalpha*cosbeta - sinalpha*sinbeta*singamma;
205 Double_t hxpsi = sinalpha*cosbeta - cosalpha*sinbeta*singamma;
206 Double_t hychi = sinalpha*cosgamma;
207 Double_t hypsi = cosalpha*cosgamma;
208 Double_t stubxchi = -cosalpha;
209 Double_t stubxpsi = sinalpha;
210 Double_t stubychi = sinalpha;
211 Double_t stubypsi = cosalpha;
212
213 if(cosalpha <= 0.707) { // x-like wire, need dist from x=0 line
214 fReadoutX = 1;
215 fReadoutCorr = 1/sinalpha;
216 } else {
217 fReadoutX = 0;
218 fReadoutCorr = 1/cosalpha;
219 }
220
221 Double_t sumsqupsi = hzpsi*hzpsi+hxpsi*hxpsi+hypsi*hypsi;
222 Double_t sumsquchi = hzchi*hzchi+hxchi*hxchi+hychi*hychi;
223 Double_t sumcross = hzpsi*hzchi + hxpsi*hxchi + hypsi*hychi;
224 Double_t denom1 = sumsqupsi*sumsquchi-sumcross*sumcross;
225 fPsi0 = (-z0*hzpsi*sumsquchi
226 +z0*hzchi*sumcross) / denom1;
227 Double_t hchi0 = (-z0*hzchi*sumsqupsi
228 +z0*hzpsi*sumcross) / denom1;
229 Double_t hphi0 = TMath::Sqrt(pow(z0+hzpsi*fPsi0+hzchi*hchi0,2)
230 + pow(hxpsi*fPsi0+hxchi*hchi0,2)
231 + pow(hypsi*fPsi0+hychi*hchi0,2) );
232 if(z0 < 0.0) hphi0 = -hphi0;
233
234 Double_t denom2 = stubxpsi*stubychi - stubxchi*stubypsi;
235
236 // Why are there 4, but only 3 used?
237 fStubCoef[0] = stubychi/(fSigma*denom2); // sin(a)/sigma
238 fStubCoef[1] = -stubxchi/(fSigma*denom2); // cos(a)/sigma
239 fStubCoef[2] = hphi0*fStubCoef[0]; // z0*sin(a)/sig
240 fStubCoef[3] = hphi0*fStubCoef[1]; // z0*cos(a)/sig
241
242 fXsp = hychi/denom2; // sin(a)
243 fYsp = -hxchi/denom2; // cos(a)
244
245 // Comput track fitting coefficients
246#define LOCRAYZT 0.0
247 fPlaneCoef[0]= hzchi; // = 0.
248 fPlaneCoef[1]=-hzchi; // = 0.
249 fPlaneCoef[2]= hychi*(z0-LOCRAYZT); // sin(a)*(z-hlocrayzt)
250 fPlaneCoef[3]= hxchi*(LOCRAYZT-z0); // cos(a)*(z-hlocrayzt)
251 fPlaneCoef[4]= hychi; // sin(a)
252 fPlaneCoef[5]=-hxchi; // cos(a)
253 fPlaneCoef[6]= hzchi*hypsi - hychi*hzpsi; // 0.
254 fPlaneCoef[7]=-hzchi*hxpsi + hxchi*hzpsi; // 0.
255 fPlaneCoef[8]= hychi*hxpsi - hxchi*hypsi; // 1.
256
257
258 fTTDConv = new THcDCLookupTTDConv(DriftMapFirstBin,fPitch/2,DriftMapBinSize,
259 NumDriftMapBins,DriftMap);
260 delete [] DriftMap;
261
262 Int_t nWires = fParent->GetNWires(fPlaneNum);
263 // For HMS, wire numbers start with one, but arrays start with zero.
264 // So wire number is index+1
265 for (int i=0; i<nWires; i++) {
266 Double_t pos = fPitch*( (fWireOrder==0?(i+1):fNWires-i)
268 Int_t readoutside = GetReadoutSide(i+1);
269 new((*fWires)[i]) THcDCWire( i+1, pos , fTzeroWire[i], fSigmaWire[i], readoutside, fTTDConv); //added fTzeroWire/fSigmaWire to be read in as fTOffset --Carlos
270 }
271
272 THaApparatus* app = GetApparatus();
273 const char* nm = "hod";
274 if( !app ||
275 !(fglHod = dynamic_cast<THcHodoscope*>(app->GetDetector(nm))) ) {
276 static const char* const here = "ReadDatabase()";
277 Warning(Here(here),"Hodoscope \"%s\" not found. "
278 "Event-by-event time offsets will NOT be used!!",nm);
279 }
280
281 return kOK;
282}
283//_____________________________________________________________________________
285{
286 // Initialize global variables and lookup table for decoder
287
288 // cout << "THcDriftChamberPlane::DefineVariables called " << GetName() << endl;
289
290 if( mode == kDefine && fIsSetup ) return kOK;
291 fIsSetup = ( mode == kDefine );
292
293 // Register variables in global list
294 RVarDef vars[] = {
295 {"raw.wirenum", "List of TDC wire number of all hits in DC",
296 "fRawHits.THcDCHit.GetWireNum()"},
297 {"wirenum", "List of TDC wire number (select first hit in TDc window",
298 "fHits.THcDCHit.GetWireNum()"},
299 {"rawnorefcorrtdc", "Raw TDC Values",
300 "fHits.THcDCHit.GetRawNoRefCorrTime()"},
301 {"rawtdc", "Raw TDC with reference time subtracted Values",
302 "fHits.THcDCHit.GetRawTime()"},
303 {"time","Drift times",
304 "fHits.THcDCHit.GetTime()"},
305 {"dist","Drift distancess",
306 "fHits.THcDCHit.GetDist()"},
307 {"pos"," Wire pos",
308 "fHits.THcDCHit.GetPos()"},
309 {"nhit", "Number of hits", "GetNHits()"},
310 {"RefTime", "TDC reference time", "fTdcRefTime"},
311 { 0 }
312 };
313
314 return DefineVarsFromList( vars, mode );
315}
316
317//_____________________________________________________________________________
319{
320 //cout << " Calling THcDriftChamberPlane::Clear " << GetName() << endl;
321 // Clears the hit lists
322 fHits->Clear();
324 fRawHits->Clear();
325}
326
327//_____________________________________________________________________________
329{
330 // Doesn't actually get called. Use Fill method instead
331 cout << " Calling THcDriftChamberPlane::Decode " << GetName() << endl;
332
333 return 0;
334}
335//_____________________________________________________________________________
337{
338
339 // HitCount();
340
341 return 0;
342}
343//
345 Double_t wire_num_calc=-1000;
346 if (fWireOrder==0) wire_num_calc = (pos+fCenter)/(fPitch)+fCentralWire;
347 if (fWireOrder==1) wire_num_calc = 1-((pos+fCenter)/(fPitch)+fCentralWire-fNWires);
348 return(wire_num_calc);
349}
350//_____________________________________________________________________________
352{
353 return 0;
354}
356{
363 fHits->Clear();
365 fRawHits->Clear();
367 Int_t nrawhits = rawhits->GetLast()+1;
368 fNRawhits=0;
369 Int_t ihit = nexthit;
370 Int_t nextHit = 0;
371 Int_t nextRawHit = 0;
372 while(ihit < nrawhits) {
373 THcRawDCHit* hit = (THcRawDCHit *) rawhits->At(ihit);
374 if(hit->fPlane > fPlaneNum) {
375 break;
376 }
377 Int_t wireNum = hit->fCounter;
378 THcDCWire* wire = GetWire(wireNum);
379 Bool_t First_Hit_In_Window = kTRUE;
381 for(UInt_t mhit=0; mhit<hit->GetRawTdcHit().GetNHits(); mhit++) {
382 fNRawhits++;
383 /* Sort into early, late and ontime */
384 Int_t rawnorefcorrtdc = hit->GetRawTdcHit().GetTimeRaw(mhit); // Get the ref time subtracted time
385 Int_t rawtdc = hit->GetRawTdcHit().GetTime(mhit); // Get the ref time subtracted time
386 Double_t time = - rawtdc*fNSperChan + fPlaneTimeZero - wire->GetTOffset(); // fNSperChan > 0 for 1877
387 new( (*fRawHits)[nextRawHit++] ) THcDCHit(wire, rawnorefcorrtdc,rawtdc, time, this);
388 if(rawtdc < fTdcWinMin) {
389 // Increment early counter (Actually late because TDC is backward)
390 } else if (rawtdc > fTdcWinMax) {
391 // Increment late count
392 } else {
393 if (First_Hit_In_Window) {
394 new( (*fFirstPassHits)[nextHit++] ) THcDCHit(wire, rawnorefcorrtdc,rawtdc, time, this);
395 First_Hit_In_Window = kFALSE;
396 }
397 }
398 }
399 ihit++;
400 }
401 return(ihit);
402}
404{
405 Double_t StartTime = 0.0;
406 Int_t nextHit=0;
407 if( fglHod ) StartTime = fglHod->GetStartTime();
408 if (StartTime == -1000) StartTime = 0.0;
409 for(Int_t ihit=0;ihit<(fFirstPassHits->GetLast()+1);ihit++) {
410 THcDCHit *thishit = (THcDCHit*) fFirstPassHits->At(ihit);
411 Double_t temptime= thishit->GetTime()-StartTime;
412 Int_t tempRawtime= thishit->GetRawTime();
413 Int_t tempRawNoRefCorrtime= thishit->GetRawNoRefCorrTime();
414 THcDCWire *tempWire = (THcDCWire*) thishit->GetWire();
415 thishit->SetTime(temptime);
416 thishit->ConvertTimeToDist();
417 if (temptime > fMin_DriftTime && temptime <fMax_DriftTime) {
418 new( (*fHits)[nextHit++] ) THcDCHit(tempWire, tempRawNoRefCorrtime, tempRawtime, temptime, this);
419
420 }
421
422 }
423 return 0;
424}
426{
427 Int_t readoutside;
428 //if new HMS
429 if (fVersion == 1) {
430 if ((fPlaneNum>=3 && fPlaneNum<=4) || (fPlaneNum>=9 && fPlaneNum<=10)) {
431 if (fReadoutTB>0) {
432 if (wirenum < 60) {
433 readoutside = 2;
434 } else {
435 readoutside = 4;
436 }
437 } else {
438 if (wirenum < 44) {
439 readoutside = 4;
440 } else {
441 readoutside = 2;
442 }
443 }
444 } else {
445 if (fReadoutTB>0) {
446 if (wirenum < 51) {
447 readoutside = 2;
448 } else if (wirenum >= 51 && wirenum <= 64) {
449 readoutside = 1;
450 } else {
451 readoutside =4;
452 }
453 } else {
454 if (wirenum < 33) {
455 readoutside = 4;
456 } else if (wirenum >=33 && wirenum<=46) {
457 readoutside = 1;
458 } else {
459 readoutside = 2;
460 }
461 }
462 }
463 } else {//appplies SHMS DC configuration
464 //check if x board
465 if ((fPlaneNum>=3 && fPlaneNum<=4) || (fPlaneNum>=9 && fPlaneNum<=10)) {
466 if (fReadoutTB>0) {
467 if (wirenum < 49) {
468 readoutside = 4;
469 } else {
470 readoutside = 2;
471 }
472 } else {
473 if (wirenum < 33) {
474 readoutside = 2;
475 } else {
476 readoutside = 4;
477 }
478 }
479 } else { //else is u board
480 if (fReadoutTB>0) {
481 if (wirenum < 41) {
482 readoutside = 4;
483 } else if (wirenum >= 41 && wirenum <= 63) {
484 readoutside = 3;
485 } else if (wirenum >=64 && wirenum <=69) {
486 readoutside = 1;
487 } else {
488 readoutside = 2;
489 }
490 } else {
491 if (wirenum < 39) {
492 readoutside = 2;
493 } else if (wirenum >=39 && wirenum<=44) {
494 readoutside = 1;
495 } else if (wirenum>=45 && wirenum<=67) {
496 readoutside = 3;
497 } else {
498 readoutside = 4;
499 }
500 }
501 }
502 }
503 return(readoutside);
504}
int Int_t
unsigned int UInt_t
const Data_t kBig
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 mode
ClassImp(VDC::AnalyticTTDConv) using namespace std
string::size_type pos
#define LOCRAYZT
R__EXTERN class THcParmList * gHcParms
Definition THcGlobals.h:11
char * Form(const char *fmt,...)
void Clear(Option_t *option="") override
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
const char * GetPrefix() const
THaDetectorBase * GetParent() const
THaApparatus * GetApparatus() const
Drift chamber wire hit info.
Definition THcDCHit.h:16
Int_t GetRawTime() const
Definition THcDCHit.h:37
virtual Double_t ConvertTimeToDist()
Definition THcDCHit.cxx:36
Int_t GetRawNoRefCorrTime() const
Definition THcDCHit.h:38
void SetTime(Double_t time)
Definition THcDCHit.h:50
THcDCWire * GetWire() const
Definition THcDCHit.h:34
Double_t GetTime() const
Definition THcDCHit.h:39
Drift time to distance conversion via lookup table.
Class representing a drift chamber wire.
Definition THcDCWire.h:13
Double_t GetTOffset() const
Definition THcDCWire.h:28
Analyze a package of horizontal drift chambers.
Definition THcDC.h:23
Class for a a single Hall C horizontal drift chamber plane.
virtual void Clear(Option_t *opt="")
Double_t CalcWireFromPos(Double_t pos)
virtual Int_t DefineVariables(EMode mode=kDefine)
virtual Int_t FineProcess(TClonesArray &tracks)
virtual Int_t Decode(const THaEvData &)
THcDCWire * GetWire(Int_t i) const
virtual Int_t CoarseProcess(TClonesArray &tracks)
virtual Int_t GetReadoutSide(Int_t wirenum)
THcDCTimeToDistConv * fTTDConv
virtual Int_t ProcessHits(TClonesArray *rawhits, Int_t nexthit)
virtual Int_t ReadDatabase(const TDatime &date)
Generic hodoscope consisting of multiple planes with multiple paddles with phototubes on both ends.
Double_t GetStartTime() 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
THcRawTdcHit & GetRawTdcHit()
Int_t fPlane
Definition THcRawHit.h:52
Int_t fCounter
Definition THcRawHit.h:53
Int_t GetRefTime() const
Gets reference time. In channels.
Int_t GetTimeRaw(UInt_t iHit=0) const
Gets raw TDC time. In channels.
Int_t GetTime(UInt_t iHit=0) const
Gets TDC time. In channels.
Bool_t HasRefTime() const
Queries whether reference time has been set.
const char * GetName() const override
TObject * At(Int_t idx) const override
Int_t GetLast() const override
virtual void Warning(const char *method, const char *msgfmt,...) const
R__ALWAYS_INLINE Bool_t IsZombie() const
double beta(double x, double y)
RVec< PromoteTypes< T0, T1 > > pow(const RVec< T0 > &v, const T1 &y)
double gamma(double x)
Double_t Sqrt(Double_t x)
Double_t Cos(Double_t)
Double_t Sin(Double_t)
STL namespace.