StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Stl3RawReaderMaker.cxx
1 //*-- Author : Dominik Flierl
2 // Christof Struck
3 //
5 // //
6 // This Maker reads l3 data as they come with the raw data stream //
7 // from the experiment and fills them into StEvent //
8 // //
10 //
11 // $Id: Stl3RawReaderMaker.cxx,v 1.21 2017/04/26 21:15:47 perev Exp $
12 //
13 // $Log: Stl3RawReaderMaker.cxx,v $
14 // Revision 1.21 2017/04/26 21:15:47 perev
15 // Hide m_Debug
16 //
17 // Revision 1.20 2008/01/23 19:10:33 fine
18 // Fix the lost L3_Reader class definition
19 //
20 // Revision 1.19 2007/04/28 17:56:24 perev
21 // Redundant StChain.h removed
22 //
23 // Revision 1.18 2004/08/09 18:25:21 kollegge
24 // Added tpc detector id when setting the number of points since StEvent can now store the number of points for different detectors.
25 //
26 // Revision 1.17 2004/03/26 11:25:23 kollegge
27 // Added another quality check of raw data quality to prevent inconsistent information in StL3EventSummary, fixes bug http://www.star.bnl.gov/rt2/Ticket/Display.html?id=359
28 //
29 // Revision 1.16 2004/03/26 00:31:59 dietel
30 // Check L3_SUMD, fixes http://www.star.bnl.gov/rt2/Ticket/Display.html?id=357
31 //
32 // Revision 1.15 2004/01/23 23:14:06 kollegge
33 // Crashs of the l3 online analysis code so far stopped the chain. Changed this to a warning message.
34 //
35 // Revision 1.14 2002/05/16 02:39:13 struck
36 // switch reco/embedding mode (m_Mode=0/1).
37 // Embedding mode skips L3 biased events (return kStErr).
38 // Reco mode fills StEvent as before.
39 //
40 // Revision 1.13 2002/02/27 20:17:58 struck
41 // adding globalTrack->setEncodedMethod() to mark tracks as l3 tracks in StEvent
42 //
43 // Revision 1.12 2002/02/26 21:40:38 struck
44 // move GetDataSet("StDAQReader") from Init() to Make(), solves a seg. fault in current dev release
45 //
46 // Revision 1.11 2002/02/13 22:36:32 struck
47 // major code clean-up for Stl3RawReaderMaker, first version of Stl3CounterMaker
48 //
49 // Revision 1.10 2001/11/14 23:30:56 struck
50 // major update: set 'unbiased'-flag, correct bugs in StGlobalTrack-filling
51 //
52 // Revision 1.9 2001/09/27 03:49:52 struck
53 // actual no. of gl3s handled flexible, max no. of gl3s and algorithms now global define-statements
54 //
55 // Revision 1.8 2001/09/25 01:42:51 struck
56 // cs: l3 vertex now put into StL3Trigger
57 //
58 // Revision 1.7 2001/08/29 20:30:40 struck
59 // and don't forget to delete your array ;-)
60 //
61 // Revision 1.6 2001/08/29 20:24:49 struck
62 // makes Solaris compiler happy
63 //
64 // Revision 1.5 2001/08/20 22:32:00 struck
65 // first version filling L3 counters and algorithm info into StEvent
66 //
67 //
68 
69 #include "Stl3RawReaderMaker.h"
70 #include "St_DataSetIter.h"
71 #include "StMessMgr.h"
72 #include "StDAQMaker/StDAQReader.h"
73 #include "StDaqLib/L3/L3_Reader.hh"
74 #include "tables/St_l3RunSummary_Table.h"
75 #include "tables/St_l3AlgorithmInfo_Table.h"
76 #include "StEventTypes.h"
77 #include "StEnumerations.h"
78 #include "Rtypes.h"
79 #include "TMath.h"
80 
81 #include "St_l3_Coordinate_Transformer.h"
82 #include "St_l3_Coordinates.h"
83 
84 
85 #define gufld gufld_
86 extern "C" {void gufld(Float_t *, Float_t *);}
87 
88 
89 ClassImp(Stl3RawReaderMaker)
90 
91 //_____________________________________________________________________________
92 Stl3RawReaderMaker::Stl3RawReaderMaker(const char *name):StMaker(name){
93  // l3RawReader constructor
94 }
95 
96 //_____________________________________________________________________________
97 Stl3RawReaderMaker::~Stl3RawReaderMaker(){
98 }
99 
100 
101 //_____________________________________________________________________________
102 Int_t Stl3RawReaderMaker::Init(){
103  // Init - is a first method the top level StChain calls to initialize all its makers
104 
105  mL3On = kFALSE;
106 
107  //SetDebug(1);
108 
109  // reset database pointer
110  mDbSet = 0;
111 
112  mNumberOfGl3Nodes = 0;
113  mNumberOfAlgorithms = 0;
114  mEventCounter = 0;
115 
116  // reset counter
117  for (int i=0; i<MaxNumberOfGl3Nodes; i++) {
118  mGlobalCounter[i].nProcessed = 0;
119  mGlobalCounter[i].nReconstructed = 0;
120  for (int j=0; j<MaxNumberOfAlgorithms; j++) {
121  mAlgorithmCounter[i][j].algId = 0;
122  mAlgorithmCounter[i][j].nProcessed = 0;
123  mAlgorithmCounter[i][j].nAccept = 0;
124  mAlgorithmCounter[i][j].nBuild = 0;
125  }
126  }
127 
128  // return
129  return StMaker::Init();
130 }
131 
132 
133 //_____________________________________________________________________________
135 {
136  //
137  // Make - this method is called in loop for each event
138  //
139 
140  // here we start
141  gMessMgr->Info("Stl3RawReaderMaker: Now we start l3RawReader Maker.");
142 
143 
144  // get the l3 daqreader
145  // Make Connection to raw data
146  DAQReaderSet = GetDataSet("StDAQReader");
147  if (!DAQReaderSet) {
148  gMessMgr->Error() << "Stl3RawReaderMaker::Make(): no DaqReader found!" << endm;
149  return kStWarn;
150  }
151 
152  StDAQReader *daqReader = (StDAQReader*)(DAQReaderSet->GetObject()) ;
153  if (daqReader) {
154  ml3reader = daqReader->getL3Reader();
155 
156  if (ml3reader) {
157 
158  // set flag L3 ON
159  mL3On = kTRUE;
160 
161  // debug output
162  if (Debug()) {
163  int sec = 23;
164  if (ml3reader->getGlobalTrackReader())
165  gMessMgr->Info() << ml3reader->getGlobalTrackReader()->getNumberOfTracks()
166  << " global tracks found." << endm;
167  if (ml3reader->getSl3ClusterReader(sec))
168  gMessMgr->Info() << ml3reader->getSl3ClusterReader(sec)->getNumberOfClusters()
169  << " sl3 clusters found in sec " << sec << endm;
170  if (ml3reader->getSl3TrackReader(sec))
171  gMessMgr->Info() << ml3reader->getSl3TrackReader(sec)->getNumberOfTracks()
172  << " sl3 tracks found in sec " << sec << endm;
173  if (ml3reader->getI960ClusterReader(sec))
174  gMessMgr->Info() << ml3reader->getI960ClusterReader(sec)->getNumberOfClusters()
175  << " i960 clusters found in sec " << sec << endm;
176  }
177 
178  switch (m_Mode) {
179  case 0:
180  // standard mode: fill StEvent
181  gMessMgr->Info("Stl3RawReaderMaker: Fill StEvent.");
182  if (fillStEvent() != 0) {
183  gMessMgr->Error("Stl3RawReaderMaker: problems filling l3 into StEvent.");
184  return kStErr;
185  }
186  break;
187  case 1:
188  // embedding mode: skip biased events,
189  // return kStErr (or better: kStSkip)
190  gMessMgr->Info("Stl3RawReaderMaker: Embedding Mode, check L3 bias.");
191  return checkL3Bias();
192  break;
193  default:
194  // unknown mode: do nothing!
195  gMessMgr->Error("Stl3RawReaderMaker: Unknown mode, return kStWarn");
196  return kStWarn;
197  }
198 
199  } // if (ml3reader)
200 
201  else {
202  // if L3 is on for this run but no data found
203  // warning: should happen only if we crashed during the run.
204  // new since FY04 run: DAQ doesn't send us all events
205  // so this is normal/expected behaviour
206  if (mL3On) {
207  gMessMgr->Warning("Stl3RawReaderMaker: L3 is ON, but no l3 data found in this event. Either DAQ has not send us this event or we crashed during online analysis.");
208  return kStWarn;
209  }
210  // if L3 is off so far
211  // it may be switched off for this run
212  // so don't crash the chain
213  else {
214  gMessMgr->Warning("Stl3RawReaderMaker: no l3 data found.");
215  // standard mode: return warning
216  if (m_Mode!=1) return kStWarn;
217  // embedding mode: return ok, since there is obviously no bias
218  else return kStOk;
219  }
220  }
221 
222  } // if (daqReader)
223 
224  // go home
225  return kStOk;
226 }
227 
228 //_____________________________________________________________________________
229 
230 Int_t Stl3RawReaderMaker::checkL3Bias()
231 {
232 
233  // get Gl3AlgorithmReader
234  Gl3AlgorithmReader* gl3Reader = ml3reader->getGl3AlgorithmReader();
235  if (!gl3Reader) {
236  gMessMgr->Error("Stl3RawReaderMaker: L3 is ON, but L3 summary data is missing!");
237  gMessMgr->Error("Stl3RawReaderMaker: ==> This event is bad, skip it!");
238  return kStErr;
239  }
240 
241  // check data
242  // L3_summary.on==0 indicates that the event crashed in
243  // L3 online and raw data contains no valid information for us
244  if (ml3reader->getL3_Summary()->on == 0) {
245  gMessMgr->Warning("Stl3RawReaderMaker: L3 crashed online on this event, no usefull information.");
246  gMessMgr->Warning("Stl3RawReaderMaker: ==> This event is WEIRD, but definitely unbiased!");
247  return kStOk;
248  }
249 
250  // number of algorithms
251  mNumberOfAlgorithms = gl3Reader->getNumberofAlgorithms();
252 
253  Algorithm_Data* algData = gl3Reader->getAlgorithmData();
254 
255  // default return value: kStErr == skip event
256  // should be kStSkip
257  int returnValue = kStErr;
258 
259  //
260  // check Pass-Thru Algorithm for hadronic triggers
261  // that is either TWReject && triggerWord=0x1...
262  // or TRUE && triggerWord=0x1...
263  //
264  // __STATUS as of NOV. 2001__
265  //
266  for (int i=0; i<mNumberOfAlgorithms; i++) {
267  if ( (algData[i].algId==10 || algData[i].algId==1)
268  && ml3reader->getL3_P()->trg_word>0x0fff
269  && ml3reader->getL3_P()->trg_word<0x2000
270  && algData[i].build==1 ) {
271 
272  if (Debug()) {
273  printf("L0 trigger word: %x\n", ml3reader->getL3_P()->trg_word);
274  printf("pass-thru algorithm Id: %i\n", algData[i].algId);
275  }
276  // this event is unbiased ==>> take it!
277  returnValue = kStOk;
278  }
279  }
280 
281  return returnValue;
282 }
283 
284 //_____________________________________________________________________________
285 
286 Int_t Stl3RawReaderMaker::fillStEvent()
287 {
288 
289  // get StEvent if not return
290  mStEvent = (StEvent *) GetInputDS("StEvent");
291  if (!mStEvent) {
292  gMessMgr->Error("Stl3RawReaderMaker: No StEvent");
293  return 1;
294  }
295 
296  // create StL3Trigger and connect it
297  myStL3Trigger = new StL3Trigger() ;
298  if (!myStL3Trigger) {
299  gMessMgr->Error("Stl3RawReaderMaker: No Stl3Trigger.");
300  return 1;
301  }
302  mStEvent->setL3Trigger(myStL3Trigger);
303 
304 
305  // check data
306  if (!ml3reader->getL3_SUMD()) {
307  gMessMgr->Warning("Stl3RawReaderMaker: No L3_SUMD bank.");
308  myStL3Trigger->setL3EventSummary(NULL);
309  return 0;
310  }
311 
312  if (!ml3reader->getL3_Summary()) {
313  gMessMgr->Warning("Stl3RawReaderMaker: No l3 summary.");
314  myStL3Trigger->setL3EventSummary(NULL);
315  return 0;
316  }
317 
318  // L3_summary.on==0 indicates that the event crashed in
319  // L3 online and raw data contains no valid information for us
320  if (ml3reader->getL3_Summary()->on == 0) {
321  gMessMgr->Warning("Stl3RawReaderMaker: L3 crashed online on this event, no usefull information.");
322  myStL3Trigger->setL3EventSummary(NULL);
323  return 0;
324  }
325 
326  // create StL3EventSummary
327  StL3EventSummary* myEventSummary = new StL3EventSummary(ml3reader->getL3_SUMD());
328  if (!myEventSummary) {
329  gMessMgr->Error("Stl3RawReaderMaker: No Stl3EventSummary.");
330  return 1;
331  }
332 
333  // connect StL3EventSummary to StL3Trigger
334  myStL3Trigger->setL3EventSummary(myEventSummary);
335 
336 
337  // get Gl3AlgorithmReader
338  Gl3AlgorithmReader* gl3Reader = ml3reader->getGl3AlgorithmReader();
339  if (!gl3Reader) {
340  gMessMgr->Error("Stl3RawReaderMaker: L3 is ON, but L3 summary data is missing!");
341  return 1;
342  }
343 
344  // get database info
345  mDbSet = GetDataBase("RunLog/l3");
346  // if we didn't get the db use maximum values as default
347  if (mDbSet) {
348  // get number of gl3 nodes
349  TTable* l3runSummaryTable = (TTable* )mDbSet->Find("l3RunSummary");
350  // database can be there, but now tables :-(
351  if (l3runSummaryTable) {
352  l3RunSummary_st* data = (l3RunSummary_st* )l3runSummaryTable->GetArray();
353  mNumberOfGl3Nodes = data->nGl3Nodes;
354  if (Debug()) {
355  gMessMgr->Info() << "database: runNumber = " << data->runNumber << endm;
356  gMessMgr->Info() << "database: nGl3Nodes = " << data->nGl3Nodes << endm;
357  }
358  // check limit
359  if (mNumberOfGl3Nodes>MaxNumberOfGl3Nodes) {
360  gMessMgr->Error() << "Stl3RawReaderMaker: number of gl3 nodes too high: db -> " << mNumberOfGl3Nodes
361  << " max -> " << MaxNumberOfGl3Nodes << endm;
362  return 1;
363  }
364  }
365  else {
366  mNumberOfGl3Nodes = MaxNumberOfGl3Nodes;
367  gMessMgr->Warning("Stl3RawReaderMaker: no table entry for this run in 'l3RunSummary' found!");
368  gMessMgr->Warning("Stl3RawReaderMaker: using default values.");
369  }
370  }
371  else {
372  mNumberOfGl3Nodes = MaxNumberOfGl3Nodes;
373  gMessMgr->Warning("Stl3RawReaderMaker: no database 'Runlog/l3' found!");
374  gMessMgr->Warning("Stl3RawReaderMaker: using default values.");
375  }
376 
377  // gl3 node which built this event
378  int gl3Id = ml3reader->getGl3Id();
379 
380  // number of algorithms
381  mNumberOfAlgorithms = gl3Reader->getNumberofAlgorithms();
382 
383  // check number of algorithms and gl3 node id
384  if (mNumberOfAlgorithms >= MaxNumberOfAlgorithms) {
385  gMessMgr->Error() << "Stl3RawReaderMaker: number of algorithms exceeds limit: found "
386  << mNumberOfAlgorithms << ", limit " << MaxNumberOfAlgorithms
387  << endm;
388  return 1;
389  }
390  if (gl3Id<=0 || gl3Id>= MaxNumberOfGl3Nodes) {
391  gMessMgr->Error() << "Stl3RawReaderMaker: gl3 id exceeds limit: found "
392  << gl3Id << ", limit " << MaxNumberOfGl3Nodes
393  << endm;
394  return 1;
395  }
396 
397  //--------------------------
398  // now do the bookkeeping
399  //--------------------------
400 
401  // increase the event counter
402  mEventCounter++;
403 
404  // first: global counter
405  mGlobalCounter[gl3Id-1].nProcessed = gl3Reader->getNumberOfProcessedEvents();
406  mGlobalCounter[gl3Id-1].nReconstructed = gl3Reader->getNumberOfReconstructedEvents();
407  // second: algorithm counter
408  Algorithm_Data* algData = gl3Reader->getAlgorithmData();
409  for (int i=0; i<mNumberOfAlgorithms; i++) {
410  mAlgorithmCounter[gl3Id-1][i].algId = algData[i].algId;
411  mAlgorithmCounter[gl3Id-1][i].nProcessed = algData[i].nProcessed;
412  mAlgorithmCounter[gl3Id-1][i].nAccept = algData[i].nAccept;
413  mAlgorithmCounter[gl3Id-1][i].nBuild = algData[i].nBuild;
414  }
415 
416  // get total counters
417  GlobalCounter totalCounter;
418  // make solaris happy
419  const int maxAlg = mNumberOfAlgorithms;
420  AlgorithmCounter* totalAlgCounter = new AlgorithmCounter[maxAlg];
421  //AlgorithmCounter totalAlgCounter[maxAlg];
422 
423 
424  totalCounter.nProcessed = 0;
425  totalCounter.nReconstructed = 0;
426  for (int i=0; i<mNumberOfAlgorithms; i++) {
427  totalAlgCounter[i].algId = 0;
428  totalAlgCounter[i].nProcessed = 0;
429  totalAlgCounter[i].nAccept = 0;
430  totalAlgCounter[i].nBuild = 0;
431  }
432 
433  for (int i=0; i<mNumberOfGl3Nodes; i++) {
434  // if one gl3 node wasn't seen so far,
435  // the counters are still undefined
436  if (mGlobalCounter[i].nProcessed==0 && mEventCounter<(2*mNumberOfGl3Nodes)) {
437  totalCounter.nProcessed = -1;
438  totalCounter.nReconstructed = -1;
439  for (int j=0; j<mNumberOfAlgorithms; j++) {
440  totalAlgCounter[j].algId = mAlgorithmCounter[gl3Id-1][j].algId;
441  totalAlgCounter[j].nProcessed = -1;
442  totalAlgCounter[j].nAccept = -1;
443  totalAlgCounter[j].nBuild = -1;
444  }
445  break;
446  }
447  // summ-up the counters of all gl3 nodes
448  // we have seen so far after #(mNumberOfGl3Nodes*2) events
449  totalCounter.nProcessed += mGlobalCounter[i].nProcessed;
450  totalCounter.nReconstructed += mGlobalCounter[i].nReconstructed;
451  for (int k=0; k<mNumberOfAlgorithms; k++) {
452  totalAlgCounter[k].algId = mAlgorithmCounter[gl3Id-1][k].algId;
453  totalAlgCounter[k].nProcessed += mAlgorithmCounter[i][k].nProcessed;
454  totalAlgCounter[k].nAccept += mAlgorithmCounter[i][k].nAccept;
455  totalAlgCounter[k].nBuild += mAlgorithmCounter[i][k].nBuild;
456  }
457  }
458 
459  //--------------------------
460  // fill StL3EventSummary
461  //--------------------------
462 
463  // fill counters and number of tracks
464  myEventSummary->setCounters(totalCounter.nProcessed, totalCounter.nReconstructed);
465  myEventSummary->setNumberOfTracks(ml3reader->getL3_Summary()->nTracks);
466  myEventSummary->setL0TriggerWord(ml3reader->getL3_P()->trg_word);
467 
468  //--------------------------
469  // fill StL3AlgorithmInfo
470  //--------------------------
471  // get database info if available
472  l3AlgorithmInfo_st* dbAlgInfo = 0;
473  if (mDbSet) {
474  TTable* l3algorithmInfoTable = (TTable* )mDbSet->Find("l3AlgorithmInfo");
475  if (l3algorithmInfoTable) {
476  dbAlgInfo = (l3AlgorithmInfo_st* )l3algorithmInfoTable->GetArray();
477  if (Debug()) {
478  for (int i=0; i<l3algorithmInfoTable->GetNRows(); i++) {
479  cout << " run \tidxAlg\talgId\tpreScale\tpostScale" << endl;
480  cout << dbAlgInfo[i].runNumber
481  << "\t" << dbAlgInfo[i].idxAlg
482  << "\t" << dbAlgInfo[i].algId
483  << "\t" << dbAlgInfo[i].preScale
484  << "\t" << dbAlgInfo[i].postScale << endl;
485  cout <<" GI: " << dbAlgInfo[i].GI1 << " " << dbAlgInfo[i].GI2 << " " << dbAlgInfo[i].GI3
486  << " " << dbAlgInfo[i].GI4 << " " << dbAlgInfo[i].GI5 << endl;
487  cout << " GF: " << dbAlgInfo[i].GF1 << " " << dbAlgInfo[i].GF2 << " " << dbAlgInfo[i].GF3
488  << " " << dbAlgInfo[i].GF4 << " " << dbAlgInfo[i].GF5 << endl;
489  cout << "------------" << endl;
490  }
491  }
492  // check entries, just to be _really_ safe (paranoia;-)
493  if (l3algorithmInfoTable->GetNRows()!=mNumberOfAlgorithms) {
494  gMessMgr->Warning() << "Stl3RawReaderMaker: Warning: database entries don't match raw data!" << endm;
495  gMessMgr->Warning() << "Stl3RawReaderMaker: db nAlgorithms: " << l3algorithmInfoTable->GetNRows()
496  << ", raw data: " << mNumberOfAlgorithms
497  << endm;
498  gMessMgr->Warning() << "Stl3RawReaderMaker: Skip database info for this event." << endm;
499  dbAlgInfo = 0;
500  }
501  }
502  else {
503  gMessMgr->Warning("Stl3RawReaderMaker: No entry for this run found in table 'l3algorithmInfo'.");
504  }
505  }
506  // now fill algorithm info in StEvent
507  for (int i=0; i<mNumberOfAlgorithms; i++) {
508  StL3AlgorithmInfo* myL3AlgorithmInfo = new StL3AlgorithmInfo(&algData[i]);
509  if (!myL3AlgorithmInfo) {
510  gMessMgr->Error("Stl3RawReaderMaker: No StL3AlgorithmInfo.");
511  return 1;
512  }
513  myL3AlgorithmInfo->setCounters(totalAlgCounter[i].nProcessed,
514  totalAlgCounter[i].nAccept,
515  totalAlgCounter[i].nBuild);
516 
517  //
518  // check Pass-Thru Algorithm for hadronic triggers
519  // that is either TWReject && triggerWord=0x1...
520  // or TRUE && triggerWord=0x1...
521  // __status as of Nov. 2001__
522  //
523  int passThruId = 0;
524  if ( (algData[i].algId==10 || algData[i].algId==1)
525  && ml3reader->getL3_P()->trg_word>0x0fff
526  && ml3reader->getL3_P()->trg_word<0x2000
527  && algData[i].build==1 ) {
528 
529  passThruId = algData[i].algId;
530  if (Debug()) {
531  cout << "pass-thru algorithm Id: " << passThruId << endl;
532  }
533 
534  // set unbiased flag in event summary
535  myEventSummary->setUnbiasedTrigger();
536  }
537 
538  // get algorithm info from database
539  // and check algId (paranoia again?)
540  if (dbAlgInfo && dbAlgInfo[i].algId==totalAlgCounter[i].algId) {
541  myL3AlgorithmInfo->setPreScale(dbAlgInfo[i].preScale);
542  myL3AlgorithmInfo->setPostScale(dbAlgInfo[i].postScale);
543  int intPara[IntParameterSize] = {dbAlgInfo[i].GI1, dbAlgInfo[i].GI2,
544  dbAlgInfo[i].GI3, dbAlgInfo[i].GI4,
545  dbAlgInfo[i].GI5};
546  float floatPara[FloatParameterSize] = {dbAlgInfo[i].GF1, dbAlgInfo[i].GF2,
547  dbAlgInfo[i].GF3, dbAlgInfo[i].GF4,
548  dbAlgInfo[i].GF5};
549  myL3AlgorithmInfo->setParameters(intPara, floatPara);
550 
551  if (dbAlgInfo[i].algId==passThruId)
552  myEventSummary->setUnbiasedTriggerPreScale(dbAlgInfo[i].preScale);
553 
554  } // if dbAlgInfo
555  myEventSummary->addAlgorithm(myL3AlgorithmInfo);
556 
557  } // for mNumberOfAlgorithms
558 
559  // delete counter array
560  delete [] totalAlgCounter;
561 
562 
563  // call filling routines
564  // global tracks and vertex
565  if (ml3reader->getGlobalTrackReader()) {
566 
567  // vertex position
568  StPrimaryVertex* myL3Vertex = new StPrimaryVertex;
569  vertex vert = ml3reader->getGlobalTrackReader()->getVertex();
570  StThreeVectorF* pos = new StThreeVectorF(vert.x, vert.y, vert.z);
571  myL3Vertex->setPosition(*pos);
572  myStL3Trigger->addPrimaryVertex(myL3Vertex);
573 
574  if (fillStEventWithL3GlobalTracks()!=0) return 1;
575 
576  }
577 
578  // i960 hits
579  //if ( ml3reader->getI960ClusterReader(1) )
580  // { if ( fillStEventWithi960Hits() != 0 ) return 1; } ;
581 
582 
583  // all right go home
584  return 0 ;
585 }
586 
587 
588 //_____________________________________________________________________________
589 Int_t Stl3RawReaderMaker::fillStEventWithL3GlobalTracks()
590 {
591  // get track nodes
592  StSPtrVecTrackNode& myTrackNodeVector = myStL3Trigger->trackNodes();
593  StSPtrVecTrackDetectorInfo& myTrackDetectorInfoVector = myStL3Trigger->trackDetectorInfo();
594 
595  // loop over rawdata tracks and fill them into StEvent
596  int numberOfTracks = ml3reader->getGlobalTrackReader()->getNumberOfTracks();
597  gMessMgr->Info() << "Stl3RawReaderMaker: Try to fill " << numberOfTracks << " tracks into StEvent." << endm;
598 
599  // get L3 raw tracks
600  globalTrack* globalL3Tracks = ml3reader->getGlobalTrackReader()->getTrackList();
601 
602  // get magnetic field
603  Float_t xval[3] = {0.,0.,0.};
604  Float_t bval[3];
605  gufld(xval,bval);
606  int signOfField = bval[2] < 0 ? -1 : 1;
607 
608  for (int trackindex=0; trackindex<numberOfTracks; trackindex++) {
609  StTrackDetectorInfo* detectorInfo = new StTrackDetectorInfo();
610  detectorInfo->setNumberOfPoints(globalL3Tracks[trackindex].nHits,kTpcId);
611  myTrackDetectorInfoVector.push_back(detectorInfo);
612 
613  // StTrackGeometry
614  StThreeVectorD* origin = new StThreeVectorD(globalL3Tracks[trackindex].r0 * TMath::Cos(globalL3Tracks[trackindex].phi0),
615  globalL3Tracks[trackindex].r0 * TMath::Sin(globalL3Tracks[trackindex].phi0),
616  globalL3Tracks[trackindex].z0);
617  StThreeVectorD* momentum = new StThreeVectorD(globalL3Tracks[trackindex].pt * TMath::Cos(globalL3Tracks[trackindex].psi),
618  globalL3Tracks[trackindex].pt * TMath::Sin(globalL3Tracks[trackindex].psi),
619  globalL3Tracks[trackindex].pt * globalL3Tracks[trackindex].tanl );
620  StHelixModel* helixModel = new StHelixModel( globalL3Tracks[trackindex].q,
621  globalL3Tracks[trackindex].psi,
622  0.0,
623  atan(globalL3Tracks[trackindex].tanl),
624  *origin,
625  *momentum, -1 );
626  // StGlobalTrack
628  globalTrack->setFlag(globalL3Tracks[trackindex].flag);
629  globalTrack->setLength(globalL3Tracks[trackindex].length);
630  globalTrack->setDetectorInfo(detectorInfo);
631  globalTrack->setGeometry(helixModel);
632  // set finder/fitting method
633  int l3word = kL3FitId + (1<<l3Standard);
634  //cout << l3word << endl; // should print 16390
635  globalTrack->setEncodedMethod(l3word);
636 
637  // helicity h = -sign(q*B)
638  short h = globalTrack->geometry()->charge()*signOfField > 0 ? -1 : 1;
639  globalTrack->geometry()->setHelicity(h);
640 
641  // curvature
642  double curv = fabs(0.0003 * bval[2] / globalL3Tracks[trackindex].pt);
643  globalTrack->geometry()->setCurvature(curv);
644 
645  // add dE/dx information
646  if (globalL3Tracks[trackindex].dedx > 0) {
647  StDedxPidTraits* myStDedxPidTraits = new StDedxPidTraits(kTpcId, kTruncatedMeanId,
648  globalL3Tracks[trackindex].ndedx,
649  globalL3Tracks[trackindex].dedx*1e6, 0);
650  globalTrack->addPidTraits(myStDedxPidTraits);
651  }
652 
653  // StTrackNode
654  StTrackNode* trackNode = new StTrackNode();
655  trackNode->addTrack(globalTrack);
656  myTrackNodeVector.push_back(trackNode);
657 
658 
659  //==================================
660  // check dca
661  //double kapa = fabs(0.0003 * bval[2] / globalL3Tracks[trackindex].pt);
662  //double lambda = atan(globalL3Tracks[trackindex].tanl);
663  //double phase = globalL3Tracks[trackindex].psi - h * TMath::Pi()/2;
664  //StHelixD* track = new StHelixD(kapa, lambda, phase, *origin, h);
665  //StThreeVectorD vertex(0,0,myStL3Trigger->primaryVertex()->position().z());
666  //if (globalL3Tracks[trackindex].nHits>20 && globalL3Tracks[trackindex].q<0)
667  //cout << " dca = " << track->distance(vertex) << ", bval = " << bval[2] << endl;
668  //==================================
669 
670 
671  }
672  // ok
673  return 0 ;
674 }
675 
676 
677 //_____________________________________________________________________________
678 Int_t Stl3RawReaderMaker::fillStEventWithi960Hits()
679 {
680  // create StTpcHitCollection and connect it to StEvent
681  StTpcHitCollection* mHitCollection = new StTpcHitCollection() ;
682  myStL3Trigger->setTpcHitCollection(mHitCollection);
683 
684  // prepare transformation from pad,time,padrow to x,y,z
685  St_l3_Coordinate_Transformer transformer ;
686  transformer.Print_parameters() ;
687  St_l3_xyz_Coordinate XYZ(0,0,0) ;
688  St_l3_ptrs_Coordinate PTRS(0,0,0,0) ;
689 
690  // loop over clusters and fill them into StEvent
691  Int_t totalcluster = 0 ;
692  for ( Int_t secindex=1 ; secindex<=24 ; secindex+=2 )
693  {
694  if ( ml3reader->getI960ClusterReader(secindex)->getClusterList() )
695  {
696  gMessMgr->Info() << "Stl3RawReaderMaker: Found some i960 clusters in sector:" << secindex <<endm ;
697  L3_Cluster* myl3cluster = ml3reader->getI960ClusterReader(secindex)->getClusterList() ;
698  Int_t numOfClusters = ml3reader->getI960ClusterReader(secindex)->getNumberOfClusters() ;
699 
700  for (Int_t clindex=0 ; clindex < numOfClusters ; clindex++)
701  {
702  totalcluster++ ;
703  // pad,time,row,q,flag
704  Double_t pad = ((Double_t)(myl3cluster[clindex].pad)) / 64 ;
705  Double_t time = ((Double_t)(myl3cluster[clindex].time)) / 64 ;
706  Int_t row = myl3cluster[clindex].padrow ;
707  Int_t q = myl3cluster[clindex].charge ;
708  Char_t flag = myl3cluster[clindex].flags ;
709 
710  // upper 4 bits are RB, lower 4 bits are MZ
711  Int_t RBMZ = myl3cluster[clindex].RB_MZ ;
712  Int_t rb = RBMZ >> 4 ;
713  Int_t mz = RBMZ & 15 ;
714 
715  // determine sector
716  Int_t sector = 0 ;
717  if (rb<=6) { sector = secindex; } else { sector = secindex+1; }
718 
719  // coordinate transformation
720  PTRS.Setptrs((Double_t) pad, (Double_t) time,(Double_t) row, (Double_t) sector) ;
721  transformer.raw_to_global(PTRS,XYZ) ;
722 
723  // some output
724  if (clindex%500==0)
725  {
726  cout << XYZ.Getx() <<"\t" << XYZ.Gety() <<"\t" << XYZ.Getz() <<"\t";
727  cout << row <<"\t" << sector << "\t" << rb << "\t" << mz << "\t" << q << "\t";
728  cout << ((Double_t)(myl3cluster[clindex].pad)) / 64 <<"\t" ;
729  cout << ((Double_t)(myl3cluster[clindex].time)) / 64 <<"\t" ;
730  cout << (Int_t)(myl3cluster[clindex].padrow) <<"\n" ;
731  }
732 
733  // Fill it
734  // position and error
735  StThreeVectorF* pos = new StThreeVectorF(XYZ.Getx(),XYZ.Gety(),XYZ.Getz()) ;
736  StThreeVectorF* poserror = new StThreeVectorF(0,0,0) ;
737  // pack sec and row : bits 4-8 = sector[1,24] and bits 9-14 = padrow[1-45]
738  ULong_t hw = 0 ;
739  ULong_t hrow = 0 ;
740  ULong_t hsec = 0 ;
741  if ( row >=1 && row <=45 ) { hrow = row << 9 ; } else { hrow=0 ; }
742  if ( sector >=1 && sector <=24 ) { hsec = sector << 4 ; } else { hsec=0 ; }
743  hw = hw | hrow ;
744  hw = hw | hsec ;
745  // track reference counter set always to 0
746  UChar_t c = 0 ;
747  // create hit
748  StTpcHit* tpcHit = new StTpcHit(*pos,*poserror,hw,q,c) ;
749  tpcHit->setFlag(flag) ;
750  // add to hit collection
751  if (tpcHit) { mHitCollection->addHit(tpcHit) ;} else { delete tpcHit; return 1;}
752  } // clusters
753  } // if ...
754  } // sectors
755  cout << "total found clusters " << totalcluster << endl;
756  // ok
757  return 0 ;
758 }
759 
760 
Int_t m_Mode
counters
Definition: StMaker.h:81
Definition: rb.hh:21
virtual Long_t GetNRows() const
Returns the number of the used rows for the wrapped table.
Definition: TTable.cxx:1388
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: TTable.h:48
Definition: Stypes.h:44
Definition: Stypes.h:41
virtual TDataSet * Find(const char *path) const
Definition: TDataSet.cxx:362