StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StSvtOnlineSeqAdjSimMaker.cxx
1 /***************************************************************************
2  *
3  * $Id: StSvtOnlineSeqAdjSimMaker.cxx,v 1.16 2010/04/15 18:37:35 baumgart Exp $
4  *
5  * Author: Petr Chaloupka
6  ***************************************************************************
7  *
8  * Description: Simulates online seqence adjusting
9  *
10  ***************************************************************************/
11 
12 
13 #include "StSvtOnlineSeqAdjSimMaker.h"
14 #include "StSvtClassLibrary/StSvtHybridPixelsC.hh"
15 #include "StSvtClassLibrary/StSvtHybridPixelsD.hh"
16 #include "StSvtClassLibrary/StSvtHybridCollection.hh"
17 #include "StSvtClassLibrary/StSvtData.hh"
18 #include "StSvtClassLibrary/StSvtHybridBadAnodes.hh"
19 #include "StSvtClassLibrary/StSvtHybridData.hh"
20 #include "StSvtClassLibrary/StSvtConfig.hh"
21 #include "StSvtClassLibrary/StSvtDaq.hh"
22 #include "StSequence.hh"
23 #include "StMCTruth.h"
24 #include "StMessMgr.h"
25 #include "StSvtConversionTable.h"
26 
27 //__________________________________________________________________________________________________
29 
31 {
32  mRawData=NULL;
33  m8bitPixelColl=NULL;
34  mPixelColl=NULL;
35  mSvtBadAnodes=NULL;
36  mCurrentPixelData=NULL;
37  mCurrent8bitPixelData=NULL;
38  mPedOffsetAdjustment=2;
39 
40  //This is because of some Makers downd the chain
41  GetConfig();
42  SetRawData();
43 }
44 
45 //____________________________________________________________________________
46 StSvtOnlineSeqAdjSimMaker::~StSvtOnlineSeqAdjSimMaker()
47 {
48 }
49 
50 //____________________________________________________________________________
51 Int_t StSvtOnlineSeqAdjSimMaker::GetConfig()
52 {
53  mConfig=NULL;
54  St_DataSet *dataSet = NULL;
55  dataSet = GetDataSet("StSvtConfig");
56  if (!dataSet)
57  {
58  gMessMgr->Warning() << " No SvtConfig data set" << endm;
59  dataSet = new St_ObjectSet("StSvtConfig");
60  AddConst(dataSet);
61  }
62  else mConfig=((StSvtConfig*)(dataSet->GetObject()));
63 
64  if (!mConfig) {
65  gMessMgr->Warning() << "SvtConfig data set is empty- seting default full configuration" << endm;
66  mConfig=new StSvtConfig();
67  mConfig->setConfiguration("FULL");
68  mConfig->setNumberOfAnodes(240);
69  mConfig->setNumberOfTimeBins(128);
70  dataSet->SetObject(mConfig);
71  }
72 
73  return kStOk;
74 }
75 
76 //__________________________________________________________________________________________________
77 Int_t StSvtOnlineSeqAdjSimMaker::GetDaqParams()
78 {
79  mDaq=NULL;
80 
81  St_DataSet* dataSet;
82  dataSet = GetDataSet("StSvtDaq");
83  if (dataSet==NULL){
84  gMessMgr->Error()<<"BIG TROUBLE:No DAQ parameters from database!!!!"<<endm;
85  return kStFatal;
86  }
87 
88  mDaq = (StSvtDaq*)dataSet->GetObject();
89  if (mDaq==NULL){
90  gMessMgr->Error()<<"BIG TROUBLE:No DAQ parameters are empty!!!!"<<endm;
91  return kStFatal;
92  }
93 
94  SetNumberTBinsToClear(mDaq->getClearedTimeBins());
95  SetExtraPixelsBefore(mDaq->getPixelsBefore());
96  SetExtraPixelsAfter(mDaq->getPixelsAfter());
97  SetPedOffset(mDaq->getPedOffset());
98  Set_n_seq_lo(mDaq->getSeqLo());
99  Set_n_seq_hi(mDaq->getSeqHi());
100  Set_thresh_lo(mDaq->getThreshLo());
101  Set_thresh_hi(mDaq->getThreshHi());
102 
103  return kStOk;
104 }
105 
106 //____________________________________________________________________________
108 Int_t StSvtOnlineSeqAdjSimMaker::GetPixelData()
109 {
110  St_DataSet* dataSet=NULL;
111  dataSet = GetDataSet("StSvtPixelData");
112  if (!dataSet) {
113  gMessMgr->Error()<<"no StSvtPixelData dataset"<<endm;
114  return kStErr;
115  }
116 
117  mPixelColl=(StSvtData*)dataSet->GetObject();
118  if (!mPixelColl){
119  gMessMgr->Error()<<"StSvtPixelData is empty"<<endm;
120  return kStErr;
121  }
122 
123  dataSet=NULL;
124 
125  dataSet = GetDataSet("StSvt8bitPixelData");
126  if (!dataSet) {
127  gMessMgr->Error()<<"no StSvt8bitPixelData dataset"<<endm;
128  return kStErr;
129  }
130  m8bitPixelColl=(StSvtData*)dataSet->GetObject();
131  if (!m8bitPixelColl){
132  gMessMgr->Error()<<"StSvt8bitPixelData is empty"<<endm;
133  return kStErr;
134  }
135 
136  return kStOk;
137 }
138 
139 //____________________________________________________________________________
141 void StSvtOnlineSeqAdjSimMaker::GetBadAnodes()
142 {
143  St_DataSet *dataSet;
144 
145  dataSet = GetDataSet("StSvtBadAnodes");
146  if( !dataSet) {
147  gMessMgr->Warning() << "StSvtOnlineSeqAdjSimMaker: No Svt Bad Anodes data set" << endm;
148  return;
149  }
150 
151  mSvtBadAnodes = (StSvtHybridCollection*)(dataSet->GetObject());
152  if( !mSvtBadAnodes) {
153  gMessMgr->Warning() << "StSvtOnlineSeqAdjSimMaker: No Svt Bad Anodes data " << endm;
154  return;
155  }
156 
157  /*if (Debug())*/ gMessMgr->Info()<<"StSvtOnlineSeqAdjSimMaker:Svt Bad Anode list found"<<endm;
158 }
159 
160 
161 //____________________________________________________________________________
162 void StSvtOnlineSeqAdjSimMaker::SetRawData()
163 { //this makes new or replaces raw data for SvtDaqMaker
164 
165  St_ObjectSet* set=(St_ObjectSet*)GetDataSet("StSvtRawData");
166 
167  if (!set) {
168  set = new St_ObjectSet("StSvtRawData");
169  AddData(set);
170  }
171  mRawData = new StSvtData(mConfig->getConfiguration());
172  set->SetObject(mRawData);
173 }
174 
175 //____________________________________________________________________________
177 {
178  return StMaker::Init();
179 }
180 
181 //____________________________________________________________________________
184 {
185  if (Debug()) gMessMgr->Info()<<"StSvtOnlineSeqAdjSimMaker::InitRun"<<endm;
186  GetConfig();
187  GetBadAnodes();
188  Int_t res;
189  if ((res=GetDaqParams())!=kStOk) return res;
190 
191  gMessMgr->Info()<< " DAQ parameters used for simulation:"<<endm;
192  gMessMgr->Info() << " PedOffSet = "<<mPedOffset<<endm;
193  gMessMgr->Info() << " thresh_lo = "<<m_thresh_lo <<endm;
194  gMessMgr->Info() << " thresh_hi = "<<m_thresh_hi <<endm;
195  gMessMgr->Info() << " n_seq_lo = "<<m_n_seq_lo <<endm;
196  gMessMgr->Info() << " n_seq_hi = "<< m_n_seq_hi <<endm;
197 
198 
199  if (Debug()) gMessMgr->Info()<<"StSvtOnlineSeqAdjSimMaker::InitRun...END"<<endm;
200  return StMaker::InitRun(runumber);
201 }
202 
203 //____________________________________________________________________________
205 {
206  return kStOK;
207 }
208 
209 //____________________________________________________________________________
210 void StSvtOnlineSeqAdjSimMaker::Clear(const char*)
211 {
212  mRawData=NULL;
213  StMaker::Clear();
214 }
215 
216 //____________________________________________________________________________
217 void StSvtOnlineSeqAdjSimMaker::SetAdjParams(int thresh_lo,int n_seq_lo,int thresh_hi,int n_seq_hi)
218 {
219 m_thresh_lo =thresh_lo;
220 m_n_seq_lo =n_seq_lo;
221 m_thresh_hi =thresh_hi;
222 m_n_seq_hi =n_seq_hi;
223 }
224 
225 //____________________________________________________________________________
228 {
229  if (Debug()) gMessMgr->Info()<<"StSvtOnlineSeqAdjSimMaker::Make"<<endm;
230  SetRawData();
231  Int_t res;
232  if ((res=GetPixelData())!= kStOk) return res;
233 
234  for(int Barrel = 1;Barrel <= mPixelColl->getNumberOfBarrels();Barrel++) {
235  for (int Ladder = 1;Ladder <= mPixelColl->getNumberOfLadders(Barrel);Ladder++) {
236  for (int Wafer = 1;Wafer <= mPixelColl->getNumberOfWafers(Barrel);Wafer++) {
237  for( int Hybrid = 1;Hybrid <= mPixelColl->getNumberOfHybrids();Hybrid++){
238 
239  mCurrentIndex = mPixelColl->getHybridIndex(Barrel,Ladder,Wafer,Hybrid);
240  if( mCurrentIndex < 0) continue;
241  //cout<<"index:" <<mCurrentIndex<<endl;
242 
243  mCurrentPixelData = (StSvtHybridPixelsD*)mPixelColl->at(mCurrentIndex);
244  mCurrent8bitPixelData = (StSvtHybridPixelsC*)m8bitPixelColl->at(mCurrentIndex);
245 
246  if(!mCurrent8bitPixelData) {
247  mCurrent8bitPixelData = new StSvtHybridPixelsC(Barrel, Ladder, Wafer, Hybrid);
248  m8bitPixelColl->put_at(mCurrent8bitPixelData,mCurrentIndex);
249  }
250 
251  if(!mCurrentPixelData) { //no data from simulation Maker
252  mCurrent8bitPixelData->Reset();
253  }
254 
255  //No we have the pixel data sructures,turn it now into the real data
256  //ie. simulate the DAQ
257 
258  Conversion10to8bit();
259  ClearMask();
260  if (mKillBadAnodes) KillBadAnodes();
261  if (mNumberTBinsToClear>0) ClearFirstTbins();
262  RawAnodes();
263  SequenceSearch();
264  WriteMask();
265  FillRawData();
266 
267  }
268  }
269  }
270  }
271 
272  if (Debug()) gMessMgr->Info()<<"StSvtOnlineSeqAdjSimMaker::Make...END"<<endm;
273  return kStOK;
274 }
275 
276 //____________________________________________________________________________
277 void StSvtOnlineSeqAdjSimMaker::Conversion10to8bit()
278 {
279  double *fromArray=mCurrentPixelData->GetArray();
280  Char_t *toArray=mCurrent8bitPixelData->GetArray();
281 
282  for (int i=0;i<mCurrentPixelData->getTotalNumberOfPixels();i++)
283  {
284  double adc=fromArray[i];
285  if (adc<=0) adc=0.;
286  if (adc>=1023) adc=1023.;
287  unsigned int adc1=(unsigned int)adc; //conversion to 10 bits from double - ?is it "compiler" reliable?
288  toArray[i]=(Char_t)StSvt10to8ConversionTable[adc1]; //conversion to 8 bits
289  }
290 }
291 
292 //____________________________________________________________________________
293 void StSvtOnlineSeqAdjSimMaker::FillRawData()
294 {
295  StSvtHybridData *hybridData;
296  hybridData = new StSvtHybridData(mCurrentPixelData->getBarrelID(), mCurrentPixelData->getLadderID(), mCurrentPixelData->getWaferID(),mCurrentPixelData->getHybridID());
297  mRawData->put_at(hybridData,mCurrentIndex);
298 
299  int anodes=0; //number of anodes with some sequences
300 
301  //raw anodes need to checked if they are not flagged bad
302  StSvtHybridBadAnodes* badAnode =NULL;
303  if (mKillBadAnodes && mSvtBadAnodes) badAnode = (StSvtHybridBadAnodes*)mSvtBadAnodes->at(mCurrentIndex);
304 
305  Char_t *mAdcArray=mCurrent8bitPixelData->GetArray(); // array of [128*240]
306 
307  StSequence tmpSeq[128]; //buffer for sequences on one anode
308  StMCTruth tmpTru[128]; //buffer for truth on one anode
309  for (int ianode=0;ianode<240;ianode++)
310  {
311  int seqCount=0; //number of sequences on current anode
312 
313  //first check for raw anodes
314  int an=ianode+1;
315  if ((an==mDaq->getSavedBlackAnodes(0))||(an==mDaq->getSavedBlackAnodes(1))||(an==mDaq->getSavedBlackAnodes(2))||(an==mDaq->getSavedBlackAnodes(3)))
316  {
317  if (badAnode->isBadAnode(an)) continue; //don't write out zeros if bad
318  tmpSeq[0].startTimeBin =0;
319  tmpSeq[0].firstAdc=(unsigned char*)(mAdcArray+ianode*128);
320  tmpSeq[0].length = 128;
321  tmpTru[0] = 0;
322  seqCount=1;
323  hybridData->setListSequences(anodes, an, seqCount, tmpSeq);
324  anodes++;
325  continue;
326  }
327 
328  int pixCount=0;
329  StMCPivotTruth pivo;
330  for(int tim = 0; tim <= 128; tim++)
331  {//loop over time bins in one anode
332  unsigned char adc;
333  if (tim==128) adc=0; // make an artificial end of time sequence
334  else adc= (unsigned char)mAdcArray[ianode*128 + tim];
335 
336  if (adc>0)
337  {
338  StMCTruth tru =mCurrentPixelData->getTrackId(ianode*128 + tim);
339  if (pixCount==0){ //starting new sequence
340  pivo.Reset();
341  tmpSeq[seqCount].startTimeBin = tim;
342  tmpSeq[seqCount].firstAdc=(unsigned char*)(mAdcArray+ianode*128 + tim);
343  }
344  if(int(tru)) pivo.Add(tru,adc);
345  pixCount++;
346  }
347  else
348  {
349  if(pixCount>0){//end of sequence
350  tmpSeq[seqCount].length = pixCount;
351  tmpTru[seqCount] = pivo.Get();
352  seqCount++;
353  pixCount=0;
354  }
355  }
356 
357 
358  }
359 
360  if(seqCount>0){ //save found sequences
361  //cout<<"found sequences:"<<seqCount<<endl;
362  hybridData->setListSequences(anodes,ianode+1, seqCount, tmpSeq);
363  hybridData->setListTruth (anodes,ianode+1, seqCount, tmpTru);
364  anodes++;
365  }
366 
367 
368  }
369 
370  hybridData->setAnodeList();
371 }
372 
373 //____________________________________________________________________________
374 void StSvtOnlineSeqAdjSimMaker::WriteSequence(int anode,int begins, int ends, int NumOfHigh)
375 {
376  //check the proper size
377  //cout<<"anode:"<<anode<<"start:"<<begins<<"ends:"<<ends<<"hi:"<<NumOfHigh<<endl;
378 
379  if (NumOfHigh<=m_n_seq_hi) return;
380  if ((ends-begins+1)<=m_n_seq_lo) return;
381 
382  //extra anodes
383  begins= begins-mExtraBefore;
384  if (begins<0) begins=0;
385 
386  ends=ends+mExtraAfter;
387  if (ends>127) ends=127;
388 
389  for (int i=begins;i<=ends;i++) mMask[anode*128 + i]=kTRUE;
390 }
391 
392 //____________________________________________________________________________
393 void StSvtOnlineSeqAdjSimMaker::SequenceSearch()
394 {
395  unsigned char adc;
396  int loCount;
397  int hiCount;
398  int SeqBegins=0;
399 
400  int HiTresh=mPedOffset+m_thresh_hi+mPedOffsetAdjustment;
401  int LoTresh=mPedOffset+m_thresh_lo+mPedOffsetAdjustment;
402 
403 
404  Char_t *mAdcArray=mCurrent8bitPixelData->GetArray(); // array of [128*240]
405 
406  for(int an = 0; an < 240; an++){
407  //get ready for new anode
408  loCount=0;hiCount=0; //just for safety
409 
410  for(int tim = 0; tim <= 128; tim++)
411  {//loop over time bins in one anode
412 
413  if (tim==128) adc=0; // make an artificial end of time sequence
414  else adc=(unsigned char)mAdcArray[an*128 + tim];
415 
416  if (adc>HiTresh) hiCount++;
417 
418  if (adc>LoTresh)
419  { //inside of sequence or at the beginning
420  if (loCount==0) SeqBegins=tim; //it is the beginning of the sequence
421  loCount++;
422  }
423  else
424  { //ouside or at the end of the sequence
425  if(loCount!=0) //end of sequence
426  {
427  WriteSequence(an,SeqBegins,tim-1,hiCount);
428  loCount=0;hiCount=0;
429  }
430  }
431 
432 
433  }
434  }
435 
436 }
437 
438 
439 
440 //____________________________________________________________________________
441 void StSvtOnlineSeqAdjSimMaker::KillBadAnodes()
442 {
443 
444  if (!mSvtBadAnodes){
445  cout<<"Warning: cannot simulate bad anodes in online sequence adjusting - no anode list"<<endl;
446  return;
447  }
448  StSvtHybridBadAnodes* BadAnode = (StSvtHybridBadAnodes*)mSvtBadAnodes->at(mCurrentIndex);
449  if (!BadAnode) return;
450 
451  Char_t *mAdcArray=mCurrent8bitPixelData->GetArray(); // array of [128*240]
452 
453  for (int an=0;an<240;an++)
454  if (BadAnode->isBadAnode(an+1))
455  { //now I've found bad anode and I'm going to delete it
456  for(int tim = 0; tim < 128; tim++) mAdcArray[an*128 + tim]=0;
457  }
458 
459 }
460 
461 
462 //____________________________________________________________________________
463 void StSvtOnlineSeqAdjSimMaker::RawAnodes()
464 {
465 
466  for (int i=0 ; i<4 ; i++)
467  {
468  int an=mDaq->getSavedBlackAnodes(i);
469  if ((an<=0)||(an>240)) continue;
470  for(int tb = 0; tb < 128; tb++) mMask[an*128 + tb]=kTRUE;
471  }
472 
473 }
474 
475 //____________________________________________________________________________
476 void StSvtOnlineSeqAdjSimMaker::ClearMask()
477 {
478  memset(mMask,0,128*240*sizeof(mMask[0]));
479 }
480 
481 //____________________________________________________________________________
482 void StSvtOnlineSeqAdjSimMaker::ClearFirstTbins()
483 {
484  Char_t *mAdcArray=mCurrent8bitPixelData->GetArray(); // array of [128*240]
485 
486  for(int tim = 0; tim < mNumberTBinsToClear; tim++){
487  for(int an = 0; an < 240; an++){
488  mAdcArray[an*128 + tim]=0;
489 
490  }
491  }
492 }
493 
494 //____________________________________________________________________________
495 void StSvtOnlineSeqAdjSimMaker::WriteMask()
496 {
497  Char_t *mAdcArray=mCurrent8bitPixelData->GetArray(); // array of [128*240]
498 
499  for(int tim = 0; tim < 128; tim++){
500  for(int an = 0; an < 240; an++){
501  if (mMask[an*128 + tim]==kFALSE) mAdcArray[an*128 + tim]=0;
502  }
503  }
504 }
505 
506 
virtual Int_t InitRun(int runumber)
Reads run dependent data from the database.
virtual void AddData(TDataSet *data, const char *dir=".data")
User methods.
Definition: StMaker.cxx:332
virtual void Clear(Option_t *option="")
User defined functions.
Definition: StMaker.cxx:634
virtual void SetObject(TObject *obj)
The depricated method (left here for the sake of the backward compatibility)
Definition: TDataSet.cxx:480
Simulates proceses in the DAQ for SVT Slow Simulator: 10 to 8 bit conversion, killing of bad anodes...
virtual Int_t Init()
inherited maker routines
virtual void SetObject(TObject *obj)
The depricated method (left here for the sake of the backward compatibility)
Definition: TObjectSet.h:59
virtual TObject * GetObject() const
The depricated method (left here for the sake of the backward compatibility)
Definition: TDataSet.cxx:428
Definition: Stypes.h:40
virtual Int_t Make()
Does all data adjusting. At the end creates data in the raw format - the same format as the real data...
Definition: Stypes.h:44
Definition: Stypes.h:41