StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StTpcRTSHitMaker.cxx
1 /***************************************************************************
2  *
3  * $Id: StTpcRTSHitMaker.cxx,v 1.59 2021/05/10 21:13:19 fisyak Exp $
4  *
5  * Author: Valeri Fine, BNL Feb 2007
6  ***************************************************************************
7  *
8  * Description: Make clusters from StTpcRawData and fill the StEvent */
9 #include <assert.h>
10 #include <stdio.h>
11 #include "StTpcHitMaker.h"
12 #include "StTpcRTSHitMaker.h"
13 
14 #include "StTpcRawData.h"
15 #include "StEvent/StTpcRawData.h"
16 #include "StEvent/StTpcHit.h"
17 #include "StEvent/StEvent.h"
18 #include "StEvent/StTpcHitCollection.h"
19 #include "StThreeVectorF.hh"
20 #include "StTpcDb/StTpcDb.h"
21 #include "StDbUtilities/StCoordinates.hh"
22 #include "StDetectorDbMaker/St_tss_tssparC.h"
23 #include "StDetectorDbMaker/St_tpcPadGainT0BC.h"
24 #include "StDetectorDbMaker/St_tpcAnodeHVavgC.h"
25 #include "StDetectorDbMaker/St_tpcMaxHitsC.h"
26 #include "StDetectorDbMaker/StDetectorDbTpcRDOMasks.h"
27 #include "StDetectorDbMaker/St_tpcPadPlanesC.h"
28 #include "StDetectorDbMaker/St_itpcPadPlanesC.h"
29 #include "StDetectorDbMaker/St_tpcPadConfigC.h"
30 #include "StDetectorDbMaker/St_tpcStatusC.h"
31 #include "StDetectorDbMaker/StPath2tpxGain.h"
32 #include "StDetectorDbMaker/StPath2itpcGain.h"
33 #include "StDetectorDbMaker/St_tpcPadConfigC.h"
34 #include "StMessMgr.h"
35 #include "StDAQMaker/StDAQReader.h"
36 #include "StRtsTable.h"
37 #include "DAQ_TPX/daq_tpx.h"
38 #include "DAQ_ITPC/daq_itpc.h"
39 #include "DAQ_READER/daq_dta.h"
40 #include "DAQ_READER/daq_det.h"
41 #include "DAQ_READER/daqReader.h"
42 #include "RTS/src/DAQ_TPX/tpxFCF_flags.h" // for FCF flag definition
43 #include "RTS/src/DAQ_TPX/daq_tpx.h"
44 #include "RTS/src/DAQ_TPX/tpxCore.h"
45 #include "RTS/src/DAQ_TPX/tpxPed.h"
46 #include "RTS/src/DAQ_TPX/tpxGain.h"
47 #include "RTS/src/DAQ_TPX/tpxFCF.h"
48 #include "RTS/src/DAQ_TPX/tpxStat.h"
49 #include "RTS/src/DAQ_ITPC/itpcFCF.h"
50 #include "RTS/src/DAQ_TPC23/tpc23_base.h"
51 #include "RTS/src/DAQ_TPC23/tpx23.h"
52 #include "RTS/src/DAQ_TPC23/itpc23.h"
53 #include "TBenchmark.h"
54 #include "TH2.h"
55 #include "TFile.h"
56 #include "TPolyMarker.h"
57 #include "TBox.h"
58 ClassImp(StTpcRTSHitMaker);
59 #define __DEBUG__
60 #ifdef __DEBUG__
61 #define PrPP(A,B) if (Debug()%10 > 1) {LOG_INFO << "StTpcRTSHitMaker::" << (#A) << "\t" << (#B) << " = \t" << (B) << endm;}
62 #else
63 #define PrPP(A,B)
64 #endif
65 StTpcRTSHitMaker *StTpcRTSHitMaker::fgStTpcRTSHitMaker = 0;
66 //________________________________________________________________________________
67 StTpcRTSHitMaker::~StTpcRTSHitMaker() {
68  SafeDelete(fTpx);
69  SafeDelete(fiTpc);
70  SafeDelete(fTpx23);
71  SafeDelete(fiTpc23);
72 }
73 //________________________________________________________________________________
74 Int_t StTpcRTSHitMaker::Init() {
75  memset(maxHits,0,sizeof(maxHits));
76  maxBin0Hits = 0;
77  bin0Hits = 0;
78  return StMaker::Init();
79 }
80 //________________________________________________________________________________
81 Int_t StTpcRTSHitMaker::from_file(daq_dta *gain_dta, const Char_t *fname) {
82  // example of gains; will use file for that
83  FILE *f = fopen(fname,"r") ;
84  if(f==0) {
85  LOG_ERROR << "Can't open gain file\t" << fname << endm;
86  return -1 ;
87  }
88  static struct g_s_t {
89  float g ;
90  float t ;
91  } g_s[25][41][121] ;
92 
93  Int_t bad_ch = 0 ;
94  Int_t all_ch = 0 ;
95  while(!feof(f)) {
96  char buff[128] ;
97  Int_t sec,rdo,port,ch,row,pad ;
98  float g, t ;
99 
100  if(fgets(buff,sizeof(buff),f)==0) continue ;
101 
102  if(buff[0]=='#') continue ;
103  if(strlen(buff)<1) continue ;
104 
105  Int_t ret = sscanf(buff,"%d %d %d %d %d %d %f %f",&sec,&rdo,&port,&ch,&row,&pad,&g,&t) ;
106  if(ret != 8) continue ;
107 
108  if(g<0.01) bad_ch++ ;
109  all_ch++ ;
110 
111  g_s[sec][row][pad].g = g ;
112  g_s[sec][row][pad].t = t ;
113  }
114  LOG_INFO << Form("From gain file %s: %d/%d bad channels",fname,bad_ch,all_ch) << endm;
115  fclose(f) ;
116  // and now load them up
117  for(Int_t s=1;s<=24;s++) { // sectors loop
118  for(Int_t r=1;r<=40;r++) { // rows loop
119  daq_det_gain *gain = (daq_det_gain *) gain_dta->request(121) ;
120  for(Int_t p=0;p<=120;p++) { // pad loop
121  gain[p].gain = g_s[s][r][p].g ;
122  gain[p].t0 = g_s[s][r][p].t ;
123  }
124  gain_dta->finalize(121,s,r) ;
125  }
126  }
127  return 1;
128 }
129 //________________________________________________________________________________
130 Int_t StTpcRTSHitMaker::FixGains(tpc23_base *tpc23, Int_t sec, Int_t row, Int_t Npads) {
131  // from tpc23_base::gains_from_cache, fix dead pads and edges
132  Int_t ok = 0;
133  tpc23_base::row_pad_t (*rp_gain)[tpc23_base::ROW_MAX+1][tpc23_base::PAD_MAX+1] = tpc23->rp_gain;
134  for (Int_t pad = 1; pad <= Npads; pad++) {
135  Float_t g = rp_gain[sec-1][row][pad].gain;
136  if(g<0.01) {
137  int p1 = pad - 1 ;
138  int p2 = pad + 1 ;
139  if(p1<1) p1 = 1 ;
140  if(p2>Npads) p2 = Npads ;
141  rp_gain[sec-1][row][pad].flags |= FCF_DEAD_EDGE ;
142  rp_gain[sec-1][row][p1].flags |= FCF_DEAD_EDGE ;
143  rp_gain[sec-1][row][p2].flags |= FCF_DEAD_EDGE ;
144  }
145  }
146  rp_gain[sec-1][row][1].flags |= FCF_ROW_EDGE ; // row edge
147  rp_gain[sec-1][row][Npads].flags |= FCF_ROW_EDGE ; // row edge
148  return ok;
149 }
150 //________________________________________________________________________________
151 Int_t StTpcRTSHitMaker::InitRun(Int_t runnumber) {
152  SetAttr("minSector",1);
153  SetAttr("maxSector",24);
154  SetAttr("minRow",1);
155  if (IAttr("Cosmics")) StTpcHitMaker::SetCosmics();
156  SetAttr("maxRow",St_tpcPadConfigC::instance()->numberOfRows(20));
157  SafeDelete(fTpx);
158  SafeDelete(fiTpc);
159  SafeDelete(fTpx23);
160  SafeDelete(fiTpc23);
161  TString fnameTPX("none");
162  TString fnameITPC("none");
163  if (IAttr("USE_GAIN_FROM_FILE")) {
164  fnameTPX = StPath2tpxGain::instance()->GetPath();
165  LOG_INFO << "StTpcRTSHitMaker::InitRun: use " << fnameTPX.Data() << " for TPX" << endm;
166  fnameITPC = StPath2itpcGain::instance()->GetPath();
167  LOG_INFO << "StTpcRTSHitMaker::InitRun: use " << fnameITPC.Data() << " for iTPC" << endm;
168  }
169  if ( IAttr("TPC23")) { // TPC23
170  LOG_INFO << "StTpcRTSHitMaker::InitRun:: use TPC23 mode" << endm;
171  Int_t log_level = 0 ;
172  // Int_t fmt_version = 0 ;
173  Int_t run_type = 3;
174  Int_t online = 0 ; // NOTE setting for Offline mode!
175  fTpx23 = new tpx23;
176  fTpx23->online = online ;
177  fTpx23->run_type = run_type ;
178  fTpx23->log_level = log_level ;
179  if (fTpx23->gains_from_cache(fnameTPX) < 0 || fnameTPX == "none") { // REQUIRED even if no gain correction
180  // Tpx Load gains from Db
181  Int_t NoPadsWithLoadedGain = 0;
182  for(Int_t sector=1;sector<=24;sector++) {
183  Int_t rowMin = 1;
184  if (St_tpcPadConfigC::instance()->iTPC(sector)) rowMin = 14;
185  for(Int_t rowO = 1; rowO < rowMin; rowO++) {
186  Int_t Npads = St_tpcPadPlanesC::instance()->padsPerRow(rowO);
187  for(Int_t pad = 0; pad <= Npads; pad++) {
188  fTpx23->rp_gain[sector-1][rowO][pad].gain = 0.; // be sure that dead pads are killed
189  fTpx23->rp_gain[sector-1][rowO][pad].t0 = -9.99;
190  fTpx23->rp_gain[sector-1][rowO][pad].flags = 64;
191  }
192  }
193  for(Int_t rowO = rowMin; rowO <= 45; rowO++) {
194  Int_t Npads = St_tpcPadPlanesC::instance()->padsPerRow(rowO);
195  Int_t padMin = 1;
196  Int_t padMax = Npads;
197  Int_t NoPadsWithLoadedGainPerRow = 0;
198  for(Int_t pad = 0; pad <= Npads; pad++) {
199  fTpx23->rp_gain[sector-1][rowO][pad].gain = 0.; // be sure that dead pads are killed
200  fTpx23->rp_gain[sector-1][rowO][pad].t0 = 0.;
201  if (rowO >= rowMin && pad >= padMin && pad <= padMax) {
202  if (St_tpcPadGainT0C::instance()->Gain(sector,rowO,pad) <= 0) continue;
203  fTpx23->rp_gain[sector-1][rowO][pad].gain = St_tpcPadGainT0C::instance()->Gain(sector,rowO,pad);
204  fTpx23->rp_gain[sector-1][rowO][pad].t0 = St_tpcPadGainT0C::instance()->T0(sector,rowO,pad);
205  NoPadsWithLoadedGainPerRow++;
206  }
207  }
208  if (NoPadsWithLoadedGainPerRow > 0) {
209  NoPadsWithLoadedGain += NoPadsWithLoadedGainPerRow;
210  FixGains(fTpx23, sector, rowO, Npads);
211  }
212  }
213  }
214  LOG_INFO << "StTpcRTSHitMaker::InitRun:: loaded gains for " << NoPadsWithLoadedGain << " pads in TPX23" << endm;
215 
216  }
217  //#define __DEBUG_GAIN__
218 #ifdef __DEBUG_GAIN__
219  LOG_INFO << "StTpcRTSHitMaker::InitRun:: Print gains for Tpx23" << endm;
220  for(Int_t sector=1;sector<=24;sector++) {
221  for(Int_t rowO = 1; rowO <= 45; rowO++) {
222  Int_t Npads = St_tpcPadPlanesC::instance()->padsPerRow(rowO);
223  for(Int_t pad = 1; pad <= Npads; pad++) {
224  cout << Form("Gain/T0 s/r/p %3i/%3i/%3i %7.3f %7.3f %i",sector,rowO,pad,fTpx23->rp_gain[sector-1][rowO][pad].gain,fTpx23->rp_gain[sector-1][rowO][pad].t0,fTpx23->rp_gain[sector-1][rowO][pad].flags) << endl;
225  }
226  }
227  }
228 #endif /* __DEBUG_GAIN__ */
229  fTpx23->run_start() ;
230  // iTPC23
231  fiTpc23 = new itpc23 ;
232  fiTpc23->online = online ;
233  fiTpc23->run_type = run_type ;
234  fiTpc23->log_level = log_level ;
235  Int_t ok = fiTpc23->gains_from_cache(fnameITPC);
236  if (ok < 0 || fnameITPC == "none") { // REQUIRED even if no gain correction
237  Int_t NoPadsWithLoadedGain = 0;
238  for(Int_t sector=1;sector<=24;sector++) {
239  if (! St_tpcPadConfigC::instance()->iTPC(sector)) continue;
240  for(Int_t row = 1; row <= 40; row++) {
241  Int_t Npads = St_itpcPadPlanesC::instance()->padsPerRow(row);
242  Int_t NoPadsWithLoadedGainPerRow = 0;
243  for(Int_t pad = 0; pad <= Npads; pad++) {
244  fiTpc23->rp_gain[sector-1][row][pad].gain = 0.; // be sure that dead pads are killed
245  fiTpc23->rp_gain[sector-1][row][pad].t0 = 0.;
246  if (pad < 1) continue; // kill pad0 just in case..
247  if (St_itpcPadGainT0C::instance()->Gain(sector,row,pad) <= 0) continue;
248  fiTpc23->rp_gain[sector-1][row][pad].gain = St_itpcPadGainT0C::instance()->Gain(sector,row,pad);
249  fiTpc23->rp_gain[sector-1][row][pad].t0 = St_itpcPadGainT0C::instance()->T0(sector,row,pad);
250  //#define __DEBUG_GAIN__
251  NoPadsWithLoadedGainPerRow++;
252  }
253  if (NoPadsWithLoadedGainPerRow > 0) {
254  NoPadsWithLoadedGain += NoPadsWithLoadedGainPerRow;
255  FixGains(fiTpc23, sector, row, Npads);
256  }
257  }
258  }
259  LOG_INFO << "StTpcRTSHitMaker::InitRun:: loaded gains for " << NoPadsWithLoadedGain << " pads in iTPC23" << endm;
260  }
261 #ifdef __DEBUG_GAIN__
262  LOG_INFO << "StTpcRTSHitMaker::InitRun:: Print gains for iTpc23" << endm;
263  for(Int_t sector=1;sector<=24;sector++) {
264  for(Int_t row = 1; row <= 40; row++) {
265  Int_t Npads = St_itpcPadPlanesC::instance()->padsPerRow(row);
266  for(Int_t pad = 1; pad <= Npads; pad++) {
267  cout << Form("Gain/T0 s/r/p %3i/%3i/%3i %7.3f %7.3f %i",sector,row,pad,fiTpc23->rp_gain[sector-1][row][pad].gain,fiTpc23->rp_gain[sector-1][row][pad].t0,fiTpc23->rp_gain[sector-1][row][pad].flags) << endl;
268  }
269  }
270  }
271 #endif /* __DEBUG_GAIN__ */
272  fiTpc23->run_start() ;
273  } else {
274  fTpx = new daq_tpx() ;
275  if (GetDate() >= 20091215) fTpx->fcf_run_compatibility = 10 ;
276  if (GetDate() >= 20191215) fTpx->fcf_run_compatibility = 22 ;
277  // change default value 2 to
278  fTpx->fcf_do_cuts = 1; // 1 means always, 2 means don't cut edges (for i.e. pulser run), 0 means don't...
279  // if (GetDate() >= 20121215) fTpx->fcf_style = 2 ; // from online/RTS/src/ESB/tpx.C new for FY13!
280  if (GetDate() <= 20090101) fminCharge = 40;
281  // do gains example; one loads them from database but I don't know how...
282  // daq_dta *dta_Tpx = fTpx->put("gain");
283  tpxGain *gain = fTpx->gain_algo;
284  // gain->do_default(sector) ; // zap to all 1...
285  Int_t ok = gain->from_file(fnameTPX, 0);
286  if (ok < 0 &&fnameTPX == "none") {
287  for(Int_t sector=1;sector<=24;sector++) {
288  // Tpx
289  Int_t rowMin = 1;
290  if (St_tpcPadConfigC::instance()->iTPC(sector)) rowMin = 14;
291  for(Int_t rowO = 1; rowO <= 45; rowO++) {
292  Int_t Npads = St_tpcPadPlanesC::instance()->padsPerRow(rowO);
293  Int_t padMin = 1;
294  Int_t padMax = Npads;
295  // daq_det_gain *gain = (daq_det_gain *) dta_Tpx->request(Npads+1); // max pad+1
296  for(Int_t pad = 0; pad <= Npads; pad++) {
297  // gain[pad].gain = 0.; // be sure that dead pads are killed
298  // gain[pad].t0 = 0.;
299  gain->set_gains(sector, rowO, pad, 0, 0);
300  if (rowO >= rowMin && pad >= padMin && pad <= padMax) {
301  if (St_tpcPadGainT0C::instance()->Gain(sector,rowO,pad) <= 0) continue;
302  // gain[pad].gain = St_tpcPadGainT0C::instance()->Gain(sector,rowO,pad);
303  // gain[pad].t0 = St_tpcPadGainT0C::instance()->T0(sector,rowO,pad);
304  gain->set_gains(sector, rowO, pad, St_tpcPadGainT0C::instance()->Gain(sector,rowO,pad), St_tpcPadGainT0C::instance()->T0(sector,rowO,pad));
305  }
306  }
307  // daq_dta::finalize(uint32_t obj_cou, Int_t sec, Int_t row, Int_t pad)
308  // dta_Tpx->finalize(Npads+1,sector,rowO);
309  }
310  }
311  }
312  // ((daq_tpx*) dta_Tpx)->InitRun(runnumber);
313  // Check presence iTPC
314  static Bool_t fNoiTPCLu = IAttr("NoiTPCLu");
315  if (! fNoiTPCLu) {
316  for(Int_t sector=1;sector<=24;sector++) {
317  if (St_tpcPadConfigC::instance()->iTPC(sector)) {
318  fiTpc = new daq_itpc() ;
319  break;
320  }
321  }
322  }
323  // iTpc
324  if (fiTpc) {
325  daq_dta * dta_iTpc = fiTpc->put("gain"); // , 0, 40, 0, miTpc_RowLen);
326  Int_t ok = from_file(dta_iTpc, fnameITPC);
327  if (ok < 0 || fnameITPC == "none") {
328  for(Int_t sector=1;sector<=24;sector++) {
329  if (! St_tpcPadConfigC::instance()->iTPC(sector)) continue;
330  for(Int_t row = 1; row <= 40; row++) {
331  Int_t Npads = St_itpcPadPlanesC::instance()->padsPerRow(row);
332  daq_det_gain *gain = (daq_det_gain *) dta_iTpc->request(Npads+1); // max pad+1
333  for(Int_t pad = 0; pad <= Npads; pad++) {
334  gain[pad].gain = 0.; // be sure that dead pads are killed
335  gain[pad].t0 = 0.;
336  if (pad < 1) continue; // kill pad0 just in case..
337  if (St_itpcPadGainT0C::instance()->Gain(sector,row,pad) <= 0) continue;
338  gain[pad].gain = St_itpcPadGainT0C::instance()->Gain(sector,row,pad);
339  gain[pad].t0 = St_itpcPadGainT0C::instance()->T0(sector,row,pad);
340  //#define __DEBUG_GAIN__
341 #ifdef __DEBUG_GAIN__
342  cout << Form("Gain/T0 s/r/p %3i/%3i/%3i %7.2f %7.2f",sector,row,pad,gain[pad].gain,gain[pad].t0) << endl;
343 #endif /* __DEBUG_GAIN__ */
344  }
345  // daq_dta::finalize(uint32_t obj_cou, Int_t sec, Int_t row, Int_t pad)
346  dta_iTpc->finalize(Npads+1,sector,row);
347  }
348  // ((daq_itpc *) dta_iTpc)->InitRun(runnumber);
349  }
350  }
351  } //
352  }
353  PrintAttr();
354 
356  // XD: restore hit occupancy protection code to avoid long-hung of production jobs
358  // Prepare scaled hit maxima
359  // No hit maxima if these DB params are 0
360  Int_t maxHitsPerSector = St_tpcMaxHitsC::instance()->maxSectorHits();
361  Int_t maxBinZeroHits = St_tpcMaxHitsC::instance()->maxBinZeroHits();
362  if ( !(maxHitsPerSector > 0 || maxBinZeroHits > 0) ) return kStOK;
363  Int_t livePads = 0;
364  Int_t totalPads = 0;
365  Float_t liveFrac = 1;
366  for(Int_t sector=1;sector<=24;sector++) {
367  Int_t liveSecPads = 0;
368  Int_t totalSecPads = 0;
369 
370  // Tpx
371  Int_t rowMin = 1;
372  if (St_tpcPadConfigC::instance()->iTPC(sector)) rowMin = 14;
373  for(Int_t row = rowMin; row <= 45; row++) {
374  Int_t Npads = St_tpcPadPlanesC::instance()->padsPerRow(row);
375  totalSecPads += Npads;
376  if (StDetectorDbTpcRDOMasks::instance()->isRowOn(sector,row) &&
377  St_tpcAnodeHVavgC::instance()->livePadrow(sector,row))
378  liveSecPads += Npads;
379  }
380 
381  //iTPC
382  if (St_tpcPadConfigC::instance()->iTPC(sector)) {
383  for(Int_t row = 1; row <= 40; row++) {
384  Int_t Npads = St_itpcPadPlanesC::instance()->padsPerRow(row);
385  totalSecPads += Npads;
386  if (StDetectorDbTpcRDOMasks::instance()->isRowOn(sector,row) &&
387  St_tpcAnodeHVavgC::instance()->livePadrow(sector,row))
388  liveSecPads += Npads;
389  }
390  }
391 
392  livePads += liveSecPads;
393  totalPads += totalSecPads;
394  if (maxHitsPerSector > 0) {
395  liveFrac = TMath::Max(0.1f, ((Float_t) liveSecPads) / (1e-15f + (Float_t) totalSecPads));
396  maxHits[sector-1] = (Int_t) (liveFrac * maxHitsPerSector);
397  if (Debug()) {LOG_INFO << "maxHits in sector " << sector << " = " << maxHits[sector-1] << endm;}
398  } else {
399  maxHits[sector-1] = 0;
400  if (Debug()) {LOG_INFO << "No maxHits in sector " << sector << endm;}
401  }
402  } // end for(sector)
403  if (maxBinZeroHits > 0) {
404  liveFrac = TMath::Max(0.1f, ((Float_t) livePads) / (1e-15f + (Float_t) totalPads));
405  maxBin0Hits = (Int_t) (liveFrac * maxBinZeroHits);
406  if (Debug()) {LOG_INFO << "maxBinZeroHits " << maxBin0Hits << endm;}
407  } else {
408  maxBin0Hits = 0;
409  if (Debug()) {LOG_INFO << "No maxBinZeroHits" << endm;}
410  }
413 
414  return kStOK;
415 }
416 //________________________________________________________________________________
417 void StTpcRTSHitMaker::PrintCld(daq_cld *cld, Int_t IdTruth, Int_t quality) {
418  if (cld) {
419  LOG_INFO << Form(" pad %f[%d:%d], tb %f[%d:%d], cha %d, fla 0x%X, Id %d, Q %d ",
420  cld->pad,
421  cld->p1,
422  cld->p2,
423  cld->tb,
424  cld->t1,
425  cld->t2,
426  cld->charge,
427  cld->flags,
428  IdTruth,
429  quality
430  ) << endm;
431  }
432 }
433 //________________________________________________________________________________
434 void StTpcRTSHitMaker::PrintAdc(daq_dta *dta) {
435  if (dta) {
436  // verify data!
437  Int_t sectorOld = -1;
438  Int_t rowOld = -1;
439  Int_t adcSum = 0;
440  while(dta && dta->iterate()) {
441  if (sectorOld != dta->sec || rowOld != dta->row) {
442  if (sectorOld > 0) {
443  LOG_INFO << Form("*** sec %2d, row %2d, Sum adc = %d",sectorOld,rowOld,adcSum) << endm;
444  }
445  sectorOld = dta->sec;
446  rowOld = dta->row;
447  adcSum = 0;
448  }
449  if (Debug() > 1) {LOG_INFO << Form("*** sec %2d, row %2d, pad %3d: %3d pixels",dta->sec,dta->row,dta->pad,dta->ncontent) << endm;}
450  for(UInt_t i=0;i<dta->ncontent;i++) {
451  if (Debug() > 1) {
452  LOG_INFO << Form(" %2d: adc %4d, tb %3d: track %4d",i,
453  dta->sim_adc[i].adc,
454  dta->sim_adc[i].tb,
455  dta->sim_adc[i].track_id
456  ) << endm;
457  }
458  adcSum += dta->sim_adc[i].adc;
459  }
460  }
461  if (sectorOld > 0) {
462  LOG_INFO << Form("*** sec %2d, row %2d, Sum adc = %d",sectorOld,rowOld,adcSum) << endm;
463  }
464  }
465 }
466 //________________________________________________________________________________
468 #ifdef __BENCHMARK__
469  TBenchmark *myBenchmark = new TBenchmark();
470  myBenchmark->Reset();
471  // myBenchmark->Start("StTpcRTSHitMaker::Make");
472 #endif
473  if (St_tpcStatusC::instance()->isDead()) {
474  LOG_WARN << "TPC status indicates it is unusable for this event. Ignoring hits." << endm;
475  return kStOK;
476  }
477  if (IAttr("TPC23")) return Make23();
478  static Short_t ADCs[__MaxNumberOfTimeBins__];
479 #ifdef __TFG__VERSION__
480  static Int_t IDTs[__MaxNumberOfTimeBins__];
481 #else
482  static UShort_t IDTs[__MaxNumberOfTimeBins__];
483 #endif
484  static StTpcCoordinateTransform transform;
485  static StThreeVectorF hard_coded_errors;
486  StEvent* rEvent = (StEvent*) GetInputDS("StEvent");
487  if (! rEvent) {
488  LOG_WARN << "There is no StEvent" << endm;
489  return kStWarn;
490  }
491  StTpcHitCollection *hitCollection = rEvent->tpcHitCollection();
492  if (! hitCollection ) {
493  hitCollection = new StTpcHitCollection();
494  rEvent->setTpcHitCollection(hitCollection);
495  }
496  TDataSet* tpcRawEvent = GetInputDS("Event");
497  if (! tpcRawEvent) {
498  LOG_WARN << "There is not Tpc Raw Event" << endm;
499  return kStWarn;
500  }
501  // if (Debug()) tpcRawEvent->ls();
502  StTpcRawData *tpcRawData = (StTpcRawData *) tpcRawEvent->GetObject();
503  if (! tpcRawData) {
504  LOG_WARN << "There is not Tpc Raw Data" << endm;
505  return kStWarn;
506  }
507  // create (or reuse) the adc_sim bank...
508  // add a bunch of adc data for a specific sector:row:pad
509  static Int_t minSector = IAttr("minSector");
510  static Int_t maxSector = IAttr("maxSector");
511  static Int_t minRow = IAttr("minRow");
512  static Int_t maxRow = IAttr("maxRow");
513  StMaker* maskMk = GetMakerInheritsFrom("StMtdTrackingMaskMaker");
514  UInt_t mask = (maskMk ? maskMk->UAttr("TpcSectorsByMtd") : ~0U); // 24 bit masking for sectors 1..24
515  bin0Hits = 0;
516  for (Int_t sector = minSector; sector <= maxSector; sector++) {
517  StTpcDigitalSector *digitalSector = tpcRawData->GetSector(sector);
518  if (! digitalSector) continue;
519  if (!((1U<<(sector-1)) & mask)) continue; // sector masking
520  UShort_t Id = 0;
521  Int_t hitsAdded = 0;
522  for (Int_t iTpcType = 1; iTpcType >= 0; iTpcType--) {// Tpx iTPC
523  daq_dta *dta = 0;
524  Int_t row1 = minRow;
525  Int_t row2 = maxRow;
526  // Check presense of iTPC and adjust row range
527  if (St_tpcPadConfigC::instance()->iTPC(sector)) {
528  // daq_tpx::get(const Char_t *bank="*", Int_t c1=-1, Int_t c2=-1, Int_t c3=-1, void *p1=0, void *p2=0)
529  // put(const Char_t *bank="*", Int_t c1=-1, Int_t c2=-1, Int_t c3=-1, void *p1=0, void *p2=0)
530  // daq_itpc::get(const Char_t *bank="*",Int_t c1=-1, Int_t c2=-1, Int_t c3=-1, void *p1=0, void *p2=0)
531  // put(const Char_t *in_bank="*", Int_t sector=-1, Int_t row=-1, Int_t pad=-1, void *p1=0, void *p2=0)
532  if (! iTpcType) { // Tpx
533  row1 = TMath::Max(row1, 41);
534  if (fTpx) dta = fTpx->put("adc_sim");
535  } else { // iTpc {
536  row2 = TMath::Min(40, row2);
537  if (fiTpc) dta = fiTpc->put("adc_sim");
538  }
539  } else { // no iTPC
540  row2 = TMath::Min(45, row2);
541  if (fTpx) dta = fTpx->put("adc_sim");
542  }
543  if (! dta) continue;
544  Int_t NoAdcs = 0;
545  for (Int_t row = row1; row <= row2; row++) {
546  if (! St_tpcPadGainT0BC::instance()->livePadrow(sector,row)) continue;
547  Int_t Npads = digitalSector->numberOfPadsInRow(row);
548  Int_t rowO = row;
549  // daq_dta *dta = 0;
550  if (! iTpcType) {
551  if (St_tpcPadConfigC::instance()->iTPC(sector) && row > 40) { // Tpx sector with iTPC
552  rowO = row - 40 + 13; // old row count
553  }
554  }
555  for(Int_t pad = 1; pad <= Npads; pad++) {
556  UInt_t ntimebins = digitalSector->numberOfTimeBins(row,pad);
557  if (! ntimebins) continue;
558  // allocate space for at least 512 pixels (timebins)
559  daq_sim_adc_tb *d = (daq_sim_adc_tb *) dta->request(__MaxNumberOfTimeBins__);
560  // add adc data for this specific sector:row:pad
561  digitalSector->getTimeAdc(row,pad,ADCs,IDTs);
562  UInt_t l = 0;
563  for (UInt_t k = 0; k < __MaxNumberOfTimeBins__; k++) {
564  if (ADCs[k]) {
565  d[l].adc = ADCs[k];
566  d[l].tb = k;
567  d[l].track_id = IDTs[k];
568  l++;
569  }
570  }
571  if (l > 0) {
572 #ifdef __BENCHMARK__
573  // myBenchmark->Start("StTpcRTSHitMaker::Make::finalize");
574 #endif
575  // daq_dta::finalize(uint32_t obj_cou, Int_t s=0, Int_t row=0, Int_t pad=0) ;
576  dta->finalize(l,sector,rowO,pad);
577 #ifdef __BENCHMARK__
578  // myBenchmark->Stop("StTpcRTSHitMaker::Make::finalize");
579 #endif
580  NoAdcs += l;
581  }
582  } // pad loop
583  } // row loop
584  if (! NoAdcs) continue;
585  daq_dta *dtaX = 0;
586  if (! iTpcType) {
587  if (fTpx) dtaX = fTpx->get("adc_sim");
588  } else {
589  if (fiTpc) dtaX = fiTpc->get("adc_sim");
590  }
591  dta = dtaX;
592  if (Debug()) {
593  PrintAdc(dta);
594  }
595  dtaX = 0;
596  if (! iTpcType) { //Tpx
597  if (! fTpx) continue;
598  if (IAttr("TpxClu2D")) {
599  dtaX = fTpx->get("cld_2d_sim");
600  } else {
601  dtaX = fTpx->get("cld_sim");;
602  }
603  } else { // iTpc
604  if (! fiTpc) continue;
605  dtaX = fiTpc->get("cld_sim");
606  }
607  if (! dtaX) continue;
608  daq_dta *dd = dtaX;
609 #ifndef __TFG__VERSION__
610  Double_t ADC2GeV = 0;
611  Int_t rowOld = -1;
612 #endif /* ! __TFG__VERSION__ */
613  static Int_t iBreak = 0;
614  while(dd && dd->iterate()) {
615  if (Debug()) {
616  LOG_INFO << Form("CLD sec %2d: row %2d: %d clusters",dd->sec, dd->row, dd->ncontent) << endm;
617  }
618  for(UInt_t i=0;i<dd->ncontent;i++) {
619  daq_cld cld;
620  Int_t IdTruth = 0;
621  UShort_t quality = 0;
622  if (! iTpcType) {
623  daq_sim_cld &dc = ((daq_sim_cld *) dd->Void)[i] ;
624  cld = dc.cld;
625  IdTruth = dc.track_id;
626  quality = dc.quality;
627  } else {
628  daq_sim_cld_x &dc_x = ((daq_sim_cld_x *) dd->Void)[i] ;
629  cld = dc_x.cld;
630  IdTruth = dc_x.track_id;
631  quality = dc_x.quality;
632  }
633  if (Debug()) {
634  PrintCld(&cld, IdTruth, quality);
635  iBreak++;
636  }
637  if (cld.p1 > cld.p2) continue;
638  if (cld.t1 > cld.t2) continue;
639  if (cld.tb >= __MaxNumberOfTimeBins__) continue;
640  if (cld.charge < fminCharge) continue;
641  if ( ! (cld.pad > 0 && cld.pad <= 182 &&
642  cld.tb >= 0 && cld.tb < 512)) continue;
643  /*tpxFCF.h
644  #define FCF_ONEPAD 1
645  #define FCF_DOUBLE_PAD 2 // offline: merged
646  #define FCF_MERGED 2
647  #define FCF_BIG_CHARGE 8
648  #define FCF_ROW_EDGE 16 // 0x10 touched end of row
649  #define FCF_BROKEN_EDGE 32 // 0x20 touches one of the mezzanine edges
650  #define FCF_DEAD_EDGE 64 // 0x40 touches a dead pad
651  */
652  if ( cld.flags && (cld.flags & ~(FCF_ONEPAD | FCF_MERGED | FCF_BIG_CHARGE))) continue;
653  Int_t rowO = dd->row;
654  Int_t padrow = rowO;
655  if (! iTpcType && St_tpcPadConfigC::instance()->iTPC(sector) && rowO > 13) padrow = rowO + 40 - 13;
656  Int_t row = padrow;
657  Float_t pad = cld.pad;
658  Int_t iRdo = StDetectorDbTpcRDOMasks::instance()->rdoForPadrow(sector,row,pad);
659  if ( ! StDetectorDbTpcRDOMasks::instance()->isOn(sector,iRdo)) continue;
660 
661  StTpcPadCoordinate Pad(dd->sec, padrow, cld.pad, cld.tb); PrPP(Make,Pad);
662  static StTpcLocalSectorCoordinate LS;
663  static StTpcLocalCoordinate L;
664  transform(Pad,LS,kFALSE,kTRUE); PrPP(Make,LS); // don't useT0, useTau
665  transform(LS,L); PrPP(Make,L);
666  UInt_t hw = 1; // detid_tpc
667  //yf if (isiTpcSector) hw += 1U << 1;
668  hw += dd->sec << 4; // (padrow/100 << 4); // sector
669  hw += padrow << 9; // (padrow%100 << 9); // padrow
670 #ifndef __TFG__VERSION__
671  if (padrow != rowOld) {
672  rowOld = padrow;
673  Double_t gain = St_tpcPadConfigC::instance()->IsRowInner(dd->sec,padrow) ?
674  St_tss_tssparC::instance()->gain_in() :
675  St_tss_tssparC::instance()->gain_out();
676  Double_t wire_coupling = St_tpcPadConfigC::instance()->IsRowInner(dd->sec,padrow) ?
677  St_tss_tssparC::instance()->wire_coupling_in() :
678  St_tss_tssparC::instance()->wire_coupling_out();
679  ADC2GeV = ((Double_t) St_tss_tssparC::instance()->ave_ion_pot() *
680  (Double_t) St_tss_tssparC::instance()->scale())/(gain*wire_coupling) ;
681  }
682  Double_t q = ADC2GeV*cld.charge;
683 #else /* used in TFG till 07/31/20 */
684  Double_t q = 0;
685 #endif /* ! __TFG__VERSION__ */
686  Id++;
687  StTpcHit *hit = StTpcHitMaker::StTpcHitFlag(L.position(),hard_coded_errors,hw,q
688  , (UChar_t ) 0 // counter
689  , IdTruth
690  , quality
691  , Id // id =0,
692  , cld.p1 // mnpad
693  , cld.p2 // mxpad
694  , cld.t1 // mntmbk
695  , cld.t2 // mxtmbk
696  , cld.pad
697  , cld.tb
698  , cld.charge
699  , cld.flags);
700  if (! hit) continue;
701  hitsAdded++;
702  if (hit->minTmbk() == 0) bin0Hits++;
703  if (Debug() > 1) hit->Print();
704  hitCollection->addHit(hit);
705  if (Debug() > 1) {cout << "Add hit #" << hitCollection->numberOfHits() << endl;}
706  }
707  }
708  } // end iTpcType loop
709  if (maxHits[sector-1] && hitsAdded > maxHits[sector-1]) {
710  LOG_ERROR << "Too many hits (" << hitsAdded << ") in one sector ("
711  << sector << "). Skipping event." << endm;
712  return kStSkip;
713  }
714  } // end sec loop
715  if (maxBin0Hits && bin0Hits > maxBin0Hits) {
716  LOG_ERROR << "Too many hits (" << bin0Hits
717  << ") starting at time bin 0. Skipping event." << endm;
718  return kStSkip;
719  }
720  if (! IAttr("NoTpxAfterBurner")) StTpcHitMaker::AfterBurner(hitCollection);
721 #ifdef __BENCHMARK__
722  myBenchmark->Stop("StTpcRTSHitMaker::Make");
723  myBenchmark->Show("StTpcRTSHitMaker::Make");
724  myBenchmark->Show("StTpcRTSHitMaker::Make::finalize");
725  delete myBenchmark;
726 #endif
727  return kStOK;
728 }
729 //________________________________________________________________________________
730 Int_t StTpcRTSHitMaker::Make23() {
731  if (! fTpx23 && ! fiTpc23) return kStErr;
732 #ifdef __BENCHMARK__
733  TBenchmark *myBenchmark = new TBenchmark();
734  myBenchmark->Reset();
735  // myBenchmark->Start("StTpcRTSHitMaker::Make");
736 #endif
737  if (St_tpcStatusC::instance()->isDead()) {
738  LOG_WARN << "TPC status indicates it is unusable for this event. Ignoring hits." << endm;
739  return kStOK;
740  }
741  static Short_t ADCs[__MaxNumberOfTimeBins__];
742 #ifdef __TFG__VERSION__
743  static Int_t IDTs[__MaxNumberOfTimeBins__];
744 #else
745  static UShort_t IDTs[__MaxNumberOfTimeBins__];
746 #endif
747  static StTpcCoordinateTransform transform;
748  static StThreeVectorF hard_coded_errors;
749  StEvent* rEvent = (StEvent*) GetInputDS("StEvent");
750  if (! rEvent) {
751  LOG_WARN << "There is no StEvent" << endm;
752  return kStWarn;
753  }
754  StTpcHitCollection *hitCollection = rEvent->tpcHitCollection();
755  if (! hitCollection ) {
756  hitCollection = new StTpcHitCollection();
757  rEvent->setTpcHitCollection(hitCollection);
758  }
759  TDataSet* tpcRawEvent = GetInputDS("Event");
760  if (! tpcRawEvent) {
761  LOG_WARN << "There is not Tpc Raw Event" << endm;
762  return kStWarn;
763  }
764  // if (Debug()) tpcRawEvent->ls();
765  StTpcRawData *tpcRawData = (StTpcRawData *) tpcRawEvent->GetObject();
766  if (! tpcRawData) {
767  LOG_WARN << "There is not Tpc Raw Data" << endm;
768  return kStWarn;
769  }
770  // create (or reuse) the adc_sim bank...
771  // add a bunch of adc data for a specific sector:row:pad
772  static Int_t minSector = IAttr("minSector");
773  static Int_t maxSector = IAttr("maxSector");
774  static Int_t minRow = IAttr("minRow");
775  static Int_t maxRow = IAttr("maxRow");
776  bin0Hits = 0;
777  for (Int_t sector = minSector; sector <= maxSector; sector++) {
778  StTpcDigitalSector *digitalSector = tpcRawData->GetSector(sector);
779  if (! digitalSector) continue;
780  UShort_t Id = 0;
781  Int_t hitsAdded = 0;
782  // daq_dta *dta = 0;
783  for (Int_t iTpcType = 1; iTpcType >= 0; iTpcType--) {// Tpx iTPC
784  tpc23_base *tpc23 = 0;
785  Int_t row1 = minRow;
786  Int_t row2 = maxRow;
787  // Check presense of iTPC and adjust row range
788  if (St_tpcPadConfigC::instance()->iTPC(sector)) {
789  // daq_tpx::get(const Char_t *bank="*", Int_t c1=-1, Int_t c2=-1, Int_t c3=-1, void *p1=0, void *p2=0)
790  // put(const Char_t *bank="*", Int_t c1=-1, Int_t c2=-1, Int_t c3=-1, void *p1=0, void *p2=0)
791  // daq_itpc::get(const Char_t *bank="*",Int_t c1=-1, Int_t c2=-1, Int_t c3=-1, void *p1=0, void *p2=0)
792  // put(const Char_t *in_bank="*", Int_t sector=-1, Int_t row=-1, Int_t pad=-1, void *p1=0, void *p2=0)
793  if (! iTpcType) { // Tpx
794  row1 = TMath::Max(row1, 41);
795  tpc23 = fTpx23;
796  } else { // iTpc {
797  row2 = TMath::Min(40, row2);
798  tpc23 = fiTpc23;
799  }
800  } else { // no iTPC
801  row2 = TMath::Min(45, row2);
802  tpc23 = fTpx23;
803  }
804  if (! tpc23) continue;
805 #ifndef __TFG__VERSION__
806  Double_t ADC2GeV = 0;
807  Int_t rowOld = -1;
808 #endif /* ! __TFG__VERSION__ */
809  for (Int_t row = row1; row <= row2; row++) {
810  if (! St_tpcPadGainT0BC::instance()->livePadrow(sector,row)) continue;
811  tpc23->sim_evt_start(sector) ; // prepare start of event
812  Int_t NoAdcs = 0;
813  Int_t Npads = digitalSector->numberOfPadsInRow(row);
814  Int_t rowO = row;
815  // daq_dta *dta = 0;
816  if (! iTpcType) {
817  if (St_tpcPadConfigC::instance()->iTPC(sector) && row > 40) { // Tpx sector with iTPC
818  rowO = row - 40 + 13; // old row countx1
819  }
820  }
821  for(Int_t pad = 1; pad <= Npads; pad++) {
822  UInt_t ntimebins = digitalSector->numberOfTimeBins(row,pad);
823  if (! ntimebins) continue;
824  // allocate space for at least 512 pixels (timebins)
825  // daq_sim_adc_tb *d = (daq_sim_adc_tb *) dta->request(__MaxNumberOfTimeBins__);
826  // add adc data for this specific sector:row:pad
827  digitalSector->getTimeAdc(row,pad,ADCs,IDTs);
828  UInt_t l = 0;
829  static UShort_t adcs[__MaxNumberOfTimeBins__];
830  memset(adcs, 0, sizeof(adcs));
831  static Int_t idts[__MaxNumberOfTimeBins__];
832  memset(idts, 0, sizeof(idts));
833  for (UInt_t k = 0; k < __MaxNumberOfTimeBins__; k++) {
834  if (ADCs[k]) {
835  adcs[k] = ADCs[k];
836  idts[k] = IDTs[k];
837  l++;
838  }
839  }
840  if (l > 0) {
841  Int_t padrow = rowO;
842  // tpc23->sim_do_pad(padrow,pad,ADCs,IDTs) ;
843  tpc23->do_ch_sim(padrow,pad,adcs,idts) ;
844  NoAdcs += l;
845  }
846  } // pad loop
847 
848  if (NoAdcs) {
849 
850  // Tonko's ------------------------------------------------------------
851 
852  // allocate output storage based upon the count of found sequences
853  const Int_t words_per_cluster = 5 ; // 5 for simulation, 2 normally
854 
855  tpc23->s2_max_words = tpc23->sequence_cou*words_per_cluster + 2000 ; // and add a bit more...
856 
857  tpc23->s2_start = (UInt_t *) malloc(tpc23->s2_max_words*4) ;
858 
859  // this actually runs the clusterfinder
860  tpc23->evt_stop() ;
861 
862 
863  if(tpc23->s2_words) { // if anything found..
864  UInt_t *p_buff = tpc23->s2_start ;
865  UInt_t *end_buff = p_buff + tpc23->s2_words ;
866  if (Debug()) {
867  cout << Form("*** sequences found %d",tpc23->sequence_cou) << endl ;
868  }
869 
870  if(tpc23->s2_words >= tpc23->s2_max_words) {
871  LOG_ERROR << "Whoa -- lots of words " << tpc23->s2_words << "/" << tpc23->s2_max_words << "\t"
872  << "\tsector " << sector << "\trow " << row << endm;
873  }
874 
875  while(p_buff < end_buff) {
876  UInt_t padrow ;
877  UInt_t version ;
878  UInt_t int_cou ;
879  Int_t ints_per_cluster ;
880  // TPX and iTPC have slightly different formats; maintained compatibility
881  if(tpc23 == fTpx23) {// det==0) { // TPX
882  padrow = *p_buff++ ;
883  int_cou = *p_buff++ ;
884 
885  version = (padrow>>16) ;
886 
887  ints_per_cluster = fTpx23 -> online? 2 : 5; // 5 for sim, 2 for real, Tommy's fix
888  }
889  else {
890  padrow = *p_buff++ ;
891  version = *p_buff++ ;
892  int_cou = *p_buff++ ;
893 
894  ints_per_cluster = (padrow>>16) ;
895 
896  }
897 
898  padrow &= 0xFFFF ;
899 
900  Int_t clusters = int_cou / ints_per_cluster ;
901 
902  for(Int_t i=0;i<clusters;i++) {
903  daq_sim_cld_x dc ;
904  tpc23->fcf_decode(p_buff,&dc,version) ;
905 
906  // nice flags printout
907  char c_flags[128] ;
908  c_flags[0] = 0 ;
909 
910  if(dc.cld.flags & FCF_ONEPAD) strcat(c_flags,"one+") ;
911  if(dc.cld.flags & FCF_MERGED) strcat(c_flags,"merge+") ;
912  if(dc.cld.flags & FCF_DEAD_EDGE) strcat(c_flags,"dead+") ;
913  if(dc.cld.flags & FCF_ROW_EDGE) strcat(c_flags,"edge+") ; if(dc.cld.flags & FCF_ONEPAD) strcat(c_flags," one") ;
914  if(dc.cld.flags & FCF_BROKEN_EDGE) strcat(c_flags,"small+") ;
915  if(dc.cld.flags & FCF_BIG_CHARGE) strcat(c_flags,"charge+") ;
916 
917  if(strlen(c_flags)) {
918  c_flags[strlen(c_flags)-1] = 0 ;
919  }
920  if (Debug()) {
921  cout << Form("row %d/%d: %f %d %d %f %d %d %d 0x%02X[%s]",row,rowO,
922  dc.cld.pad,dc.cld.p1,dc.cld.p2,
923  dc.cld.tb,dc.cld.t1,dc.cld.t2,
924  dc.cld.charge,
925  dc.cld.flags,c_flags) << endl;
926  cout << Form(" track_id %u, quality %d, pixels %d, max_adc %d",
927  dc.reserved[0],
928  dc.quality,
929  dc.pixels, //<<<<<<<<<< Add to TpcHit ?
930  dc.max_adc) << endl;; //<<<<<<<<<<
931  }
932 
933  p_buff += ints_per_cluster ;
934  // Tonko's end ------------------------------------------------------------
935  Int_t IdTruth = dc.reserved[0];
936  UShort_t quality = dc.quality;
937  // Fill TpcHit
938  if (dc.cld.p1 > dc.cld.p2) continue;
939  if (dc.cld.t1 > dc.cld.t2) continue;
940  if (dc.cld.tb >= __MaxNumberOfTimeBins__) continue;
941  if (dc.cld.charge < fminCharge) continue;
942  if ( ! (dc.cld.pad > 0 && dc.cld.pad <= 182 &&
943  dc.cld.tb >= 0 && dc.cld.tb < 512)) continue;
944  /*tpxFCF.h
945  #define FCF_ONEPAD 1
946  #define FCF_DOUBLE_PAD 2 // offline: merged
947  #define FCF_MERGED 2
948  #define FCF_BIG_CHARGE 8
949  #define FCF_ROW_EDGE 16 // 0x10 touched end of row
950  #define FCF_BROKEN_EDGE 32 // 0x20 touches one of the mezzanine edges
951  #define FCF_DEAD_EDGE 64 // 0x40 touches a dead pad
952  */
953  if ( dc.cld.flags && (dc.cld.flags & ~(FCF_ONEPAD | FCF_MERGED | FCF_BIG_CHARGE))) continue;
954  Float_t pad = dc.cld.pad;
955  Int_t iRdo = StDetectorDbTpcRDOMasks::instance()->rdoForPadrow(sector,row,pad);
956  if ( ! StDetectorDbTpcRDOMasks::instance()->isOn(sector,iRdo)) continue;
957 
958  StTpcPadCoordinate Pad(sector, row, dc.cld.pad, dc.cld.tb); PrPP(Make,Pad);
959  static StTpcLocalSectorCoordinate LS;
960  static StTpcLocalCoordinate L;
961  transform(Pad,LS,kFALSE,kTRUE); PrPP(Make,LS); // don't useT0, useTau
962  transform(LS,L); PrPP(Make,L);
963  UInt_t hw = 1; // detid_tpc
964  //yf if (isiTpcSector) hw += 1U << 1;
965  hw += sector << 4; // (row/100 << 4); // sector
966  hw += row << 9; // (row%100 << 9); // row
967 #ifndef __TFG__VERSION__
968  if (row != rowOld) {
969  rowOld = row;
970  Double_t gain = St_tpcPadConfigC::instance()->IsRowInner(sector,row) ?
971  St_tss_tssparC::instance()->gain_in() :
972  St_tss_tssparC::instance()->gain_out();
973  Double_t wire_coupling = St_tpcPadConfigC::instance()->IsRowInner(sector,row) ?
974  St_tss_tssparC::instance()->wire_coupling_in() :
975  St_tss_tssparC::instance()->wire_coupling_out();
976  ADC2GeV = ((Double_t) St_tss_tssparC::instance()->ave_ion_pot() *
977  (Double_t) St_tss_tssparC::instance()->scale())/(gain*wire_coupling) ;
978  }
979  Double_t q = ADC2GeV*dc.cld.charge;
980 #else /* used in TFG till 07/31/20 */
981  Double_t q = 0;
982 #endif /* ! __TFG__VERSION__ */
983  Id++;
984  StTpcHit *hit = StTpcHitMaker::StTpcHitFlag(L.position(),hard_coded_errors,hw,q
985  , (UChar_t ) 0 // counter
986  , IdTruth
987  , quality
988  , Id // id =0,
989  , dc.cld.p1 // mnpad
990  , dc.cld.p2 // mxpad
991  , dc.cld.t1 // mntmbk
992  , dc.cld.t2 // mxtmbk
993  , dc.cld.pad
994  , dc.cld.tb
995  , dc.cld.charge
996  , dc.cld.flags);
997  if (! hit) continue;
998  hitsAdded++;
999  if (hit->minTmbk() == 0) bin0Hits++;
1000  if (Debug() > 1) hit->Print();
1001  hitCollection->addHit(hit);
1002  if (Debug() > 1) {cout << "Add hit #" << hitCollection->numberOfHits() << endl;}
1003  }
1004  }
1005  tpc23->run_stop() ; // dumps some statistics...
1006 
1007  // free allocated storage
1008  free(tpc23->s2_start) ;
1009  tpc23->s2_start = 0 ;
1010  } // clusters loop
1011  } // row loop
1012 
1013  } // end of row
1014 
1015 
1016  } // end iTpcType loop
1017  } // end of sector
1018 
1019 
1020 
1021  if (! IAttr("NoTpxAfterBurner")) StTpcHitMaker::AfterBurner(hitCollection);
1022 #ifdef __BENCHMARK__
1023  myBenchmark->Stop("StTpcRTSHitMaker::Make");
1024  myBenchmark->Show("StTpcRTSHitMaker::Make");
1025  myBenchmark->Show("StTpcRTSHitMaker::Make::finalize");
1026  delete myBenchmark;
1027 #endif
1028  return kStOK;
1029 }
1030 //________________________________________________________________________________
1031 TH2F *StTpcRTSHitMaker::PlotSecRow(Int_t sec, Int_t row, Int_t flags) {
1032  /*
1033  Usage:
1034  root.exe ' bfc.C(1,"P2019a,StiCA,evout,NoHistos,noTags,noRunco,PicoVtxVpdOrDefault,TpxRaw,TPC23,TpxRaw,TPC23,USE_GAIN_FROM_FILE,tpxO,NoTpxAfterBurner","/RTS/TONKO/data/st_physics_adc_20192001_raw_5500002.daq","","Plot.root")'
1035  TH2F *plot = StTpcRTSHitMaker::instance()->PlotSecRow(2,43,-1);
1036  */
1037  TH2F *plot = 0;
1038  StEvent* pEvent = (StEvent*) StMaker::GetChain()->GetInputDS("StEvent");
1039  if (!pEvent) { cout << "Can't find StEvent" << endl; return plot;}
1040  StTpcHitCollection* TpcHitCollection = pEvent->tpcHitCollection();
1041  if (! TpcHitCollection) { cout << "No TPC Hit Collection" << endl; return plot;}
1042  StTpcSectorHitCollection* sectorCollection = TpcHitCollection->sector(sec-1);
1043  if (! sectorCollection) { cout << "No TPC Hit Collection for sector " << sec << endl; return plot;}
1044  StTpcPadrowHitCollection *rowCollection = sectorCollection->padrow(row-1);
1045  if (! rowCollection) { cout << "No TPC Hit Collection for sector " << sec << " and row " << row << endl; return plot;}
1046  TFile *f = StMaker::GetChain()->GetTFile();
1047  if (! f) {cout << "No TFile is opened" << endl; return plot;}
1048  f->cd();
1049  TString PlotName(Form("s%02ir%03i",sec,row));
1050  plot = (TH2F *) f->Get(PlotName);
1051  if (plot) {
1052  plot->Reset();
1053  plot->SetTitle(Form("ADC versus pad and time buckets for sec = %i and row = %i in event %i",sec,row, pEvent->id()));
1054  } else {
1055  Int_t npad = St_tpcPadConfigC::instance()->padsPerRow(sec,row);
1056  plot = new TH2F(PlotName,Form("ADC versus pad and time buckets for sec = %i and row = %i in event %i",sec,row, pEvent->id()), 360,-0.5,359.5,npad,0.5,npad+0.5);
1057  plot->SetStats(1);
1058  }
1059  TDataSet* tpcRawEvent = StMaker::GetChain()->GetInputDS("Event");
1060  if (! tpcRawEvent) {
1061  LOG_WARN << "There is not Tpc Raw Event" << endm;
1062  } else {
1063  StTpcRawData *tpcRawData = (StTpcRawData *) tpcRawEvent->GetObject();
1064  if (! tpcRawData) {
1065  LOG_WARN << "There is not Tpc Raw Data" << endm;
1066  } else {
1067  StTpcDigitalSector *digitalSector = tpcRawData->GetSector(sec);
1068  if (! digitalSector) {
1069  LOG_WARN << "There is not digital sector" << sec << endm;
1070  } else {
1071  Int_t Npads = digitalSector->numberOfPadsInRow(row);
1072  for(Int_t pad = 1; pad <= Npads; pad++) {
1073  static Short_t ADCs[__MaxNumberOfTimeBins__];
1074 #ifdef __TFG__VERSION__
1075  static Int_t IDTs[__MaxNumberOfTimeBins__];
1076 #else
1077  static UShort_t IDTs[__MaxNumberOfTimeBins__];
1078 #endif
1079  digitalSector->getTimeAdc(row,pad,ADCs,IDTs);
1080  for (UInt_t k = 0; k < __MaxNumberOfTimeBins__; k++) {
1081  if (ADCs[k]) plot->Fill(k,pad,ADCs[k]);
1082  }
1083  }
1084  }
1085  }
1086  }
1087  StSPtrVecTpcHit &hits = rowCollection->hits();
1088  Long_t NoHits = hits.size();
1089  for (Long64_t k = 0; k < NoHits; k++) {
1090  const StTpcHit *tpcHit = static_cast<const StTpcHit *> (hits[k]);
1091  Int_t color = 1;
1092  Int_t style = 43;
1093  Double_t offset = 0.0;
1094  Int_t flag = tpcHit->flag();
1095  if (flag & 256) {color = 2; style = 29; offset = 0.1; flag &= ~0x100;} // online kTpxO
1096  if (flags > -1 && flag < flags) continue;
1097  Double_t tb = tpcHit->timeBucket();
1098  Double_t pad = tpcHit->pad();
1099  TPolyMarker *pm = new TPolyMarker(1,&tb, &pad);
1100  pm->SetMarkerStyle(style);
1101  pm->SetMarkerColor(color);
1102  pm->SetMarkerSize(2);
1103  plot->GetListOfFunctions()->Add(pm);
1104  TBox *box = new TBox(tpcHit->minTmbk()-0.5 + offset,tpcHit->minPad()-0.5 + offset, tpcHit->maxTmbk()+0.5 + offset,tpcHit->maxPad()+0.5 + offset);
1105  box->SetLineWidth(4);
1106  box->SetLineColor(color);
1107  box->SetFillStyle(0);
1108  box->SetFillColor(0);
1109  plot->GetListOfFunctions()->Add(box);
1110  }
1111  plot->Draw("colz");
1112  return plot;
1113 }
Definition: itpc23.h:8
void version(std::ostream &os=std::cout)
print HepMC version
Definition: Version.h:27
Definition: Stypes.h:49
virtual TObject * GetObject() const
The depricated method (left here for the sake of the backward compatibility)
Definition: TDataSet.cxx:428
Definition: Stypes.h:42
Definition: Stypes.h:40
Definition: tpx23.h:48
Definition: Stypes.h:44