Hall C ROOT/C++ Analyzer (hcana)
Loading...
Searching...
No Matches
THcCherenkov.cxx
Go to the documentation of this file.
1
8#include "THcCherenkov.h"
9#include "THcHodoscope.h"
10#include "TClonesArray.h"
11#include "THcSignalHit.h"
12#include "THaEvData.h"
13#include "THaDetMap.h"
14#include "THcDetectorMap.h"
15#include "THcGlobals.h"
16#include "THaCutList.h"
17#include "THcParmList.h"
18#include "THcHitList.h"
19#include "THaApparatus.h"
20#include "VarDef.h"
21#include "VarType.h"
22#include "THaTrack.h"
23#include "TClonesArray.h"
24#include "TMath.h"
25#include "THaTrackProj.h"
27
28#include <algorithm>
29#include <cstring>
30#include <cstdio>
31#include <cstdlib>
32#include <iostream>
33#include <string>
34#include <iomanip>
35
36
37using namespace std;
38
39using std::cout;
40using std::cin;
41using std::endl;
42using std::setw;
43using std::setprecision;
44
45//_____________________________________________________________________________
46THcCherenkov::THcCherenkov( const char* name, const char* description,
47 THaApparatus* apparatus ) :
48 THaNonTrackingDetector(name,description,apparatus)
49{
50 // Normal constructor with name and description
55 frAdcPed = new TClonesArray("THcSignalHit", MaxNumCerPmt*MaxNumAdcPulse);
60 // Normal constructor with name and description
69
70 fNumAdcHits = vector<Int_t> (MaxNumCerPmt, 0.0);
71 fNumGoodAdcHits = vector<Int_t> (MaxNumCerPmt, 0.0);
72 fNumTracksMatched = vector<Int_t> (MaxNumCerPmt, 0.0);
73 fNumTracksFired = vector<Int_t> (MaxNumCerPmt, 0.0);
74 fNpe = vector<Double_t> (MaxNumCerPmt, 0.0);
75 fGoodAdcPed = vector<Double_t> (MaxNumCerPmt, 0.0);
76 fGoodAdcMult = vector<Double_t> (MaxNumCerPmt, 0.0);
77 fGoodAdcHitUsed = vector<Double_t> (MaxNumCerPmt, 0.0);
78 fGoodAdcPulseInt = vector<Double_t> (MaxNumCerPmt, 0.0);
79 fGoodAdcPulseIntRaw = vector<Double_t> (MaxNumCerPmt, 0.0);
80 fGoodAdcPulseAmp = vector<Double_t> (MaxNumCerPmt, 0.0);
81 fGoodAdcPulseTime = vector<Double_t> (MaxNumCerPmt, 0.0);
82 fGoodAdcTdcDiffTime = vector<Double_t> (MaxNumCerPmt, 0.0);
83
84 InitArrays();
85}
86
87//_____________________________________________________________________________
90{
91 // Constructor
92 frAdcPedRaw = NULL;
93 frAdcPulseIntRaw = NULL;
94 frAdcPulseAmpRaw = NULL;
95 frAdcPulseTimeRaw = NULL;
96 frAdcPed = NULL;
97 frAdcPulseInt = NULL;
98 frAdcPulseAmp = NULL;
99 frAdcPulseTime = NULL;
100 fAdcErrorFlag = NULL;
101 frAdcSampPedRaw = NULL;
105 frAdcSampPed = NULL;
106 frAdcSampPulseInt = NULL;
107 frAdcSampPulseAmp = NULL;
108 frAdcSampPulseTime = NULL;
109
110 InitArrays();
111}
112
113//_____________________________________________________________________________
115{
116 // Destructor
117 delete frAdcPedRaw; frAdcPedRaw = NULL;
118 delete frAdcPulseIntRaw; frAdcPulseIntRaw = NULL;
119 delete frAdcPulseAmpRaw; frAdcPulseAmpRaw = NULL;
121 delete frAdcPed; frAdcPed = NULL;
122 delete frAdcPulseInt; frAdcPulseInt = NULL;
123 delete frAdcPulseAmp; frAdcPulseAmp = NULL;
124 delete frAdcPulseTime; frAdcPulseTime = NULL;
125 delete fAdcErrorFlag; fAdcErrorFlag = NULL;
126 delete frAdcSampPedRaw; frAdcSampPedRaw = NULL;
130 delete frAdcSampPed; frAdcSampPed = NULL;
134
135 DeleteArrays();
136}
137
138//_____________________________________________________________________________
140{
141 fGain = NULL;
142 fPedSum = NULL;
143 fPedSum2 = NULL;
144 fPedLimit = NULL;
145 fPedMean = NULL;
146 fPedCount = NULL;
147 fPed = NULL;
148 fThresh = NULL;
149}
150//_____________________________________________________________________________
152{
153 delete [] fGain; fGain = NULL;
154
155 // 6 Gev variables
156 delete [] fPedSum; fPedSum = NULL;
157 delete [] fPedSum2; fPedSum2 = NULL;
158 delete [] fPedLimit; fPedLimit = NULL;
159 delete [] fPedMean; fPedMean = NULL;
160 delete [] fPedCount; fPedCount = NULL;
161 delete [] fPed; fPed = NULL;
162 delete [] fThresh; fThresh = NULL;
163
164 delete [] fPedDefault; fPedDefault = 0;
167 delete [] fRegionValue; fRegionValue = 0;
168}
169
170//_____________________________________________________________________________
172{
173 // cout << "THcCherenkov::Init for: " << GetName() << endl;
174
175 string EngineDID = string(GetApparatus()->GetName()).substr(0, 1) + GetName();
176 std::transform(EngineDID.begin(), EngineDID.end(), EngineDID.begin(), ::toupper);
177 if(gHcDetectorMap->FillMap(fDetMap, EngineDID.c_str()) < 0) {
178 static const char* const here = "Init()";
179 Error(Here(here), "Error filling detectormap for %s.", EngineDID.c_str());
180 return kInitError;
181 }
182
183 EStatus status;
184 if((status = THaNonTrackingDetector::Init( date )))
185 return fStatus=status;
186
187 // Should probably put this in ReadDatabase as we will know the
188 // maximum number of hits after setting up the detector map
189 InitHitList(fDetMap, "THcCherenkovHit", fDetMap->GetTotNumChan()+1,
190 0, fADC_RefTimeCut);
191
193 if( !app ||
194 !(fglHod = dynamic_cast<THcHodoscope*>(app->GetDetector("hod"))) ) {
195 static const char* const here = "ReadDatabase()";
196 Warning(Here(here),"Hodoscope \"%s\" not found. ","hod");
197 }
198
199 fPresentP = 0;
200 THaVar* vpresent = gHaVars->Find(Form("%s.present",GetApparatus()->GetName()));
201 if(vpresent) {
202 fPresentP = (Bool_t *) vpresent->GetValuePointer();
203 }
204 return fStatus = kOK;
205}
206
207//_____________________________________________________________________________
209{
210 // This function is called by THaDetectorBase::Init() once at the beginning
211 // of the analysis.
212
213 // cout << "THcCherenkov::ReadDatabase for: " << GetName() << endl; // Ahmed
214
215 string prefix = string(GetApparatus()->GetName()).substr(0, 1) + GetName();
216 std::transform(prefix.begin(), prefix.end(), prefix.begin(), ::tolower);
217
218 CreateMissReportParms(prefix.c_str());
219
220 fNRegions = 4; // Defualt if not set in paramter file
221
222 DBRequest list_1[] = {
223 {"_tot_pmts", &fNelem, kInt},
224 {"_num_regions", &fNRegions, kInt},
225 {0}
226 };
227
228 gHcParms->LoadParmValues(list_1, prefix.c_str());
229
230 Bool_t optional = true;
231
232 cout << "Created Cherenkov detector " << GetApparatus()->GetName() << "."
233 << GetName() << " with " << fNelem << " PMTs" << endl;
234
235 // 6 GeV pedestal paramters
236 fPedLimit = new Int_t[fNelem];
237 fGain = new Double_t[fNelem];
238 fPedMean = new Double_t[fNelem];
242 // Region parameters
247
248 DBRequest list[]={
249 {"_ped_limit", fPedLimit, kInt, (UInt_t) fNelem, optional},
250 {"_adc_to_npe", fGain, kDouble, (UInt_t) fNelem},
251 {"_SampThreshold", &fSampThreshold, kDouble,0,1},
252 {"_SampNSA", &fSampNSA, kInt,0,1},
253 {"_SampNSAT", &fSampNSAT, kInt,0,1},
254 {"_SampNSB", &fSampNSB, kInt,0,1},
255 {"_OutputSampWaveform", &fOutputSampWaveform, kInt,0,1},
256 {"_UseSampWaveform", &fUseSampWaveform, kInt,0,1},
257 {"_red_chi2_min", &fRedChi2Min, kDouble},
258 {"_red_chi2_max", &fRedChi2Max, kDouble},
259 {"_beta_min", &fBetaMin, kDouble},
260 {"_beta_max", &fBetaMax, kDouble},
261 {"_enorm_min", &fENormMin, kDouble},
262 {"_enorm_max", &fENormMax, kDouble},
263 {"_dp_min", &fDpMin, kDouble},
264 {"_dp_max", &fDpMax, kDouble},
265 {"_mirror_zpos", &fMirrorZPos, kDouble},
266 {"_npe_thresh", &fNpeThresh, kDouble},
267 {"_debug_adc", &fDebugAdc, kInt, 0, 1},
268 {"_adcTimeWindowMin", fAdcTimeWindowMin, kDouble,(UInt_t) fNelem,1},
269 {"_adcTimeWindowMax", fAdcTimeWindowMax, kDouble, (UInt_t) fNelem,1},
270 {"_PedDefault", fPedDefault, kInt, (UInt_t) fNelem,1},
271 {"_adc_tdc_offset", &fAdcTdcOffset, kDouble, 0, 1},
272 {"_region", &fRegionValue[0], kDouble, (UInt_t) fRegionsValueMax},
273 {"_adcrefcut", &fADC_RefTimeCut, kInt, 0, 1},
274 {0}
275 };
276 fSampThreshold = 5.;
277 fSampNSA = 0; // use value stored in event 125 info
278 fSampNSB = 0; // use value stored in event 125 info
279 fSampNSAT = 2; // default value in THcRawHit::SetF250Params
280 cout << " fSampThreshold 1 = " << fSampThreshold << endl;
281 fOutputSampWaveform = 0; // 0= no output , 1 = output Sample Waveform
282 fUseSampWaveform = 0; // 0= do not use , 1 = use Sample Waveform
283 for (Int_t i=0;i<fNelem;i++) {
284 fAdcTimeWindowMin[i]=-1000.;
285 fAdcTimeWindowMax[i]=1000.;
286 fPedDefault[i]=0;
287 }
288 fDebugAdc = 0; // Set ADC debug parameter to false unless set in parameter file
289 fAdcTdcOffset = 0.0;
290 fADC_RefTimeCut = 0;
291
292 gHcParms->LoadParmValues((DBRequest*)&list, prefix.c_str());
293 cout << " fSampThreshold = " << fSampThreshold << endl;
294
295 // if (fDebugAdc) cout << "Cherenkov ADC Debug Flag Set To TRUE" << endl;
296
297 fIsInit = true;
298
299 // cout << "Track Matching Parameters for: " << GetName() << endl;
300 // for (Int_t iregion = 0; iregion < fNRegions; iregion++) {
301 // cout << "Region = " << iregion + 1 << endl;
302 // for (Int_t ivalue = 0; ivalue < 8; ivalue++)
303 // cout << fRegionValue[GetIndex(iregion, ivalue)] << " ";
304 // cout << endl;
305 // }
306
307 // Create arrays to hold pedestal results
309
310 return kOK;
311}
312
313//_____________________________________________________________________________
315{
316 // Initialize global variables for histogramming and tree
317 // cout << "THcCherenkov::DefineVariables called for: " << GetName() << endl;
318
319 if( mode == kDefine && fIsSetup ) return kOK;
320 fIsSetup = ( mode == kDefine );
321
322 // Register variables in global list
323
324 if (fDebugAdc) {
325 RVarDef vars[] = {
326 {"numAdcHits", "Number of ADC Hits Per PMT", "fNumAdcHits"}, // Cherenkov occupancy
327 {"totNumAdcHits", "Total Number of ADC Hits", "fTotNumAdcHits"}, // Cherenkov multiplicity
328 {"adcPedRaw", "Raw ADC pedestals", "frAdcPedRaw.THcSignalHit.GetData()"},
329 {"adcPulseIntRaw", "Raw ADC pulse integrals", "frAdcPulseIntRaw.THcSignalHit.GetData()"},
330 {"adcPulseAmpRaw", "Raw ADC pulse amplitudes", "frAdcPulseAmpRaw.THcSignalHit.GetData()"},
331 {"adcPulseTimeRaw", "Raw ADC pulse times", "frAdcPulseTimeRaw.THcSignalHit.GetData()"},
332 {"adcPed", "ADC pedestals", "frAdcPed.THcSignalHit.GetData()"},
333 {"adcPulseInt", "ADC pulse integrals", "frAdcPulseInt.THcSignalHit.GetData()"},
334 {"adcPulseAmp", "ADC pulse amplitudes", "frAdcPulseAmp.THcSignalHit.GetData()"},
335 {"adcPulseTime", "ADC pulse times", "frAdcPulseTime.THcSignalHit.GetData()"},
336 {"adcSampPedRaw", "Raw ADCSAMP pedestals", "frAdcSampPedRaw.THcSignalHit.GetData()"},
337 {"adcSampPulseIntRaw", "Raw ADCSAMP pulse integrals", "frAdcSampPulseIntRaw.THcSignalHit.GetData()"},
338 {"adcSampPulseAmpRaw", "Raw ADCSAMP pulse amplitudes", "frAdcSampPulseAmpRaw.THcSignalHit.GetData()"},
339 {"adcSampPulseTimeRaw", "Raw ADCSAMP pulse times", "frAdcSampPulseTimeRaw.THcSignalHit.GetData()"},
340 {"adcSampPed", "ADCSAMP pedestals", "frAdcSampPed.THcSignalHit.GetData()"},
341 {"adcSampPulseInt", "ADCSAMP pulse integrals", "frAdcSampPulseInt.THcSignalHit.GetData()"},
342 {"adcSampPulseAmp", "ADCSAMP pulse amplitudes", "frAdcSampPulseAmp.THcSignalHit.GetData()"},
343 {"adcSampPulseTime", "ADCSAMP pulse times", "frAdcSampPulseTime.THcSignalHit.GetData()"},
344 { 0 }
345 };
346 DefineVarsFromList( vars, mode);
347 } //end debug statement
348
349 if (fOutputSampWaveform==1) {
350 RVarDef vars[] = {
351 {"adcSampWaveform", "FADC Sample Waveform", "fSampWaveform"},
352 { 0 }
353 };
354 DefineVarsFromList( vars, mode);
355 }
356
357 RVarDef vars[] = {
358 {"adcCounter", "ADC counter numbers", "frAdcPulseIntRaw.THcSignalHit.GetPaddleNumber()"},
359 {"adcSampleCounter", "ADC SAMP counter numbers", "frAdcSampPulseIntRaw.THcSignalHit.GetPaddleNumber()"},
360 {"adcErrorFlag", "Error Flag for When FPGA Fails", "fAdcErrorFlag.THcSignalHit.GetData()"},
361
362 {"numGoodAdcHits", "Number of Good ADC Hits Per PMT", "fNumGoodAdcHits"}, // Cherenkov occupancy
363 {"totNumGoodAdcHits", "Total Number of Good ADC Hits", "fTotNumGoodAdcHits"}, // Cherenkov multiplicity
364
365 {"numTracksMatched", "Number of Tracks Matched Per Region", "fNumTracksMatched"},
366 {"numTracksFired", "Number of Tracks that Fired Per Region", "fNumTracksFired"},
367 {"totNumTracksMatched", "Total Number of Tracks Matched Per Region", "fTotNumTracksMatched"},
368 {"totNumTracksFired", "Total Number of Tracks that Fired", "fTotNumTracksFired"},
369
370 {"xAtCer", "Track X at Cherenkov mirror", "fXAtCer"},
371 {"yAtCer", "Track Y at Cherenkov mirror", "fYAtCer"},
372
373 {"npe", "Number of PEs", "fNpe"},
374 {"npeSum", "Total Number of PEs", "fNpeSum"},
375
376 {"goodAdcPed", "Good ADC pedestals", "fGoodAdcPed"},
377 {"goodAdcMult", "Good ADC Multiplicity", "fGoodAdcMult"},
378 {"goodAdcHitUsed", "Good ADC Hit Used", "fGoodAdcHitUsed"},
379 {"goodAdcPulseInt", "Good ADC pulse integrals", "fGoodAdcPulseInt"},
380 {"goodAdcPulseIntRaw", "Good ADC raw pulse integrals", "fGoodAdcPulseIntRaw"},
381 {"goodAdcPulseAmp", "Good ADC pulse amplitudes", "fGoodAdcPulseAmp"},
382 {"goodAdcPulseTime", "Good ADC pulse times", "fGoodAdcPulseTime"},
383 {"goodAdcTdcDiffTime", "Good Hodo Start - ADC pulse times", "fGoodAdcTdcDiffTime"},
384 {"RefTime", "Raw ADC RefTime (chan) ", "fRefTime"}, // Raw reference time
385 { 0 }
386 };
387
388 return DefineVarsFromList(vars, mode);
389}
390//_____________________________________________________________________________
391inline
393{
394 // Clear the hit lists
395 fNhits = 0;
396 fTotNumAdcHits = 0;
400
401 fXAtCer = 0.0;
402 fYAtCer = 0.0;
403
404 fNpeSum = 0.0;
406
411
412 frAdcPed->Clear();
417
422
427
428 for (UInt_t ielem = 0; ielem < fNumAdcHits.size(); ielem++)
429 fNumAdcHits.at(ielem) = 0;
430 for (UInt_t ielem = 0; ielem < fNumGoodAdcHits.size(); ielem++)
431 fNumGoodAdcHits.at(ielem) = 0;
432 for (UInt_t ielem = 0; ielem < fNumTracksMatched.size(); ielem++)
433 fNumTracksMatched.at(ielem) = 0;
434 for (UInt_t ielem = 0; ielem < fNumTracksFired.size(); ielem++)
435 fNumTracksFired.at(ielem) = 0;
436 for (UInt_t ielem = 0; ielem < fGoodAdcPed.size(); ielem++) {
437 fGoodAdcPed.at(ielem) = 0.0;
438 fGoodAdcMult.at(ielem) = 0.0;
439 fGoodAdcHitUsed.at(ielem) = 0.0;
440 fGoodAdcPulseInt.at(ielem) = 0.0;
441 fGoodAdcPulseIntRaw.at(ielem) = 0.0;
442 fGoodAdcPulseAmp.at(ielem) = 0.0;
443 fGoodAdcPulseTime.at(ielem) = kBig;
444 fGoodAdcTdcDiffTime.at(ielem) = kBig;
445 fNpe.at(ielem) = 0.0;
446 }
447
448}
449
450//_____________________________________________________________________________
452{
453 // Get the Hall C style hitlist (fRawHitList) for this event
454 Bool_t present = kTRUE; // Suppress reference time warnings
455 if(fPresentP) { // if this spectrometer not part of trigger
456 present = *fPresentP;
457 }
458 // THcHallCSpectrometer *app = dynamic_cast<THcHallCSpectrometer*>(GetApparatus());
459 // cout << "Cerenkov Event num = " << evdata.GetEvNum() << " spec = " << app->GetName() << " det = " << GetName()<< endl;
460 fNhits = DecodeToHitList(evdata, !present);
461 if(gHaCuts->Result("Pedestal_event")) {
463 fAnalyzePedestals = 1; // Analyze pedestals first normal events
464 return(0);
465 }
466
469 fAnalyzePedestals = 0; // Don't analyze pedestals next event
470 }
471
472 Int_t ihit = 0;
473 UInt_t nrAdcHits = 0;
474 UInt_t nrSampAdcHits = 0;
475 fSampWaveform.clear();
476
477 while(ihit < fNhits) {
478
480 Int_t npmt = hit->fCounter;
481 THcRawAdcHit& rawAdcHit = hit->GetRawAdcHitPos();
482 if ((rawAdcHit.GetNPulses() >0 || rawAdcHit.GetNSamples() >0) && rawAdcHit.HasRefTime()) {
483 fRefTime=rawAdcHit.GetRefTime() ;
484 }
485 if ( fUseSampWaveform == 0 ) {
486 for (UInt_t thit = 0; thit < rawAdcHit.GetNPulses(); thit++) {
487
488 ((THcSignalHit*) frAdcPedRaw->ConstructedAt(nrAdcHits))->Set(npmt, rawAdcHit.GetPedRaw());
489 ((THcSignalHit*) frAdcPed->ConstructedAt(nrAdcHits))->Set(npmt, rawAdcHit.GetPed());
490
491 ((THcSignalHit*) frAdcPulseIntRaw->ConstructedAt(nrAdcHits))->Set(npmt, rawAdcHit.GetPulseIntRaw(thit));
492 ((THcSignalHit*) frAdcPulseInt->ConstructedAt(nrAdcHits))->Set(npmt, rawAdcHit.GetPulseInt(thit));
493
494 ((THcSignalHit*) frAdcPulseAmpRaw->ConstructedAt(nrAdcHits))->Set(npmt, rawAdcHit.GetPulseAmpRaw(thit));
495 ((THcSignalHit*) frAdcPulseAmp->ConstructedAt(nrAdcHits))->Set(npmt, rawAdcHit.GetPulseAmp(thit));
496
497 ((THcSignalHit*) frAdcPulseTimeRaw->ConstructedAt(nrAdcHits))->Set(npmt, rawAdcHit.GetPulseTimeRaw(thit));
498 ((THcSignalHit*) frAdcPulseTime->ConstructedAt(nrAdcHits))->Set(npmt, rawAdcHit.GetPulseTime(thit)+fAdcTdcOffset);
499
500 if (rawAdcHit.GetPulseAmpRaw(thit) > 0) ((THcSignalHit*) fAdcErrorFlag->ConstructedAt(nrAdcHits))->Set(npmt, 0);
501 if (rawAdcHit.GetPulseAmpRaw(thit) <= 0) ((THcSignalHit*) fAdcErrorFlag->ConstructedAt(nrAdcHits))->Set(npmt, 1);
502 if (rawAdcHit.GetPulseAmpRaw(thit) <= 0 && rawAdcHit.GetNSamples() >0) ((THcSignalHit*) fAdcErrorFlag->ConstructedAt(nrAdcHits))->Set(npmt, 2);
503
504 if (rawAdcHit.GetPulseAmpRaw(thit) <= 0) {
505 Double_t PeakPedRatio= rawAdcHit.GetF250_PeakPedestalRatio();
506 Int_t NPedSamples= rawAdcHit.GetF250_NPedestalSamples();
507 Double_t AdcToC = rawAdcHit.GetAdcTopC();
508 Double_t AdcToV = rawAdcHit.GetAdcTomV();
509 if (fPedDefault[npmt-1] !=0) {
510 Double_t tPulseInt = AdcToC*(rawAdcHit.GetPulseIntRaw(thit) - fPedDefault[npmt-1]*PeakPedRatio);
511 ((THcSignalHit*) frAdcPulseInt->ConstructedAt(nrAdcHits))->Set(npmt, tPulseInt);
512 ((THcSignalHit*) frAdcPedRaw->ConstructedAt(nrAdcHits))->Set(npmt, fPedDefault[npmt-1]);
513 ((THcSignalHit*) frAdcPed->ConstructedAt(nrAdcHits))->Set(npmt, float(fPedDefault[npmt-1])/float(NPedSamples)*AdcToV);
514
515 }
516 ((THcSignalHit*) frAdcPulseAmp->ConstructedAt(nrAdcHits))->Set(npmt, 0.);
517
518 }
519
520
521 ++nrAdcHits;
523 fNumAdcHits.at(npmt-1) = npmt;
524 }
525 }
526 //
527 if (rawAdcHit.GetNSamples() >0 ) {
529 if (fSampNSA == 0) fSampNSA=rawAdcHit.GetF250_NSA();
530 if (fSampNSB == 0) fSampNSB=rawAdcHit.GetF250_NSB();
531 rawAdcHit.SetF250Params(fSampNSA,fSampNSB,4); // Set NPED =4
532 if (fSampNSAT != 2) rawAdcHit.SetSampNSAT(fSampNSAT);
533 rawAdcHit.SetSampIntTimePedestalPeak();
534 fSampWaveform.push_back(float(npmt));
535 fSampWaveform.push_back(float(rawAdcHit.GetNSamples()));
536 for (UInt_t thit = 0; thit < rawAdcHit.GetNSamples(); thit++) {
537 fSampWaveform.push_back(rawAdcHit.GetSample(thit)); // ped subtracted sample (mV)
538 }
539 for (UInt_t thit = 0; thit < rawAdcHit.GetNSampPulses(); thit++) {
540 ((THcSignalHit*) frAdcSampPedRaw->ConstructedAt(nrSampAdcHits))->Set(npmt, rawAdcHit.GetSampPedRaw());
541 ((THcSignalHit*) frAdcSampPed->ConstructedAt(nrSampAdcHits))->Set(npmt, rawAdcHit.GetSampPed());
542
543 ((THcSignalHit*) frAdcSampPulseIntRaw->ConstructedAt(nrSampAdcHits))->Set(npmt, rawAdcHit.GetSampPulseIntRaw(thit));
544 ((THcSignalHit*) frAdcSampPulseInt->ConstructedAt(nrSampAdcHits))->Set(npmt, rawAdcHit.GetSampPulseInt(thit));
545
546 ((THcSignalHit*) frAdcSampPulseAmpRaw->ConstructedAt(nrSampAdcHits))->Set(npmt, rawAdcHit.GetSampPulseAmpRaw(thit));
547 ((THcSignalHit*) frAdcSampPulseAmp->ConstructedAt(nrSampAdcHits))->Set(npmt, rawAdcHit.GetSampPulseAmp(thit));
548 ((THcSignalHit*) frAdcSampPulseTimeRaw->ConstructedAt(nrSampAdcHits))->Set(npmt, rawAdcHit.GetSampPulseTimeRaw(thit));
549 ((THcSignalHit*) frAdcSampPulseTime->ConstructedAt(nrSampAdcHits))->Set(npmt, rawAdcHit.GetSampPulseTime(thit)+fAdcTdcOffset);
550 //
551 if ( rawAdcHit.GetNPulses() ==0 || fUseSampWaveform ==1 ) {
552 ((THcSignalHit*) frAdcPedRaw->ConstructedAt(nrAdcHits))->Set(npmt, rawAdcHit.GetSampPedRaw());
553 ((THcSignalHit*) frAdcPed->ConstructedAt(nrAdcHits))->Set(npmt, rawAdcHit.GetSampPed());
554
555 ((THcSignalHit*) frAdcPulseIntRaw->ConstructedAt(nrAdcHits))->Set(npmt,rawAdcHit.GetSampPulseIntRaw(thit));
556 ((THcSignalHit*) frAdcPulseInt->ConstructedAt(nrAdcHits))->Set(npmt,rawAdcHit.GetSampPulseInt(thit));
557
558 ((THcSignalHit*) frAdcPulseAmpRaw->ConstructedAt(nrAdcHits))->Set(npmt,rawAdcHit.GetSampPulseAmpRaw(thit));
559 ((THcSignalHit*) frAdcPulseAmp->ConstructedAt(nrAdcHits))->Set(npmt,rawAdcHit.GetSampPulseAmp(thit) );
560
561 ((THcSignalHit*) frAdcPulseTimeRaw->ConstructedAt(nrAdcHits))->Set(npmt,rawAdcHit.GetSampPulseTimeRaw(thit) );
562 ((THcSignalHit*) frAdcPulseTime->ConstructedAt(nrAdcHits))->Set(npmt, rawAdcHit.GetSampPulseTime(thit)+fAdcTdcOffset);
563 ((THcSignalHit*) fAdcErrorFlag->ConstructedAt(nrAdcHits))->Set(npmt, 3);
564 if (fUseSampWaveform ==1) ((THcSignalHit*) fAdcErrorFlag->ConstructedAt(nrAdcHits))->Set(npmt, 0);
565 ++nrAdcHits;
567 fNumAdcHits.at(npmt-1) = npmt;
568 }
569 ++nrSampAdcHits;
570 }
571 }
572 // Int_t t=0;
573 // std::cin >> t;
574
575 //
576 ihit++;
577 }
578 return ihit;
579}
580
581//_____________________________________________________________________________
583{
584 return(0);
585}
586
587//_____________________________________________________________________________
589{
590 Double_t StartTime = 0.0;
591 if( fglHod ) StartTime = fglHod->GetStartTime();
592 Double_t OffsetTime = 0.0;
593 if( fglHod ) OffsetTime = fglHod->GetOffsetTime();
594 for(Int_t ipmt = 0; ipmt < fNelem; ipmt++) {
595 fAdcPulseAmpTest[ipmt] = -1000.;
596 fAdcGoodElem[ipmt]=-1;
597 }
598 //
599 for(Int_t ielem = 0; ielem < frAdcPulseInt->GetEntries(); ielem++) {
600 Int_t npmt = ((THcSignalHit*) frAdcPulseInt->ConstructedAt(ielem))->GetPaddleNumber() - 1;
601 Double_t pulseTime = ((THcSignalHit*) frAdcPulseTime->ConstructedAt(ielem))->GetData();
602 Double_t pulseAmp = ((THcSignalHit*) frAdcPulseAmp->ConstructedAt(ielem))->GetData();
603 Int_t errorFlag = ((THcSignalHit*) fAdcErrorFlag->ConstructedAt(ielem))->GetData();
604 Double_t adctdcdiffTime = StartTime-pulseTime+OffsetTime;
605 Bool_t pulseTimeCut = adctdcdiffTime > fAdcTimeWindowMin[npmt] && adctdcdiffTime < fAdcTimeWindowMax[npmt];
606 fGoodAdcMult.at(npmt) += 1;
607 if (errorFlag == 0 || errorFlag == 3) {
608 if (pulseTimeCut && pulseAmp > fAdcPulseAmpTest[npmt]) {
609 fAdcGoodElem[npmt]=ielem;
610 fAdcPulseAmpTest[npmt] = pulseAmp;
611 }
612 } else if (errorFlag == 1) {
613 fAdcGoodElem[npmt]=ielem;
614 }
615 }
616 // Loop over the npmt
617 for(Int_t npmt = 0; npmt < fNelem; npmt++) {
618 Int_t ielem = fAdcGoodElem[npmt];
619 if (ielem != -1) {
620 Double_t pulsePed = ((THcSignalHit*) frAdcPed->ConstructedAt(ielem))->GetData();
621 Double_t pulseInt = ((THcSignalHit*) frAdcPulseInt->ConstructedAt(ielem))->GetData();
622 Double_t pulseIntRaw = ((THcSignalHit*) frAdcPulseIntRaw->ConstructedAt(ielem))->GetData();
623 Double_t pulseAmp = ((THcSignalHit*) frAdcPulseAmp->ConstructedAt(ielem))->GetData();
624 Double_t pulseTime = ((THcSignalHit*) frAdcPulseTime->ConstructedAt(ielem))->GetData();
625 Double_t adctdcdiffTime = StartTime-pulseTime+OffsetTime;
626 fGoodAdcPed.at(npmt) = pulsePed;
627 fGoodAdcHitUsed.at(npmt) = ielem+1;
628 fGoodAdcPulseInt.at(npmt) = pulseInt;
629 fGoodAdcPulseIntRaw.at(npmt) = pulseIntRaw;
630 fGoodAdcPulseAmp.at(npmt) = pulseAmp;
631 fGoodAdcPulseTime.at(npmt) = pulseTime;
632 fGoodAdcTdcDiffTime.at(npmt) = adctdcdiffTime;
633
634 fNpe.at(npmt) = fGain[npmt]*pulseInt;
635 fNpeSum += fNpe.at(npmt);
636
638 fNumGoodAdcHits.at(npmt) = npmt + 1;
639 }
640 }
641 return 0;
642}
643
644//_____________________________________________________________________________
646{
647
648 Int_t nTracks = tracks.GetLast() + 1;
649
650 for (Int_t itrack = 0; itrack < nTracks; itrack++) {
651
652 THaTrack* track = dynamic_cast<THaTrack*> (tracks[itrack]);
653 if (track->GetIndex() != 0) continue; // Select the best track
654
655 Double_t trackChi2 = track->GetChi2();
656 Int_t trackNDoF = track->GetNDoF();
657 Double_t trackRedChi2 = trackChi2/trackNDoF;
658 Double_t trackBeta = track->GetBeta();
659 Double_t trackEnergy = track->GetEnergy();
660 Double_t trackMom = track->GetP();
661 Double_t trackDp = track->GetDp();
662 Double_t trackENorm = trackEnergy/trackMom;
663 Double_t trackXfp = track->GetX();
664 Double_t trackYfp = track->GetY();
665 Double_t trackTheta = track->GetTheta();
666 Double_t trackPhi = track->GetPhi();
667
668 Bool_t trackRedChi2Cut = trackRedChi2 > fRedChi2Min && trackRedChi2 < fRedChi2Max;
669 Bool_t trackBetaCut = trackBeta > fBetaMin && trackBeta < fBetaMax;
670 Bool_t trackENormCut = trackENorm > fENormMin && trackENorm < fENormMax;
671 Bool_t trackDpCut = trackDp > fDpMin && trackDp < fDpMax;
672 fXAtCer = trackXfp + trackTheta * fMirrorZPos;
673 fYAtCer = trackYfp + trackPhi * fMirrorZPos;
674
675 if (trackRedChi2Cut && trackBetaCut && trackENormCut && trackDpCut) {
676
677 // Project the track to the Cherenkov mirror planes
678
679 // cout << "Cherenkov Detector: " << GetName() << " has fNRegions = " << fNRegions << endl;
680 // cout << "nTracks = " << nTracks << "\t" << "trackChi2 = " << trackChi2
681 // << "\t" << "trackNDof = " << trackNDoF << "\t" << "trackRedChi2 = " << trackRedChi2 << endl;
682 // cout << "trackBeta = " << trackBeta << "\t" << "trackEnergy = " << trackEnergy << "\t"
683 // << "trackMom = " << trackMom << "\t" << "trackENorm = " << trackENorm << endl;
684 // cout << "trackXfp = " << trackXfp << "\t" << "trackYfp = " << trackYfp << "\t"
685 // << "trackTheta = " << trackTheta << "\t" << "trackPhi = " << trackPhi << endl;
686 // cout << "fMirrorZPos = " << fMirrorZPos << "\t" << "fXAtCer = " << fXAtCer << "\t" << "fYAtCer = " << fYAtCer << endl;
687 // cout << "=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:=:" << endl;
688
689 for (Int_t iregion = 0; iregion < fNRegions; iregion++) {
690
691 if ((TMath::Abs(fRegionValue[GetIndex(iregion, 0)] - fXAtCer) < fRegionValue[GetIndex(iregion, 4)]) &&
692 (TMath::Abs(fRegionValue[GetIndex(iregion, 1)] - fYAtCer) < fRegionValue[GetIndex(iregion, 5)]) &&
693 (TMath::Abs(fRegionValue[GetIndex(iregion, 2)] - trackTheta) < fRegionValue[GetIndex(iregion, 6)]) &&
694 (TMath::Abs(fRegionValue[GetIndex(iregion, 3)] - trackPhi) < fRegionValue[GetIndex(iregion, 7)])) {
695
697 fNumTracksMatched.at(iregion) = iregion + 1;
698
699 if (fNpeSum > fNpeThresh) {
701 fNumTracksFired.at(iregion) = iregion + 1;
702 } // NPE threshold cut
703 } // Regional cuts
704 } // Loop over regions
705 } // Tracking cuts
706 } // Track loop
707
708 return 0;
709}
710
711//_____________________________________________________________________________
713{
715 fMinPeds = 500; // In engine, this is set in parameter file
716 fPedSum = new Int_t [fNelem];
717 fPedSum2 = new Int_t [fNelem];
718 fPedCount = new Int_t [fNelem];
719 fPed = new Double_t [fNelem];
720 fThresh = new Double_t [fNelem];
721 for(Int_t i = 0; i < fNelem; i++) {
722 fPedSum[i] = 0;
723 fPedSum2[i] = 0;
724 fPedCount[i] = 0;
725 }
726}
727
728//_____________________________________________________________________________
730{
731 // Extract data from the hit list, accumulating into arrays for
732 // calculating pedestals
733
734 Int_t nrawhits = rawhits->GetLast()+1;
735
736 Int_t ihit = 0;
737 while(ihit < nrawhits) {
738 THcCherenkovHit* hit = (THcCherenkovHit *) rawhits->At(ihit);
739
740 Int_t element = hit->fCounter - 1;
741 Int_t nadc = hit->GetRawAdcHitPos().GetPulseIntRaw();
742 if(nadc <= fPedLimit[element]) {
743 fPedSum[element] += nadc;
744 fPedSum2[element] += nadc*nadc;
745 fPedCount[element]++;
746 if(fPedCount[element] == fMinPeds/5) {
747 fPedLimit[element] = 100 + fPedSum[element]/fPedCount[element];
748 }
749 }
750 ihit++;
751 }
752
754
755 return;
756}
757
758//_____________________________________________________________________________
760{
761 // Use the accumulated pedestal data to calculate pedestals
762 // Later add check to see if pedestals have drifted ("Danger Will Robinson!")
763 // cout << "Plane: " << fPlaneNum << endl;
764 for(Int_t i=0; i<fNelem;i++) {
765
766 // PMT tubes
767 fPed[i] = ((Double_t) fPedSum[i]) / TMath::Max(1, fPedCount[i]);
768 fThresh[i] = fPed[i] + 15;
769
770 // Just a copy for now, but allow the possibility that fXXXPedMean is set
771 // in a parameter file and only overwritten if there is a sufficient number of
772 // pedestal events. (So that pedestals are sensible even if the pedestal events were
773 // not acquired.)
774 if(fMinPeds > 0) {
775 if(fPedCount[i] > fMinPeds) {
776 fPedMean[i] = fPed[i];
777 }
778 }
779 }
780 // cout << " " << endl;
781
782}
783
784//_____________________________________________________________________________
786{
787 return fNRegions * nValue + nRegion;
788}
789
790
791//_____________________________________________________________________________
792void THcCherenkov::Print(const Option_t* opt) const
793{
795 // Print out the pedestals
796 cout << endl;
797 cout << "Cherenkov Pedestals" << endl;
798 // Ahmed
799 cout << "No. ADC" << endl;
800 for(Int_t i=0; i<fNelem; i++)
801 cout << " " << i << " " << fPed[i] << endl;
802 cout << endl;
803}
804
805//_____________________________________________________________________________
810//_____________________________________________________________________________
812{
813 MissReport(Form("%s.%s", GetApparatus()->GetName(), GetName()));
814 return 0;
815}
int Int_t
unsigned int UInt_t
const Data_t kBig
bool Bool_t
double Double_t
const Bool_t kTRUE
const char Option_t
Option_t Option_t TPoint TPoint const char mode
R__EXTERN class THaVarList * gHaVars
R__EXTERN class THaCutList * gHaCuts
ClassImp(VDC::AnalyticTTDConv) using namespace std
R__EXTERN class THcParmList * gHcParms
Definition THcGlobals.h:11
R__EXTERN class THcDetectorMap * gHcDetectorMap
Definition THcGlobals.h:12
char * Form(const char *fmt,...)
void Clear(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="")
virtual const char * Here(const char *) const
virtual void Print(Option_t *opt="") const
virtual Int_t Result(const char *cutname="", EWarnMode mode=kWarn)
UInt_t GetTotNumChan() const
THaDetMap * fDetMap
THaApparatus * GetApparatus() const
Double_t GetChi2() const
Int_t GetNDoF() const
Double_t GetX() const
Double_t GetPhi() const
Double_t GetEnergy() const
Double_t GetDp() const
Double_t GetTheta() const
Double_t GetBeta() const
Double_t GetY() const
Int_t GetIndex() const
Double_t GetP() const
virtual THaVar * Find(const char *name) const
const void * GetValuePointer() const
Class representing a Cherenkov PMT hit.
Class for gas Cherenkov detectors.
Double_t fENormMax
Double_t fRedChi2Min
virtual void CalculatePedestals()
vector< Double_t > fGoodAdcPulseTime
Double_t fAdcTdcOffset
TClonesArray * frAdcPulseTimeRaw
Int_t fNPedestalEvents
Double_t fRefTime
Double_t GetCerNPE()
vector< Double_t > fGoodAdcPulseIntRaw
vector< Int_t > fNumGoodAdcHits
TClonesArray * frAdcSampPulseInt
Double_t * fPed
Int_t * fPedSum
Double_t fBetaMax
Double_t * fThresh
virtual void Print(const Option_t *opt) const
Int_t fAnalyzePedestals
Double_t fYAtCer
virtual Int_t Decode(const THaEvData &)
Int_t fTotNumTracksMatched
Double_t * fRegionValue
virtual void Clear(Option_t *opt="")
Int_t * fPedCount
virtual Int_t ApplyCorrections(void)
Int_t fTotNumTracksFired
Double_t fXAtCer
virtual void InitializePedestals()
Double_t fNpeThresh
Int_t End(THaRunBase *run)
virtual Int_t FineProcess(TClonesArray &tracks)
Int_t fTotNumAdcHits
TClonesArray * frAdcPulseAmp
Double_t fDpMax
vector< Double_t > fGoodAdcTdcDiffTime
virtual Int_t CoarseProcess(TClonesArray &tracks)
TClonesArray * frAdcSampPedRaw
vector< Double_t > fSampWaveform
virtual void AccumulatePedestals(TClonesArray *rawhits)
Int_t fADC_RefTimeCut
Int_t GetIndex(Int_t nRegion, Int_t nValue)
vector< Int_t > fNumAdcHits
vector< Double_t > fGoodAdcMult
vector< Double_t > fGoodAdcPulseInt
Double_t fENormMin
TClonesArray * frAdcSampPulseTime
Int_t * fPedLimit
TClonesArray * frAdcPedRaw
TClonesArray * frAdcPulseIntRaw
Double_t fDpMin
Double_t * fAdcPulseAmpTest
vector< Double_t > fNpe
Int_t fUseSampWaveform
static const Int_t MaxNumCerPmt
TClonesArray * frAdcPulseTime
TClonesArray * frAdcSampPulseAmpRaw
TClonesArray * frAdcPed
Double_t * fGain
vector< Int_t > fNumTracksMatched
Int_t fTotNumGoodAdcHits
Bool_t * fPresentP
Int_t * fPedSum2
Double_t fSampThreshold
virtual ~THcCherenkov()
Double_t * fAdcTimeWindowMax
Double_t * fPedMean
Double_t fRedChi2Max
TClonesArray * frAdcSampPulseIntRaw
TClonesArray * frAdcSampPulseTimeRaw
Int_t fRegionsValueMax
Double_t fMirrorZPos
Double_t * fAdcTimeWindowMin
vector< Double_t > fGoodAdcHitUsed
THcHodoscope * fglHod
TClonesArray * frAdcSampPulseAmp
static const Int_t MaxNumAdcPulse
Int_t * fPedDefault
Int_t * fAdcGoodElem
Double_t fNpeSum
TClonesArray * fAdcErrorFlag
TClonesArray * frAdcSampPed
virtual Int_t ReadDatabase(const TDatime &date)
vector< Double_t > fGoodAdcPed
vector< Int_t > fNumTracksFired
vector< Double_t > fGoodAdcPulseAmp
Double_t fBetaMin
TClonesArray * frAdcPulseInt
virtual Int_t DefineVariables(EMode mode=kDefine)
TClonesArray * frAdcPulseAmpRaw
Int_t fOutputSampWaveform
virtual Int_t FillMap(THaDetMap *detmap, const char *detectorname)
A standard Hall C spectrometer apparatus.
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.
void MissReport(const char *name) const
Generic hodoscope consisting of multiple planes with multiple paddles with phototubes on both ends.
Double_t GetStartTime() const
Double_t GetOffsetTime() const
Int_t LoadParmValues(const DBRequest *list, const char *prefix="")
Retrieve parameter values from the parameter cache.
Class representing a single raw ADC hit.
Definition THcRawAdcHit.h:7
Int_t GetSampPulseIntRaw(UInt_t iPulse=0) const
UInt_t GetNSampPulses() const
Double_t GetSampPulseInt(UInt_t iPulse=0) const
Int_t GetPedRaw() const
Gets raw signal pedestal. In channels.
Int_t GetSampPulseAmpRaw(UInt_t iPulse=0) const
Double_t GetSampPulseTime(UInt_t iPulse=0) const
UInt_t GetNPulses() const
Gets number of set pulses.
Int_t GetF250_NSB() const
UInt_t GetNSamples() const
Gets number of set samples.
Int_t GetSampPedRaw() const
Bool_t HasRefTime() const
void SetSampNSAT(Int_t nsat)
static constexpr Double_t GetAdcTopC()
void SetSampThreshold(Double_t thres)
Double_t GetPulseAmp(UInt_t iPulse=0) const
Gets pedestal subtracted pulse amplitude. In channels.
Int_t GetRefTime() const
Double_t GetSample(UInt_t iSample=0) const
void SetSampIntTimePedestalPeak()
Int_t GetPulseAmpRaw(UInt_t iPulse=0) const
Gets raw pulse amplitude. In channels.
Double_t GetPed() const
Gets sample pedestal. In channels.
Double_t GetF250_PeakPedestalRatio() const
Int_t GetPulseTimeRaw(UInt_t iPulse=0) const
Gets raw pulse time. In subsamples.
Double_t GetSampPulseAmp(UInt_t iPulse=0) const
Double_t GetPulseTime(UInt_t iPulse=0) const
Int_t GetF250_NPedestalSamples() const
Int_t GetPulseIntRaw(UInt_t iPulse=0) const
Gets raw pulse integral. In channels.
Int_t GetSampPulseTimeRaw(UInt_t iPulse=0) const
Double_t GetSampPed() const
Int_t GetF250_NSA() const
static constexpr Double_t GetAdcTomV()
void SetF250Params(Int_t NSA, Int_t NSB, Int_t NPED)
Sets F250 parameters used for pedestal subtraction.
Double_t GetPulseInt(UInt_t iPulse=0) const
Gets pedestal subtracted pulse integral. In channels.
Int_t fCounter
Definition THcRawHit.h:53
THcRawAdcHit & GetRawAdcHitPos()
const char * GetName() const override
Int_t GetEntries() const override
TObject * At(Int_t idx) const override
Int_t GetLast() const override
virtual void Warning(const char *method, const char *msgfmt,...) const
virtual void Error(const char *method, const char *msgfmt,...) const
Double_t Abs(Double_t d)
Double_t Max(Double_t a, Double_t b)
STL namespace.
void tracks()