StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StSsdDaqMaker.cxx
1 // $Id: StSsdDaqMaker.cxx,v 1.22 2017/04/26 20:12:06 perev Exp $
2 //
3 // $log$
4 //
5 // Id: StSsdDaqMaker.cxx,v 1.5 2005/04/22 16:12:02 lmartin Exp $
6 //
7 // Log: StSsdDaqMaker.cxx,v $
8 // Revision 1.5 2005/04/22 16:12:02 lmartin
9 // bug in the ladder index fixed
10 //
11 // Revision 1.4 2005/04/21 14:59:47 lmartin
12 // useless print removed
13 //
14 // Revision 1.3 2005/04/21 14:55:17 lmartin
15 // bug in the offset correction fixed
16 //
17 // Revision 1.2 2005/04/21 09:50:28 lmartin
18 // Hardware offset corrected for specific ladders
19 //
20 // Revision 1.1 2005/04/15 15:11:24 lmartin
21 // StSsdDaqMaker
22 //
23 
24 //*-- Authors : Lilian Martin
25 //*-- : Joerg Reinnarth
26 //*-- : Jonathan Bouchet
27 
28 #include "StSsdDaqMaker.h"
29 #include "TDataSetIter.h"
30 #include "StMessMgr.h"
31 #include "StDAQMaker/StDAQReader.h"
32 #include "StDAQMaker/StSSDReader.h"
33 #include "ssdLadderMap.h"
34 #include "tables/St_spa_strip_Table.h"
35 #include "tables/St_ssdConfiguration_Table.h"
36 #include "tables/St_ssdPedStrip_Table.h"
37 #include "StSsdUtil/StSsdConfig.hh"
38 #include "TH1.h"
39 #include "TH2.h"
40 #include "TFile.h"
41 #include "TNtuple.h"
42 ClassImp(StSsdDaqMaker)
43 
44 //_____________________________________________________________________________
46 
53 StSsdDaqMaker::StSsdDaqMaker(const char *name):StMaker(name){
54  //
55 }
56 
57 
58 //_____________________________________________________________________________
60 
70  //
71 }
72 
73 
74 //_____________________________________________________________________________
77  // Now we try to get the configuration information (which ladders are active)
78  // table using the StSsdUtil and StSsdDbMaker...
79 
80 
81  LOG_INFO << "Read now Databases " << endm;
82  // Create SsdPedestal histograms
83  if (IAttr(".histos")) {
84  occupancy_wafer = new TH2S("occupancy_wafer","occupancy per wafer",40,0,40,20,0,20);
85  occupancy_chip = new TH2S("occupancy_chip","occupancy per chip",40,0,40,99,0,99);
86  noise_chip = new TH2S("noise_chip","mean noise per chip",40,0,40,99,0,99);
87  noise_wafer = new TH2S("noise_wafer","mean noise per wafer",40,0,40,20,0,20);
88  noise_chip_P = new TH2S("noise_chip_P","mean noise per chip for the P Side ",20,0,20,96,0,96);
89  noise_chip_N = new TH2S("noise_chip_N","mean noise per chip for the N Side",20,0,20,96,0,96);
90  pedestal_chip = new TH2S("pedestal_chip","pedestal per chip",40,0,40,99,0,99);
91  ped_zero_ladP = new TH2S("ped_zero_ladP","map of p-side wafers where strips have ped = 0",20,0,20,15,0,15);
92  ped_zero_ladN = new TH2S("ped_zero_ladN","map of n-side wafers where strips have ped = 0",20,0,20,15,0,15);
93  ped_high_ladP = new TH2S("ped_high_ladP","map of p-side wafers where strips have ped = 255",20,0,20,15,0,15);
94  ped_high_ladN = new TH2S("ped_high_ladN","map of n-side wafers where strips have ped = 255",20,0,20,15,0,15);
95  noise_zero_ladP = new TH2S("noise_zero_ladP","map of p-side wafers where strips have noise = 0",20,0,20,15,0,15);
96  noise_zero_ladN = new TH2S("noise_zero_ladN","map of n-side wafers where strips have noise = 0",20,0,20,15,0,15);
97  noise_high_ladP = new TH2S("noise_high_ladP","map of p-side wafers where strips have noise = 255",20,0,20,15,0,15);
98  noise_high_ladN = new TH2S("noise_high_ladN","map of n-side wafers where strips have noise = 255",20,0,20,15,0,15);
99  occupancy = new TH1F("deadStrips","number of dead strips",40,0,40);
100  occupancy_wafer->SetXTitle("Ladder");
101  occupancy_wafer->SetYTitle("Wafer");
102  occupancy_chip->SetXTitle("Ladder");
103  occupancy_chip->SetYTitle("Chip");
104  noise_chip->SetXTitle("Ladder");
105  noise_chip->SetYTitle("Chip");
106  noise_chip_P->SetXTitle("Ladder");
107  noise_chip_N->SetXTitle("Ladder");
108  noise_chip_P->SetYTitle("Chip");
109  noise_chip_N->SetYTitle("Chip");
110  noise_wafer->SetXTitle("Ladder");
111  noise_wafer->SetYTitle("Wafer");
112  pedestal_chip->SetXTitle("Ladder");
113  pedestal_chip->SetYTitle("Chip");
114  occupancy->SetXTitle("Ladder");
115  occupancy->SetYTitle("#");
116  ped_zero_ladP->SetXTitle("Ladder");
117  ped_zero_ladP->SetYTitle("Wafer");
118  ped_zero_ladN->SetXTitle("Ladder");
119  ped_zero_ladN->SetYTitle("Wafer");
120  ped_high_ladP->SetXTitle("Ladder");
121  ped_high_ladP->SetYTitle("Wafer");
122  ped_high_ladN->SetXTitle("Ladder");
123  ped_high_ladN->SetYTitle("Wafer");
124  noise_zero_ladP->SetXTitle("Ladder");
125  noise_zero_ladP->SetYTitle("Wafer");
126  noise_zero_ladN->SetXTitle("Ladder");
127  noise_zero_ladN->SetYTitle("Wafer");
128  noise_high_ladP->SetXTitle("Ladder");
129  noise_high_ladP->SetYTitle("Wafer");
130  noise_high_ladN->SetXTitle("Ladder");
131  noise_high_ladN->SetYTitle("Wafer");
132  mPedOut = 0;//default : we do not write the Tuple
133  if(mPedOut)
134  {
135  DeclareNTuple();
136  }
137  }
138  LOG_INFO << "Init() - Done " << endm;
139 
140  return StMaker::Init();
141 }
142 
143 //_____________________________________________________________________________
144 Int_t StSsdDaqMaker::InitRun(int runumber)
145 {
146  LOG_INFO <<"InitRun(int runumber) - Read now Databases"<<endm;
147  Int_t run = (runumber/1000000)-1;
148 
149  St_ssdConfiguration *configuration = (St_ssdConfiguration*) GetDataBase("Geometry/ssd/ssdConfiguration");
150  if (!configuration){
151  LOG_ERROR << "InitRun("<<runumber<<") - ERROR - ssdConfiguration==0"<<endm;
152  return 0;
153  }
154  ssdConfiguration_st *config = (ssdConfiguration_st*) configuration->GetTable() ;
155  if (!config){
156  LOG_ERROR <<"InitRun("<<runumber<<") - ERROR - config==0"<<endm;
157  return 0;
158  }
159 
160  mConfig = new StSsdConfig();
161 
162  int totLadderPresent = 0;
163 
164  for (int ladder = 1; ladder<=config->nMaxLadders;ladder++)
165  {
166  LOG_INFO<< " on sector = " << config->ladderIsPresent[ladder-1];
167  if (config->ladderIsPresent[ladder-1] != 0)
168  totLadderPresent++;
169  mConfig->setLadderIsActive(ladder,config->ladderIsPresent[ladder-1]);
170  }
171  PrintConfiguration(run,config);
172  mConfig->setNumberOfLadders(totLadderPresent);
173  mConfig->setNumberOfWafers(config->nMaxWafers/config->nMaxLadders);
174  mConfig->setNumberOfHybrids(2);
175  mConfig->setTotalNumberOfHybrids(2*16*totLadderPresent);
176  mConfig->setTotalNumberOfLadders(config->nMaxLadders);
177  mConfig->setNumberOfStrips(768);
178 
179  mConfig->setConfiguration();
180 
181  LOG_INFO << "_____________________________" << endm;
182  LOG_INFO << " Via Datababase......." << endm;
183  LOG_INFO << ".......numberOfSectors = " << config->nMaxSectors << endm;
184  LOG_INFO << ".......numberOfLadders = " << totLadderPresent << endm;
185  LOG_INFO << " .numberOfWafersPerLadder = " << config->nMaxWafers/config->nMaxLadders << endm;
186  LOG_INFO << "_____________________________" << endm;
187  LOG_INFO << " InitRun() - Done "<<endm;
188  return kStOk;
189 
190 }
191 
192 //_____________________________________________________________________________
193 
194 void StSsdDaqMaker::DeclareNTuple(){
195  pFile = new TFile("PedestalFile.root","RECREATE");
196  string varlist = "side:ladder:wafer:strip:pedestal:noise:chip:id_wafer:tot_strip";
197  pTuple = new TNtuple("PedestalNTuple","Pedestal Ntuple",varlist.c_str());
198  LOG_INFO << "StSsdDaqMaker::DeclareNtuple - Done - "<<endm;
199 }
200 
201 //_____________________________________________________________________________
202 // Make - this method is called in loop for each event
203 //
204 // The real data are saved in the spa_strip table
205 // The pedestal data are saved in the ssdPedStrip table
206 //
207 //
208 //
209 //
211  int strip_number,id_wafer,id_side,count,my_channel=-9999;
212  int ladderCountN[20]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0} ;
213  int ladderCountP[20]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0} ;
214  int data,pedestal,noise,channel,newchannel,ladder; char EastWest;
215  int maxChannel;
216  int nStrip =0;
217  int iLad = 0;
218  int iWaf = 0;
219  int chip = 0;
220  // SSD parameters independant from its configuration (half or full barrel)
221  // mConfig->getTotalNumberOfLadders()=20;
222  // mConfig->getNumberOfStrips()=768;
223  // mConfig->getNumberOfWafers()=16;
224  for (Int_t i=0;i<9;i++){
225  PedestalNTuple[i] = 0;
226  }
227  maxChannel=mConfig->getNumberOfStrips()*mConfig->getNumberOfWafers();
228  TDataSet *daq = GetDataSet("StDAQReader");
229  if (!daq) {
230  LOG_WARN << "Make : StDAQReader Dataset not found - Skipping the event" << endm;
231  return kStWarn;
232  }
233  StDAQReader *daqReader = (StDAQReader*) (daq->GetObject());
234  if (!daqReader) {
235  LOG_WARN << "Make : StDAQReader Object not found - Skipping the event" << endm;
236  return kStWarn;
237  }
238  if (!(daqReader->SSDPresent())) {
239  LOG_WARN << "Make : The SSD was not in the data stream - Skipping the event" << endm;
240  return kStWarn;
241  }
242  StSSDReader *stssdreader = daqReader->getSSDReader();
243  if(!stssdreader)
244  {
245  LOG_WARN << "Make : SSD reader not found - Skipping the event" << endm;
246  return kStWarn;
247  }
248  // creating the spa_strip and spa_ped_strip tables
249  St_spa_strip *spa_strip = (St_spa_strip *) GetData("spa_strip");
250  if(!spa_strip)
251  {
252  spa_strip = new St_spa_strip("spa_strip",100000);
253  AddData(spa_strip);
254  }
255  St_ssdPedStrip *ssdPedStrip = (St_ssdPedStrip *) GetData("ssdPedStrip");
256  if(!ssdPedStrip)
257  {
258  ssdPedStrip = new St_ssdPedStrip("ssdPedStrip",100000);
259  AddData(ssdPedStrip);
260  }
261  spa_strip_st out_strip;
262  ssdPedStrip_st out_ped_strip;
263  count=1;
264  // looping on the sides ...
265 static const char EW[4]={'E','W','p','n'};
266  for (id_side=0;id_side<2;id_side++)
267  {
268  EastWest=EW[id_side];
269  LOG_DEBUG <<" process "<<EW[id_side+2]<<"-side ladders "<< endm;
270  for (ladder=0;ladder<mConfig->getTotalNumberOfLadders();ladder++)
271  {
272  if (mConfig->getLadderIsActive(ladder+1)>0)
273  {
274  if ((ladder+1)== 4 || (ladder+1)== 6 || (ladder+1)==10 ||
275  (ladder+1)==11 || (ladder+1)==13 || (ladder+1)==15 ||
276  (ladder+1)==17 )
277  maxChannel=mConfig->getNumberOfStrips()*mConfig->getNumberOfWafers()-1;
278  else maxChannel=mConfig->getNumberOfStrips()*mConfig->getNumberOfWafers();
279  for (newchannel=0;newchannel<maxChannel;newchannel++)
280  {
281  if ((ladder+1)== 4 || (ladder+1)== 6 || (ladder+1)==10 ||
282  (ladder+1)==11 || (ladder+1)==13 || (ladder+1)==15 ||
283  (ladder+1)==17 )
284  channel= newchannel+1;
285  else channel=newchannel;
286  if(stssdreader->getSsdData(ladder+1,EastWest,newchannel,data,pedestal,noise)==0)
287  {
288  //We are looking at a physics run
289  // filling the out_strip structure...
290  //the strip id is coded as follow:
291  // id_strip=10000*(10*strip_number+id_side)+id_wafer
292  //strip number=1-mConfig->getNumberOfStrips()
293  //id_side=0 for p side, 1 for n side
294  if (data>0) {
295  if (id_side==1) my_channel=maxChannel-1-channel;
296  else my_channel=channel;
297  // the ssd mapping tables are inverted at the moment so we have to scan
298  // them to get the correct channel.
299  if (id_side==1) {
300  for (int kk=0;kk<maxChannel;kk++) {
301  if (ssd_ladder_mapN[kk]==channel) my_channel=kk;
302  }
303  }
304  else {
305  for (int kk=0;kk<maxChannel;kk++) {
306  if (ssd_ladder_mapP[kk]==channel) my_channel=kk;
307  }
308  }
309  assert(my_channel>=0);
310  strip_number=my_channel-(my_channel/mConfig->getNumberOfStrips())*mConfig->getNumberOfStrips()+1;
311  if (id_side==0)
312  id_wafer=7000+100*(mConfig->getNumberOfWafers()-(my_channel/mConfig->getNumberOfStrips()))+ladder+1;
313  else
314  id_wafer=7000+100*((my_channel/mConfig->getNumberOfStrips())+1)+ladder+1;
315  out_strip.id = count;
316  out_strip.adc_count = data;
317  out_strip.id_strip = 10000*(10*strip_number+id_side)+id_wafer;
318  out_strip.id_mchit[0] = 0;
319  out_strip.id_mchit[1] = 0;
320  out_strip.id_mchit[2] = 0;
321  out_strip.id_mchit[3] = 0;
322  out_strip.id_mchit[4] = 0;
323  spa_strip->AddAt(&out_strip);
324  if (id_side ==0) ladderCountP[ladder]++;
325  else ladderCountN[ladder]++;
326  count++;
327  } // end else for if data > 0
328  }
329  else // else for if(stssdreader->getSsdData(ladder,EastWest,channel,data,pedestal,noise)==0)
330  {
331  //We are looking at a pedestal run
332  // filling the out_strip structure...
333  //the strip id is coded as follow:
334  // id_strip=10000*(10*strip_number+id_side)+id_wafer
335  //strip number=1-mConfig->getNumberOfStrips()
336  //id_side=0 for p side, 1 for n side
337  if (pedestal>=0) {
338  // if (id_side==1) my_channel=maxChannel-1-channel;
339  // else my_channel=channel;
340  // the ssd mapping tables are inverted at the moment so we have to scan
341  // them to get the correct channel.
342  if (id_side==1) {
343  for (int kk=0;kk<maxChannel;kk++) {
344  if (ssd_ladder_mapN[kk]==channel) my_channel=kk;
345  }
346  // my_channel=ssd_ladder_mapN[channel];
347  }
348  else {
349  for (int kk=0;kk<maxChannel;kk++) {
350  if (ssd_ladder_mapP[kk]==channel) my_channel=kk;
351  }
352  // my_channel=ssd_ladder_mapP[channel];
353  }
354  strip_number=my_channel-(my_channel/mConfig->getNumberOfStrips())*mConfig->getNumberOfStrips()+1;
355  if (id_side==0)
356  id_wafer=7000+100*(mConfig->getNumberOfWafers()-(my_channel/mConfig->getNumberOfStrips()))+ladder+1;
357  else
358  id_wafer=7000+100*((my_channel/mConfig->getNumberOfStrips())+1)+ladder+1;
359  out_ped_strip.id = count;
360  out_ped_strip.id_strip = 10000*(10*strip_number+id_side)+id_wafer;
361  out_ped_strip.noise = noise;
362  out_ped_strip.pedestal = pedestal;
363  nStrip=(int)(out_ped_strip.id_strip/100000.);
364  iWaf = (int)((id_wafer - 7*1000)/100 - 1);
365  iLad = (int)(id_wafer - 7*1000 - (iWaf+1)*100 - 1);
366  nStrip=nStrip-1; // to have the good number of chip [0;95]
367  chip=(int)((nStrip+(768*(iWaf)))/128.);
368  ssdPedStrip->AddAt(&out_ped_strip);
369  if (id_side ==0) {
370  ladderCountP[ladder]++;
371  occupancy_wafer->Fill(2*iLad,iWaf,1);
372  occupancy_chip->Fill(2*iLad,chip,1);
373  noise_chip->Fill(2*iLad,chip,(out_ped_strip.noise/(16.)));
374  noise_wafer->Fill(2*iLad,iWaf,(out_ped_strip.noise/(16.)));
375  noise_chip_P->Fill(iLad,chip,(out_ped_strip.noise/(16.)));
376  pedestal_chip->Fill(2*iLad,chip,out_ped_strip.pedestal);
377  if(pedestal==0)ped_zero_ladP->Fill(iLad,iWaf);
378  if(pedestal==255)ped_high_ladP->Fill(iLad,iWaf);
379  if(noise==0)noise_zero_ladP->Fill(iLad,iWaf);
380  if(noise==255)noise_high_ladP->Fill(iLad,iWaf);
381 
382  //store the pedestal and noise if mPedOut is on
383  if(mPedOut){
384  PedestalNTuple[0]=0;
385  PedestalNTuple[1]=iLad;
386  PedestalNTuple[2]=iWaf;
387  PedestalNTuple[3]=nStrip+1;
388  PedestalNTuple[4]=out_ped_strip.pedestal;
389  PedestalNTuple[5]=out_ped_strip.noise;
390  PedestalNTuple[6]=(nStrip/128)+1;
391  PedestalNTuple[7]=id_wafer;
392  PedestalNTuple[8]=nStrip+(iWaf-1)*768;
393  pTuple->Fill(PedestalNTuple);
394  }
395  }
396  else {
397  ladderCountN[ladder]++;
398  occupancy_wafer->Fill((2*iLad)+1,iWaf,1);
399  occupancy_chip->Fill((2*iLad)+1,chip,1);
400  noise_chip->Fill((2*iLad)+1,chip,(out_ped_strip.noise/(16.)));
401  noise_wafer->Fill((2*iLad)+1,iWaf,(out_ped_strip.noise/(16.)));
402  noise_chip_N->Fill(iLad,chip,(out_ped_strip.noise/(16.)));
403  pedestal_chip->Fill((2*iLad)+1,chip,out_ped_strip.pedestal);
404  if(pedestal==0)ped_zero_ladN->Fill(iLad,iWaf);
405  if(pedestal==255)ped_high_ladN->Fill(iLad,iWaf);
406  if(noise==0)noise_zero_ladN->Fill(iLad,iWaf);
407  if(noise==255)noise_high_ladN->Fill(iLad,iWaf);
408  if(mPedOut){
409  PedestalNTuple[0]=1;
410  PedestalNTuple[1]=iLad;
411  PedestalNTuple[2]=iWaf;
412  PedestalNTuple[3]=nStrip-1;
413  PedestalNTuple[4]=out_ped_strip.pedestal;
414  PedestalNTuple[5]=out_ped_strip.noise;
415  PedestalNTuple[6]=(nStrip/128)+1;
416  PedestalNTuple[7]=id_wafer;
417  PedestalNTuple[8]=nStrip+(iWaf-1)*768;
418  pTuple->Fill(PedestalNTuple);
419  }
420  }
421  count++;
422  for (int i=0;i<20;i++)
423  {
424  occupancy->SetFillColor(2);
425  occupancy->Fill((2*i),(12288-ladderCountP[i]));
426  occupancy->Fill((2*i)+1,(12288-ladderCountN[i]));
427  }
428  //LOG_INFO << "Make()/ssdPedStrip->NRows= "<<ssdPedStrip->GetNRows()<<endm;
429  } // end if pedestal > 0
430  } // end if(stssdreader->getSsdData(ladder,EastWest,channel,data,pedestal,noise)==0)
431  } // end for (channel=0;channel<maxChannel;channel++)
432  } // end if (mConfig->getLadderIsActive(ladder)>0)
433  } // end for (ladder=0;ladder<mConfig->mConfig->getTotalNumberOfLadders();ladder++)
434  }
435  // end for (id_side=0;id_side<2;id_side++)
436 
437  LOG_DEBUG <<"Make()/Number of raw data in the SSD ";
438  LOG_DEBUG << "Make()/Active Ladders: ";
439  for (int i=0;i<mConfig->getTotalNumberOfLadders();i++){
440  if (mConfig->getLadderIsActive(i+1)>0) {
441  LOG_DEBUG.width(5);
442  LOG_DEBUG <<i+1<<" ";
443  } }
444  LOG_DEBUG <<endm;
445  LOG_DEBUG << "Make()/Counts (p-side): ";
446  for (int i=0;i<mConfig->getTotalNumberOfLadders();i++){
447  if (mConfig->getLadderIsActive(i+1)>0) {
448  LOG_DEBUG.width(5);
449  LOG_DEBUG <<ladderCountP[i]<<" ";
450  } }
451  LOG_DEBUG << endm;
452  LOG_DEBUG << "Make()/Counts (n-side): ";
453  for (int i=0;i<mConfig->getTotalNumberOfLadders();i++){
454  if (mConfig->getLadderIsActive(i+1)>0) {
455  LOG_DEBUG.width(5);
456  LOG_DEBUG <<ladderCountN[i]<<" ";
457  } }
458  LOG_DEBUG<<endm;
459  if((spa_strip->GetNRows()==0)&&(ssdPedStrip && ssdPedStrip->GetNRows()!=0)){
460  LOG_INFO << "Make()/ Read Pedestal & Noise"<< endm;
461  LOG_INFO << "Make()/ssdPedStrip->NRows= "<<ssdPedStrip->GetNRows()<<endm;
462  }
463  else
464  if((! ssdPedStrip || ssdPedStrip->GetNRows()==0) && (spa_strip->GetNRows()!=0)){
465  LOG_INFO << "Make()/ Read Signal from Physics Run"<< endm;
466  LOG_INFO << "Make()/ spa_strip->NRows= "<<spa_strip->GetNRows()<<endm;
467  }
468  return kStOK;
469 }
470 
471 //_____________________________________________________________________________
473 {
474  LOG_INFO << Form("Finish()") << endm;
475  if(mPedOut){
476  LOG_INFO << "Write Pedestal tuple"<< endm;
477  pFile->Write();
478  pFile->Close();
479  }
480  return kStOK;
481 }
482 //_____________________________________________________________________________
483 void StSsdDaqMaker::PrintConfiguration(Int_t runumber,ssdConfiguration_st *config)
484 {
485  switch(runumber){
486  case 4 : {
487  LOG_INFO <<"Configuration of ladders for run IV" <<endm;
488  break;
489  }
490  case 5 : {
491  LOG_INFO <<"Configuration of ladders for run V" << endm;
492  break;
493  }
494  case 6 : {
495  LOG_INFO <<"Configuration of ladders for run VI"<< endm;
496  break;
497  }
498  case 7 : {
499  LOG_INFO <<"Configuration of ladders for run VII" << endm;
500  break;
501  }
502  }
503  Int_t i =0;
504  Int_t totladderPresent =0;
505  LOG_INFO << "PrintLadderSummary:ladder id :";
506  for (i=1;i<=config->nMaxLadders;i++){
507  LOG_INFO.width(3);
508  LOG_INFO << i;
509  }
510  LOG_INFO <<endm;
511  LOG_INFO << "PrintLadderSummary:Active Ladders on sectors: ";
512  for (i=1;i<=config->nMaxLadders;i++){
513  LOG_INFO.width(3);
514  LOG_INFO <<mConfig->getLadderIsActive(i);
515  if(mConfig->getLadderIsActive(i)>0)totladderPresent++;
516 
517  }
518  LOG_INFO << endm;
519  LOG_INFO << "totLadderActive = "<<totladderPresent<<endm;
520 }
virtual void AddData(TDataSet *data, const char *dir=".data")
User methods.
Definition: StMaker.cxx:332
virtual Int_t Init()
Init - is a first method the top level StChain calls to initialize all its makers.
virtual Int_t Finish()
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
virtual ~StSsdDaqMaker()
This is StSsdDaqMaker destructor.
virtual Int_t Make()
Class to read SSD data from DAQ files.
Definition: StSsdDaqMaker.h:72
Definition: Stypes.h:41