StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StAnalysisMaker.cxx
1 //
2 // This is a STAR typical comment header. You should modify
3 // it to reflect your changes.
4 // As a minimum it should contain the name of the author, the
5 // date it was written/modified, and a short description of what
6 // the class is meant to do. The cvs strings $X$ (where X=Id, Log)
7 // are not needed when you do not intend to put the file under
8 // cvs control. Remove them.
9 //
25 //
26 // Include header files. What has to be included strongly depends
27 // on your implementation. StEventTypes.h contains all includes
28 // you need to use StEvent.
29 //
30 //#define __TPC_LOCAL_COORDINATES__
31 #include "StAnalysisMaker.h"
32 #include "StEventTypes.h"
33 #include "StMessMgr.h"
34 #include "StDcaGeometry.h"
35 #if ROOT_VERSION_CODE < 334081
36 #include "TArrayL.h"
37 #else
38 #include "TArrayL64.h"
39 #endif
40 #include "TClassTable.h"
41 #include "TNtuple.h"
42 #include "StThreeVectorF.hh"
43 #include "StDetectorName.h"
44 #ifdef __TPC_LOCAL_COORDINATES__
45 #include "StDbUtilities/StTpcCoordinateTransform.hh"
46 #include "StDbUtilities/StCoordinates.hh"
47 #endif /* __TPC_LOCAL_COORDINATES__ */
48 //
49 // The following line defines a static string. Currently it contains
50 // the cvs Id. The compiler will put the string (literally) in the
51 // object file. It thus ends up in the shared library.
52 // The UNIX command 'strings' allowsto print all printable character in
53 // a non-text file. This way one can check the version of the file
54 // contained in a given shared library. If you do not intend to put
55 // the file under cvs control (likely) you can remove the line.
56 
57 //
58 // Proptotypes of little functions which perform
59 // specific analysis tasks. You'll find them
60 // in the same directory as StAnalysisMaker.cxx.
61 // You most likely will not need them but they can serve
62 // as an example for your own functions.
63 //
64 //
65 // This is needed to make your maker work in root4star.
66 // It can be place anywhere in the file. Note that this
67 // is a macro, that's why the ';' is missing.
68 //
69 ClassImp(StAnalysisMaker);
70 
72 StAnalysisMaker::StAnalysisMaker(const Char_t *name) : StMaker(name)
73 {
74  mEventCounter = 0;
75 }
77  //
78  // A good place for printout and to summarize
79  // the run.
80  //
81  gMessMgr->Info() << "StAnalysisMaker::Finish() "
82  << "Processed " << mEventCounter << " events." << endm;
83 
84  return kStOK;
85 }
86 
92  mEventCounter++; // increase counter
93 
94  //
95  // Get pointer to StEvent
96  //
97  StEvent* event;
98  event = (StEvent *) GetInputDS("StEvent");
99  if (!event){
100  gMessMgr->Warning() << "StAnalysisMaker::Make : No StEvent" << endm;
101  return kStOK; // if no event, we're done
102  }
103 
104  //
105  // The following is only needed since the
106  // QA folks use this maker for their QA runs.
107  // You do not need this.
108  //
109  summarizeEvent(event, mEventCounter);
110 
111  //
112  // See if this event survives the event filter.
113  // If not we stop here right away.
114  //
115  if (!accept(event)){
116  gMessMgr->Warning() << "StAnalysisMaker::Make : Event was not accepted" << endm;
117  return kStOK;
118  }
119  return kStOK;
120 }
121 
122 bool StAnalysisMaker::accept(StEvent* event)
123 {
124  //
125  // This is a kind of very simple event filter.
126  // We select only events with a valid event vertex,
127  // i.e. event->primaryVertex() returns a non-zero pointer.
128  //
129  return event->primaryVertex();
130 }
131 
132 bool StAnalysisMaker::accept(StTrack* track)
133 {
134  //
135  // This is a kind of very simple track filter.
136  // We only check for positive flags.
137  // Note that this method works for global and
138  // primary tracks since we deal with the base
139  // class only (StTrack).
140  //
141  return track && track->flag() >= 0;
142 }
143 //________________________________________________________________________________
144 void StAnalysisMaker::PrintStEvent(TString opt) {
145  // opt = vpgl3 => "v" print vertex, "p" and primary tracks, "g" print global tracks, "l3" global tracks
146  StEvent* pEvent = (StEvent*) StMaker::GetChain()->GetInputDS("StEvent");
147  if (!pEvent) return;
148  cout << "Event: Run "<< pEvent->runId() << " Event No: " << pEvent->id() << endl;
149  UInt_t NpVX = pEvent->numberOfPrimaryVertices();
150  if (NpVX) {
151  if (opt.Contains("v",TString::kIgnoreCase)) {
152  for (UInt_t i = 0; i < NpVX; i++) {
153  const StPrimaryVertex *vx = pEvent->primaryVertex(i);
154  vx->Print(Form("Vertex: %3i ",i));
155 #ifdef StTrackMassFit_hh
156  const StTrackMassFit *pf = vx->parent();
157  if (pf) cout << *pf << endl;
158 #endif
159  if (opt.Contains("p",TString::kIgnoreCase)) {
160  UInt_t nDaughters = vx->numberOfDaughters();
161  for (UInt_t j = 0; j < nDaughters; j++) {
162  const StTrack* track = vx->daughter(j);
163  if (! track) continue;
164  cout << *((const StPrimaryTrack *)track) << endl;
165 #ifdef StTrackMassFit_hh
166  const StVertex* vxEnd = track->endVertex();
167  if (vxEnd) cout << *vxEnd << endl;
168 #endif
169  }
170 #ifdef StTrackMassFit_hh
171  UInt_t nMassFits = vx->numberOfMassFits();
172  for (UInt_t j = 0; j < nMassFits; j++) {
173  const StTrackMassFit * track = vx->massFit(j);
174  if (! track) continue;
175  cout << *track << endl;
176  }
177 #endif
178  }
179  }
180  }
181  } else {
182  cout << "Event: Vertex Not Found" << endl;
183  }
184  if (opt.Contains("g",TString::kIgnoreCase)) {
185  StSPtrVecTrackNode& trackNode = pEvent->trackNodes();
186  UInt_t nTracks = trackNode.size();
187  StTrackNode *node = 0;
188 #ifndef StTrackMassFit_hh
189  cout << " Global tracks " << nTracks << endl;
190 #else
191  cout << nTracks << " Track nodes" << endl;
192 #endif
193  for (UInt_t i=0; i < nTracks; i++) {
194  node = trackNode[i]; if (!node) continue;
195 #ifdef StTrackMassFit_hh
196  cout << *node << endl;
197 #else
198  UInt_t nentries = node->entries();
199  for (UInt_t j = 0; j < nentries; j++) {
200  StTrack *track = node->track(j);
201  if (! track) continue;
202  if (track->type() == global) {
203  StGlobalTrack* gTrack = (StGlobalTrack* ) track;
204  cout << *gTrack << endl;
205  } else if (track->type() == primary) {
206  StPrimaryTrack* pTrack = (StPrimaryTrack* ) track;
207  cout << *pTrack << endl;
208  }
209  }
210 #endif
211  }
212  }
213  if (opt.Contains("l3",TString::kIgnoreCase)) {
214  if (pEvent->l3Trigger()) {
215  StSPtrVecTrackNode& trackNode = pEvent->l3Trigger()->trackNodes();
216  UInt_t nTracks = trackNode.size();
217  StTrackNode *node = 0;
218  cout << " L3 global tracks " << nTracks << endl;
219  for (UInt_t i=0; i < nTracks; i++) {
220  node = trackNode[i]; if (!node) continue;
221  StGlobalTrack* gTrack = dynamic_cast<StGlobalTrack*>(node->track(global));
222  if (gTrack) cout << *gTrack << endl;
223  }
224  }
225  }
226 }
227 //________________________________________________________________________________
228 void StAnalysisMaker::PrintGlobalTrack(Int_t itk) {
229  StEvent* pEvent = (StEvent*) StMaker::GetChain()->GetInputDS("StEvent");
230  if (!pEvent) return;
231  cout << "Event: Run "<< pEvent->runId() << " Event No: " << pEvent->id() << endl;
232  StSPtrVecTrackNode& trackNode = pEvent->trackNodes();
233  UInt_t nTracks = trackNode.size();
234  StTrackNode *node = 0;
235  cout << " Global tracks " << endl;
236  for (UInt_t i = 0; i < nTracks; i++) {
237  node = trackNode[i]; if (!node) continue;
238  StGlobalTrack* gTrack = static_cast<StGlobalTrack*>(node->track(global));
239  if (itk != 0 && gTrack->key() != itk) continue;
240  cout << *gTrack << endl;
241  if (! gTrack->detectorInfo()) {cout << "=============== detectorInfo is missing" << endl; continue;}
242  StPtrVecHit hvec = gTrack->detectorInfo()->hits();
243  for (UInt_t j=0; j<hvec.size(); j++) {// hit loop
244  if (hvec[j]->detector() == kTpcId) {
245  StTpcHit *tpcHit = static_cast<StTpcHit *> (hvec[j]);
246  if (! tpcHit) continue;
247  cout << *tpcHit << endl;
248  } else {
249  cout << *hvec[j] << endl;
250  }
251  }
252  }
253 }
254 //________________________________________________________________________________
255 void StAnalysisMaker::PrintVertex(Int_t ivx) {
256  // opt = vpg => "v" print vertex, "p" and primary tracks, "g" print global tracks
257  StEvent* pEvent = (StEvent*) StMaker::GetChain()->GetInputDS("StEvent");
258  if (!pEvent) return;
259  cout << "Event: Run "<< pEvent->runId() << " Event No: " << pEvent->id() << endl;
260  UInt_t NpVX = pEvent->numberOfPrimaryVertices();
261  if (NpVX) {
262  for (Int_t i = 0; i < NpVX; i++) {
263  if (ivx >= 0 && i != ivx) continue;
264  const StPrimaryVertex *vx = pEvent->primaryVertex(i);
265  vx->Print(Form("Vertex: %3i ",i));
266  UInt_t nDaughters = vx->numberOfDaughters();
267  for (UInt_t j = 0; j < nDaughters; j++) {
268  StPrimaryTrack* pTrack = (StPrimaryTrack*) vx->daughter(j);
269  if (! pTrack) continue;
270  cout << *pTrack << endl;
271  if (! pTrack->detectorInfo()) {cout << "=============== detectorInfo is missing" << endl; continue;}
272  StPtrVecHit hvec = pTrack->detectorInfo()->hits();
273  for (UInt_t j=0; j<hvec.size(); j++) {// hit loop
274  if (hvec[j]->detector() == kTpcId) {
275  StTpcHit *tpcHit = static_cast<StTpcHit *> (hvec[j]);
276  if (! tpcHit) continue;
277  cout << *tpcHit << endl;
278  } else {
279  cout << *hvec[j] << endl;
280  }
281  }
282  }
283  }
284  } else {
285  cout << "Event: Vertex Not Found" << endl;
286  }
287 }
288 //________________________________________________________________________________
289 void StAnalysisMaker::PrintTpcHits(Int_t sector, Int_t row, Int_t plot, Int_t IdTruth) {
290  // plot = 1 => All hits;
291  // plot = 2 => prompt hits only |z| > 190
292  struct BPoint_t {
293  Float_t sector,row,x,y,z,q,adc,pad,timebucket,IdTruth,xL,yL,zL;
294  };
295  static const Char_t *vname = "sector:row:x:y:z:q:adc:pad:timebucket:IdTruth"
296 #ifdef __TPC_LOCAL_COORDINATES__
297  ":xL:yL:zL"
298 #endif /* __TPC_LOCAL_COORDINATES__ */
299  ;
300  BPoint_t BPoint;
301  static TNtuple *Nt = 0;
302  if (plot && Nt == 0) {
303  TFile *tf = StMaker::GetTopChain()->GetTFile();
304  if (tf) {tf->cd(); Nt = new TNtuple("TpcHit","TpcHit",vname);}
305  }
306  StEvent* pEvent = (StEvent*) StMaker::GetChain()->GetInputDS("StEvent");
307  if (!pEvent) { cout << "Can't find StEvent" << endl; return;}
308  // StSPtrVecTrackNode& trackNode = pEvent->trackNodes();
309  Int_t TotalNoOfTpcHits = 0;
310  StTpcHitCollection* TpcHitCollection = pEvent->tpcHitCollection();
311  if (! TpcHitCollection) { cout << "No TPC Hit Collection" << endl; return;}
312  UInt_t numberOfSectors = TpcHitCollection->numberOfSectors();
313  for (UInt_t i = 0; i< numberOfSectors; i++) {
314  if (sector == 0 || (Int_t) i+1 == sector) {
315  StTpcSectorHitCollection* sectorCollection = TpcHitCollection->sector(i);
316  if (sectorCollection) {
317  Int_t numberOfPadrows = sectorCollection->numberOfPadrows();
318  // Int_t noHits = 0;
319  for (int j = 0; j< numberOfPadrows; j++) {
320  if (row == 0 || j+1 == row) {
321  StTpcPadrowHitCollection *rowCollection = sectorCollection->padrow(j);
322  if (rowCollection) {
323  StSPtrVecTpcHit &hits = rowCollection->hits();
324 #if ROOT_VERSION_CODE < 334081
325  Long_t NoHits = hits.size();
326  TArrayL idxT(NoHits); Long_t *idx = idxT.GetArray();
327 #else
328  Long64_t NoHits = hits.size();
329  TArrayL64 idxT(NoHits); Long64_t *idx = idxT.GetArray();
330 #endif
331  if (! NoHits) continue;
332  TotalNoOfTpcHits += NoHits;
333 #if 1
334  TArrayD dT(NoHits); Double_t *d = dT.GetArray();
335  for (Long64_t k = 0; k < NoHits; k++) {
336  const StTpcHit *tpcHit = static_cast<const StTpcHit *> (hits[k]);
337 #if 0
338  const StThreeVectorF& xyz = tpcHit->position();
339  d[k] = xyz.z();
340 #else
341  d[k] = tpcHit->id();
342 #endif
343  }
344  idx[0] = 0;
345  if (NoHits > 1) TMath::Sort(NoHits,d,idx,kFALSE);
346 #endif
347  for (Long64_t k = 0; k < NoHits; k++) {
348 #if 1
349  Int_t l = idx[k];
350 #else
351  Int_t l = k;
352 #endif
353  StTpcHit *tpcHit = static_cast<StTpcHit *> (hits[l]);
354  if (! tpcHit) continue;
355  if (IdTruth >= 0 && tpcHit->idTruth() != IdTruth) continue;
356  if (! plot) tpcHit->Print();
357  else {
358  if (Nt) {
359  const StThreeVectorF& xyz = tpcHit->position();
360  if (plot == 2 && TMath::Abs(xyz.z()) < 195.0) continue;
361 #ifdef __TPC_LOCAL_COORDINATES__
363  StGlobalCoordinate glob(xyz);
365  tran(glob,lTpc,i+1,j+1);
366  BPoint.xL = lTpc.position().x();
367  BPoint.yL = lTpc.position().y();
368  BPoint.zL = lTpc.position().z();
369 #endif /* __TPC_LOCAL_COORDINATES__ */
370  BPoint.sector = i+1;
371  BPoint.row = j+1;
372  BPoint.x = xyz.x();
373  BPoint.y = xyz.y();
374  BPoint.z = xyz.z();
375  BPoint.q = 1.e6*tpcHit->charge();
376  BPoint.adc = tpcHit->adc();
377  BPoint.pad = tpcHit->pad();
378  BPoint.timebucket = tpcHit->timeBucket();
379  BPoint.IdTruth = tpcHit->idTruth();
380  Nt->Fill(&BPoint.sector);
381  }
382  }
383  }
384  }
385  }
386  }
387  }
388  }
389  // break;
390  }
391  cout << "TotalNoOfTpcHits = " << TotalNoOfTpcHits << endl;
392 }
393 //________________________________________________________________________________
394 void StAnalysisMaker::PrintEmcHits(Int_t det, Int_t mod, const Option_t *opt) {
395  TString Opt(opt);
396  StEvent* pEvent = (StEvent*) StMaker::GetChain()->GetInputDS("StEvent");
397  if (!pEvent) { cout << "Can't find StEvent" << endl; return;}
398  StEmcCollection* emccol=(StEmcCollection*) pEvent->emcCollection();
399  if (! emccol) { cout << "No Emc Hit Collection" << endl; return;}
400  //cout <<"Filling hits and clusters \n";
401  Int_t d1 = 0, d2 = 7;
402  if (det >= 0 && det <= 7) {d1 = d2 = det;}
403  for(Int_t d = d1; d <= d2; d++) {
404  StDetectorId id = static_cast<StDetectorId>(d+kBarrelEmcTowerId);
405  if (id != kBarrelEmcTowerId && id != kEndcapEmcTowerId &&
406  ! Opt.Contains("Pre",TString::kIgnoreCase)) continue;
407  const StEmcDetector* detector=emccol->detector(id);
408  if(detector) {
409  Int_t maxMod = 121;
410  if (d > 3) maxMod = 14;
411  Int_t j1 = 1;
412  if (mod > 0 and mod < maxMod) {j1 = maxMod = mod;}
413  //cout <<"Filling hits for detetor "<<EmcDet<<endl;
414  if (Opt.Contains("Adc",TString::kIgnoreCase)) {
415  for(Int_t j = j1; j < maxMod; j++) {
416  const StEmcModule* module = detector->module(j);
417  if(module) {
418  const StSPtrVecEmcRawHit& rawHit=module->hits();
419  Int_t nhits = (Int_t) rawHit.size();
420  for(Int_t k = 0; k < nhits; k++) if (rawHit[k]->energy() > 0) cout << DetectorName(id) << "\t" << *rawHit[k] << endl;
421  }
422  }
423  }
424  if (Opt.Contains("Clu",TString::kIgnoreCase)) {
425  const StEmcClusterCollection *cl = detector->cluster();
426  if (cl) {
427  Int_t NoCls = cl->numberOfClusters();
428  if (NoCls) {
429  const StSPtrVecEmcCluster& clusters = cl->clusters();
430  for (Int_t i = 0; i < NoCls; i++) {
431  if (clusters[i]->energy() > 0) cout << DetectorName(id) << "\t" << *clusters[i] << endl;
432  }
433  }
434  }
435  }
436  }
437  }
438  if (Opt.Contains("Point",TString::kIgnoreCase)) {
439  const StSPtrVecEmcPoint& bp = emccol->barrelPoints();
440  const StSPtrVecEmcPoint& ep = emccol->endcapPoints();
441  for (Int_t i = 0; i < 2; i++) {// barrel & endcap
442  const StSPtrVecEmcPoint& p = (i == 0) ? bp : ep;
443  Int_t np = (Int_t) p.size();
444  if (np) {
445  cout << "Found " << np << " Points in ";
446  if (! i) cout << "Barrel";
447  else cout << "Encap";
448  cout << endl;
449  for (Int_t j = 0; j < np; j++) {
450  cout << *p[j] << endl;
451  }
452  }
453  }
454  }
455 }
456 //________________________________________________________________________________
457 void StAnalysisMaker::PrintSvtHits() {
458  UInt_t i,j,k,l;
459  // Double_t zPrim = 0;
460  StEvent* pEvent = (StEvent*) StMaker::GetChain()->GetInputDS("StEvent");
461  if (!pEvent) return;
462  // if (pEvent->numberOfPrimaryVertices() != 1) return;
463  StPrimaryVertex *primaryVertex = pEvent->primaryVertex();
464  if ( primaryVertex) {
465  const StThreeVectorF &primXYZ = primaryVertex->position();
466  // cout << "primaryVertex " << primXYZ << endl;
467  cout << "primaryVertex \t" << primXYZ.x() << "\t" << primXYZ.y() << "\t" << primXYZ.z() << endl;
468  }
469  Int_t TotalNoOfSvtHits = 0;
470  StSvtHitCollection* SvtHitCollection = pEvent->svtHitCollection();
471  if (! SvtHitCollection) { cout << "No SVT Hit Collection" << endl; return;}
472  UInt_t numberOfBarrels = SvtHitCollection->numberOfBarrels();
473  // Int_t vers = gClassTable->GetID("StSvtHit");
474  for ( i = 0; i< numberOfBarrels; i++) {
475  StSvtBarrelHitCollection* barrelCollection = SvtHitCollection->barrel(i);
476  if (barrelCollection) {
477  UInt_t numberOfLadders = barrelCollection->numberOfLadders();
478  // UInt_t noHits = 0;
479  for (j = 0; j< numberOfLadders; j++) {
480  StSvtLadderHitCollection *ladderCollection = barrelCollection->ladder(j);
481  if (ladderCollection) {
482  UInt_t numberOfWafers = ladderCollection->numberOfWafers();
483  for (k = 0; k < numberOfWafers; k++) {
484  StSvtWaferHitCollection* waferCollection = ladderCollection->wafer(k);
485  StSPtrVecSvtHit &hits = waferCollection->hits();
486  UInt_t NoHits = hits.size();
487  for (l = 0; l < NoHits; l++) {
488  StSvtHit *hit = hits[l];
489  if (hit) {
490  // cout << *((StHit *) hit) << endl;
491  TotalNoOfSvtHits++;
492  hit->Print();
493  }
494  }
495  }
496  }
497  }
498  }
499  }
500  cout << "Total no. of Svt Hits " << TotalNoOfSvtHits << endl;
501 }
502 //________________________________________________________________________________
503 void StAnalysisMaker::PrintSsdHits() {
504  UInt_t i,k,l;
505  // Double_t zPrim = 0;
506  StEvent* pEvent = (StEvent*) StMaker::GetChain()->GetInputDS("StEvent");
507  if (!pEvent) return;
508  // if (pEvent->numberOfPrimaryVertices() != 1) return;
509  StPrimaryVertex *primaryVertex = pEvent->primaryVertex();
510  if ( primaryVertex) {
511  const StThreeVectorF &primXYZ = primaryVertex->position();
512  // cout << "primaryVertex " << primXYZ << endl;
513  cout << "primaryVertex \t" << primXYZ.x() << "\t" << primXYZ.y() << "\t" << primXYZ.z() << endl;
514  }
515  // Int_t TotalNoOfSsdHits = 0;
516  StSsdHitCollection* SsdHitCollection = pEvent->ssdHitCollection();
517  if (! SsdHitCollection) { cout << "No SSD Hit Collection" << endl; return;}
518  UInt_t numberOfLadders = SsdHitCollection->numberOfLadders();
519  // Int_t vers = gClassTable->GetID("StSsdHit");
520  for ( i = 0; i< numberOfLadders; i++) {
521  StSsdLadderHitCollection* ladderCollection = SsdHitCollection->ladder(i);
522  if (ladderCollection) {
523  UInt_t numberOfWafers = ladderCollection->numberOfWafers();
524  for (k = 0; k < numberOfWafers; k++) {
525  StSsdWaferHitCollection* waferCollection = ladderCollection->wafer(k);
526  StSPtrVecSsdHit &hits = waferCollection->hits();
527  UInt_t NoHits = hits.size();
528  for (l = 0; l < NoHits; l++) {
529  StSsdHit *hit = hits[l];
530  if (hit) {
531  hit->Print("");
532  }
533  }
534  }
535  }
536  }
537 }
538 //________________________________________________________________________________
539 void StAnalysisMaker::PrintSstHits() {
540  UInt_t i,k,l;
541  // Double_t zPrim = 0;
542  StEvent* pEvent = (StEvent*) StMaker::GetChain()->GetInputDS("StEvent");
543  if (!pEvent) return;
544  // if (pEvent->numberOfPrimaryVertices() != 1) return;
545  StPrimaryVertex *primaryVertex = pEvent->primaryVertex();
546  if ( primaryVertex) {
547  const StThreeVectorF &primXYZ = primaryVertex->position();
548  // cout << "primaryVertex " << primXYZ << endl;
549  cout << "primaryVertex \t" << primXYZ.x() << "\t" << primXYZ.y() << "\t" << primXYZ.z() << endl;
550  }
551  // Int_t TotalNoOfSstHits = 0;
552  StSstHitCollection* SstHitCollection = pEvent->sstHitCollection();
553  if (! SstHitCollection) { cout << "No SST Hit Collection" << endl; return;}
554  UInt_t numberOfLadders = SstHitCollection->numberOfLadders();
555  // Int_t vers = gClassTable->GetID("StSstHit");
556  for ( i = 0; i< numberOfLadders; i++) {
557  StSstLadderHitCollection* ladderCollection = SstHitCollection->ladder(i);
558  if (ladderCollection) {
559  UInt_t numberOfWafers = ladderCollection->numberOfWafers();
560  for (k = 0; k < numberOfWafers; k++) {
561  StSstWaferHitCollection* waferCollection = ladderCollection->wafer(k);
562  StSPtrVecSstHit &hits = waferCollection->hits();
563  UInt_t NoHits = hits.size();
564  for (l = 0; l < NoHits; l++) {
565  StSstHit *hit = hits[l];
566  if (hit) {
567  hit->Print("");
568  }
569  }
570  }
571  }
572  }
573 }
574 //________________________________________________________________________________
575 void StAnalysisMaker::PrintToFHits() {
576  // Double_t zPrim = 0;
577  StEvent* pEvent = (StEvent*) StMaker::GetChain()->GetInputDS("StEvent");
578  if (!pEvent) return;
579  const StBTofCollection* tof = pEvent->btofCollection();
580  if (! tof) {LOG_QA << "No BToF collection" << endm; return;}
581  else {LOG_QA << "BToF collection";}
582  if (tof->tofHeader() && tof->tofHeader()->vpdVz() > -250) {
583  LOG_QA << " VpdZ:" << Form("%7.2f",tof->tofHeader()->vpdVz());
584  }
585  LOG_QA << endm;
586  const StSPtrVecBTofHit& tofHits = tof->tofHits();
587  for(size_t i=0;i<tofHits.size();i++) { //loop on hits in modules
588  StBTofHit *aHit = tofHits[i];
589  if(!aHit) continue;
590  LOG_QA << *aHit << endm;
591  }
592 }
593 //________________________________________________________________________________
594 void StAnalysisMaker::PrintRnDHits() {
595  UInt_t i=0,k=0,l;
596  // Double_t zPrim = 0;
597  StEvent* pEvent = (StEvent*) StMaker::GetChain()->GetInputDS("StEvent");
598  if (!pEvent) return;
599  // if (pEvent->numberOfPrimaryVertices() != 1) return;
600  StPrimaryVertex *primaryVertex = pEvent->primaryVertex();
601  if ( primaryVertex) {
602  const StThreeVectorF &primXYZ = primaryVertex->position();
603  // cout << "primaryVertex " << primXYZ << endl;
604  cout << "primaryVertex \t" << primXYZ.x() << "\t" << primXYZ.y() << "\t" << primXYZ.z() << endl;
605  }
606  // Int_t TotalNoOfRnDHits = 0;
607  StRnDHitCollection* RnDHitCollection = pEvent->rndHitCollection();
608  if (! RnDHitCollection) { cout << "No RND Hit Collection" << endl; return;}
609  StSPtrVecRnDHit &hits = RnDHitCollection->hits();
610  UInt_t NoHits = hits.size();
611  for (l = 0; l < NoHits; l++) {
612  StRnDHit *hit = hits[l];
613  if (hit) {
614  // cout << *((StHit *) hit) << endl;
615  const StThreeVectorF &P = hit->position();
616  printf("l:%2i w:%2i",i+1,k+1);
617  printf(" x: %8.3f y: %8.3f z: %8.3f ", P.x(), P.y(), P.z());
618  printf("l:%2i w:%2i",
619  hit->ladder(), hit->wafer());
620  printf(" Id: %4i Q: %4i",hit->idTruth(), hit->qaTruth());
621  printf(" Flag: %4i Fit: %3i",hit->flag(), hit->usedInFit());
622  printf("\n");
623  }
624  }
625 }
626 //________________________________________________________________________________
627 void StAnalysisMaker::summarizeEvent(StEvent *event, Int_t mEventCounter) {
628  if (! event) event = (StEvent*) StMaker::GetChain()->GetInputDS("StEvent");
629  static const UInt_t NoFitPointCutForGoodTrack = StVertex::NoFitPointCutForGoodTrack();
630  LOG_QA << "+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-" << endm
631  << "StAnalysisMaker, Reading Event: " << mEventCounter
632  << " Type: " << event->type()
633  << " Run: " << event->runId()
634  << " EventId: " << event->id();
635  const StTriggerIdCollection* triggerCol = event->triggerIdCollection();
636  if (triggerCol) {
637  const StTriggerId* nominal = triggerCol->nominal();
638  if (nominal) {
639  UInt_t maxTriggers = nominal->maxTriggerIds();
640  LOG_QA << " TriggerIds: ";
641  for (UInt_t i = 0; i < maxTriggers; i++) {
642  if (nominal->triggerId(i)) {LOG_QA << nominal->triggerId(i) << "|";}
643  }
644  }
645  }
646  LOG_QA << endm;
647  StSPtrVecTrackNode& trackNode = event->trackNodes();
648  UInt_t nTracks = trackNode.size();
649  StTrackNode *node = 0;
650  UInt_t nGoodTracks = 0;
651  UInt_t nGoodFtpcTracks = 0;
652  UInt_t nBeamBackTracks = 0;
653  UInt_t nGoodBeamBackTracks = 0;
654  UInt_t nShortTrackForEEmc = 0;
655  UInt_t nShortTrackForETOF = 0;
656  UInt_t pcTracks = 0; // PostCrossingTrack
657  UInt_t promptTracks = 0; // tracks with prompt hits
658  UInt_t crossMembrane = 0;
659  UInt_t nToFMatched = 0;
660  UInt_t nBEmcMatched = 0;
661  UInt_t nEEmcMatched = 0;
662  UInt_t nWestTpcOnly = 0;
663  UInt_t nEastTpcOnly = 0;
664  StGlobalTrack* gTrack = 0;
665  for (UInt_t i=0; i < nTracks; i++) {
666  node = trackNode[i]; if (!node) continue;
667  gTrack = static_cast<StGlobalTrack*>(node->track(global));
668  if (! gTrack) continue;
669  if (gTrack->flag() < 0) continue;
670  if (TMath::Abs(gTrack->flag())%100 == 11) nShortTrackForEEmc++;
671  if (TMath::Abs(gTrack->flag())%100 == 12) nShortTrackForETOF++;
672  if (gTrack->flag()/100 == 9) {
673  nBeamBackTracks++;
674  if (! gTrack->bad()) nGoodBeamBackTracks++;
675  }
676  if (gTrack->flag() >= 700 && gTrack->flag() < 900) nGoodFtpcTracks++;
677  if (gTrack->isPostXTrack()) pcTracks++;
678  if (gTrack->isPromptTrack()) promptTracks++;
679  if (gTrack->isMembraneCrossingTrack()) crossMembrane++;
680  if (gTrack->isToFMatched()) nToFMatched++;
681  if (gTrack->isBemcMatched()) nBEmcMatched++;
682  if (gTrack->isEemcMatched() ) nEEmcMatched++;
683  if (gTrack->fitTraits().numberOfFitPoints() < NoFitPointCutForGoodTrack) continue;
684  if (gTrack->isWestTpcOnly()) nWestTpcOnly++;
685  if (gTrack->isEastTpcOnly()) nEastTpcOnly++;
686  nGoodTracks++;
687  }
688  LOG_QA << "# track nodes: \t"
689  << nTracks << ": good globals with NFitP>="<< NoFitPointCutForGoodTrack << ": " << nGoodTracks;
690  if (nGoodFtpcTracks) {LOG_QA << ": Ftpc tracks : " << nGoodFtpcTracks;}
691  LOG_QA << endm;
692  if (nBeamBackTracks) {LOG_QA << "BeamBack tracks: " << nBeamBackTracks << ": good ones: " << nGoodBeamBackTracks;}
693  if (nShortTrackForEEmc) {LOG_QA << ": Short tracks pointing to EEMC : " << nShortTrackForEEmc;}
694  if (nShortTrackForETOF) {LOG_QA << ": Short tracks pointing to ETOF : " << nShortTrackForETOF;}
695  if (nBeamBackTracks || nShortTrackForEEmc || nShortTrackForETOF) {LOG_QA << endm;}
696  LOG_QA << "post (C)rossing tracks :" << pcTracks << ": (P)rompt:" << promptTracks
697  << ": (X) membrane :" << crossMembrane
698  << "(T)of/ctb matches:" << nToFMatched << " :Emc matches(B/E): " << nBEmcMatched << "/" << nEEmcMatched
699  << " :Only W:" << nWestTpcOnly << " E:" << nEastTpcOnly;
700  if (event->btofCollection()) {
701  if (event->btofCollection()->tofHeader() && event->btofCollection()->tofHeader()->vpdVz() > -250){
702  LOG_QA << " VpdZ:" << Form("%7.2f",event->btofCollection()->tofHeader()->vpdVz());
703  }
704  }
705  if (event->triggerData()) {
706  LOG_QA << ": ZdcZ:" << Form("%7.2f",event->triggerData()->zdcVertexZ());
707  }
708  LOG_QA << endm;
709  // Report for jobTracking Db
710  if (nTracks) {
711  // LOG_QA << "SequenceValue=" << mEventCounter
712  LOG_QA
713  << "StageID='3'"
714  << ",MessageKey=" << "'nodes all'"
715  << ",MessageValue='" << nTracks
716  << "'" << endm;
717  }
718 
719  if (nGoodTracks) {
720  // LOG_QA << "SequenceValue=" << mEventCounter
721  LOG_QA
722  << "StageID='3'"
723  << ",MessageKey=" << "'nodes good'"
724  << ",MessageValue='" << nGoodTracks
725  << "'" << endm;
726  }
727 
728  StPrimaryVertex *pVertex=0;
729  Int_t NoVertexPos = 0; // no. of vertices with positive rank
730  for (Int_t ipr=0;(pVertex=event->primaryVertex(ipr));ipr++) {
731  if (pVertex->ranking() > 0) NoVertexPos++;
732  }
733  LOG_QA << "StageID='3'" << ",MessageKey=" << "'No. of Vertices with positive rank'" << ",MessageValue='" << NoVertexPos << "'" << endm;
734  for (Int_t ipr=0;(pVertex=event->primaryVertex(ipr));ipr++) {
735 #ifdef StTrackMassFit_hh
736  Int_t key = pVertex->key();
737  if (key <= 0) pVertex->setKey(ipr);
738  LOG_QA << *pVertex << endm;
739 #else
740  LOG_QA << Form("#V[%3i]",ipr) << *pVertex << endm;
741 #endif
742  // Report for jobTracking Db (non-zero entry only)
743  if (pVertex->numberOfDaughters()) {
744  // LOG_QA << "SequenceValue=" << mEventCounter
745  LOG_QA
746  << "StageID='3'"
747  << ",MessageKey=" << "'primary all'"
748  << ",MessageValue='" << pVertex->numberOfDaughters()
749  << "'" << endm;
750  }
751  if (pVertex->numberOfGoodTracks()) {
752  // LOG_QA << "SequenceValue=" << mEventCounter
753  LOG_QA
754  << "StageID='3'"
755  << ",MessageKey=" << "'primary good'"
756  << ",MessageValue='" << pVertex->numberOfGoodTracks()
757  << "'" << endm;
758  }
759  }// end prim vtx
760  if (event->v0Vertices() .size()) {
761  LOG_QA << "# V0 vertices: " << event->v0Vertices().size() << endm;
762  StSPtrVecV0Vertex& v0Vertices = event->v0Vertices();
763  Int_t nv0 = v0Vertices.size();
764  for (Int_t iv0=0;iv0<nv0;iv0++) {
765  StV0Vertex *v0Vertex = v0Vertices[iv0];
766  if (! v0Vertex) continue;
767 #ifdef StTrackMassFit_hh
768  Int_t key = v0Vertex->key();
769  if (key <= 0) key = iv0;
770  LOG_QA << Form("#V[%3i]",key) << *v0Vertex << endm;
771 #else
772  LOG_QA << *v0Vertex << endm;
773 #endif
774  }
775  }// end prim vtx
776 
777  if (event->xiVertices() .size()) {
778  LOG_QA << "# Xi vertices: " << event->xiVertices().size() << endm;
779  }
780  if (event->kinkVertices().size()) { LOG_QA << "# Kink vertices: "
781  << event->kinkVertices().size() << endm;
782  }
783  // Report for jobTracking Db (non-zero entry only)
784  if (event->v0Vertices() .size()) {
785  // LOG_QA << "SequenceValue=" << mEventCounter
786  LOG_QA
787  << "StageID='3'"
788  << ",MessageKey=" << "'V0Vertices', " << "MessageValue=" << event->v0Vertices() .size() << endm;
789  }
790  if (event->xiVertices() .size()) {
791  // LOG_QA << "SequenceValue=" << mEventCounter
792  LOG_QA
793  << "StageID='3'"
794  << ",MessageKey=" << "'XiVertices', " << "MessageValue="<< event->xiVertices() .size() << endm;
795  }
796 
797  if (event->kinkVertices().size()) {
798  // LOG_QA << "SequenceValue=" << mEventCounter
799  LOG_QA
800  << "StageID='3'"
801  << ",MessageKey=" << "'KinkVertices'," << "MessageValue="<< event->kinkVertices().size() << endm;
802  }
803 
804  UInt_t TotalNoOfTpcHits = 0, noBadTpcHits = 0, noTpcHitsUsedInFit = 0;
805  StTpcHitCollection* TpcHitCollection = event->tpcHitCollection();
806  if (TpcHitCollection) {
807  UInt_t numberOfSectors = TpcHitCollection->numberOfSectors();
808  for (UInt_t i = 0; i< numberOfSectors; i++) {
809  StTpcSectorHitCollection* sectorCollection = TpcHitCollection->sector(i);
810  if (sectorCollection) {
811  Int_t numberOfPadrows = sectorCollection->numberOfPadrows();
812  for (Int_t j = 0; j< numberOfPadrows; j++) {
813  StTpcPadrowHitCollection *rowCollection = sectorCollection->padrow(j);
814  if (rowCollection) {
815  StSPtrVecTpcHit &hits = rowCollection->hits();
816  UInt_t NoHits = hits.size();
817  for (UInt_t k = 0; k < NoHits; k++) {
818  StTpcHit *tpcHit = static_cast<StTpcHit *> (hits[k]);
819  if (tpcHit) {
820  TotalNoOfTpcHits++;
821  if ( tpcHit->flag()) noBadTpcHits++;
822  if (tpcHit->usedInFit()) noTpcHitsUsedInFit++;
823  }
824  }
825  }
826  }
827  }
828  }
829  }
830  if (TotalNoOfTpcHits) {
831  LOG_QA << "# TPC hits: " << TotalNoOfTpcHits
832  << ":\tBad ones (! flag): " << noBadTpcHits
833  << ":\tUsed in Fit: " << noTpcHitsUsedInFit << endm;
834  }
835  UInt_t TotalNoOfSvtHits = 0, noBadSvtHits = 0, noSvtHitsUsedInFit = 0;
836  StSvtHitCollection* svthits = event->svtHitCollection();
837  if (svthits) {
838  StSvtHit* hit;
839  for (UInt_t barrel=0; barrel<svthits->numberOfBarrels(); ++barrel) {
840  StSvtBarrelHitCollection* barrelhits = svthits->barrel(barrel);
841  if (!barrelhits) continue;
842  for (UInt_t ladder=0; ladder<barrelhits->numberOfLadders(); ++ladder) {
843  StSvtLadderHitCollection* ladderhits = barrelhits->ladder(ladder);
844  if (!ladderhits) continue;
845  for (UInt_t wafer=0; wafer<ladderhits->numberOfWafers(); ++wafer) {
846  StSvtWaferHitCollection* waferhits = ladderhits->wafer(wafer);
847  if (!waferhits) continue;
848  const StSPtrVecSvtHit& hits = waferhits->hits();
849  for (const_StSvtHitIterator it=hits.begin(); it!=hits.end(); ++it) {
850  hit = static_cast<StSvtHit*>(*it);
851  if (!hit) continue;
852  TotalNoOfSvtHits++;
853  if (hit->flag() >3) noBadSvtHits++;
854  if (hit->usedInFit()) noSvtHitsUsedInFit++;
855  }
856  }
857  }
858  }
859  }
860  if (TotalNoOfSvtHits) {
861  LOG_QA << "# SVT hits: " << TotalNoOfSvtHits
862  << ":\tBad ones(flag >3): " << noBadSvtHits
863  << ":\tUsed in Fit: " << noSvtHitsUsedInFit << endm;
864  }
865  UInt_t TotalNoOfPxlHits = 0, noBadPxlHits = 0, noPxlHitsUsedInFit = 0;
866  StPxlHitCollection* pxlhits = event->pxlHitCollection();
867  if (pxlhits) {
868  StPxlHit* hit;
869  for (UInt_t sector=0; sector<pxlhits->numberOfSectors(); ++sector) {
870  StPxlSectorHitCollection* sectorhits = pxlhits->sector(sector);
871  if (!sectorhits) continue;
872  for (UInt_t ladder=0; ladder<sectorhits->numberOfLadders(); ++ladder) {
873  StPxlLadderHitCollection* ladderhits = sectorhits->ladder(ladder);
874  if (!ladderhits) continue;
875  for (UInt_t sensor=0; sensor<ladderhits->numberOfSensors(); ++sensor) {
876  StPxlSensorHitCollection* sensorhits = ladderhits->sensor(sensor);
877  if (!sensorhits) continue;
878  const StSPtrVecPxlHit& hits = sensorhits->hits();
879  for (const_StPxlHitIterator it=hits.begin(); it!=hits.end(); ++it) {
880  hit = static_cast<StPxlHit*>(*it);
881  if (!hit) continue;
882  TotalNoOfPxlHits++;
883  if (hit->flag() >3) noBadPxlHits++;
884  if (hit->usedInFit()) noPxlHitsUsedInFit++;
885  }
886  }
887  }
888  }
889  }
890  if (TotalNoOfPxlHits) {
891  LOG_QA << "# PXL hits: " << TotalNoOfPxlHits
892  << ":\tBad ones(flag >3): " << noBadPxlHits
893  << ":\tUsed in Fit: " << noPxlHitsUsedInFit << endm;
894  }
895  UInt_t TotalNoOfIstHits = 0, noBadIstHits = 0, noIstHitsUsedInFit = 0;
896  StIstHitCollection* isthits = event->istHitCollection();
897  if (isthits) {
898  StIstHit* hit;
899  for (Int_t ladder=0; ladder<kIstNumLadders; ++ladder) {
900  StIstLadderHitCollection* ladderhits = isthits->ladder(ladder);
901  if (!ladderhits) continue;
902  for (Int_t sensor=0; sensor<kIstNumSensorsPerLadder; ++sensor) {
903  StIstSensorHitCollection* sensorhits = ladderhits->sensor(sensor);
904  if (!sensorhits) continue;
905  const StSPtrVecIstHit& hits = sensorhits->hits();
906  for (const_StIstHitIterator it=hits.begin(); it!=hits.end(); ++it) {
907  hit = static_cast<StIstHit*>(*it);
908  if (!hit) continue;
909  TotalNoOfIstHits++;
910  if (hit->flag() >3) noBadIstHits++;
911  if (hit->usedInFit()) noIstHitsUsedInFit++;
912  }
913  }
914  }
915  }
916  if (TotalNoOfIstHits) {
917  LOG_QA << "# IST hits: " << TotalNoOfIstHits
918  << ":\tBad ones(flag >3): " << noBadIstHits
919  << ":\tUsed in Fit: " << noIstHitsUsedInFit << endm;
920  }
921 
922  UInt_t TotalNoOfSsdHits = 0, noBadSsdHits = 0, noSsdHitsUsedInFit = 0;
923  StSsdHitCollection* ssdhits = event->ssdHitCollection();
924  if (ssdhits) {
925  StSsdHit* hit;
926  for (UInt_t ladder=0; ladder<ssdhits->numberOfLadders(); ++ladder) {
927  StSsdLadderHitCollection* ladderhits = ssdhits->ladder(ladder);
928  if (!ladderhits) continue;
929  for (UInt_t wafer=0; wafer<ladderhits->numberOfWafers(); ++wafer) {
930  StSsdWaferHitCollection* waferhits = ladderhits->wafer(wafer);
931  if (!waferhits) continue;
932  const StSPtrVecSsdHit& hits = waferhits->hits();
933  for (const_StSsdHitIterator it=hits.begin(); it!=hits.end(); ++it) {
934  hit = static_cast<StSsdHit*>(*it);
935  if (!hit) continue;
936  TotalNoOfSsdHits++;
937  if (hit->flag() >3) noBadSsdHits++;
938  if (hit->usedInFit()) noSsdHitsUsedInFit++;
939  }
940  }
941  }
942  }
943  if (TotalNoOfSsdHits) {
944  LOG_QA << "# SSD hits: " << TotalNoOfSsdHits
945  << ":\tBad ones(flag>3): " << noBadSsdHits
946  << ":\tUsed in Fit: " << noSsdHitsUsedInFit << endm;
947  }
948  UInt_t TotalNoOfSstHits = 0, noBadSstHits = 0, noSstHitsUsedInFit = 0;
949  StSstHitCollection* ssthits = event->sstHitCollection();
950  if (ssthits) {
951  StSstHit* hit;
952  for (UInt_t ladder=0; ladder<ssthits->numberOfLadders(); ++ladder) {
953  StSstLadderHitCollection* ladderhits = ssthits->ladder(ladder);
954  if (!ladderhits) continue;
955  for (UInt_t wafer=0; wafer<ladderhits->numberOfWafers(); ++wafer) {
956  StSstWaferHitCollection* waferhits = ladderhits->wafer(wafer);
957  if (!waferhits) continue;
958  const StSPtrVecSstHit& hits = waferhits->hits();
959  for (const_StSstHitIterator it=hits.begin(); it!=hits.end(); ++it) {
960  hit = static_cast<StSstHit*>(*it);
961  if (!hit) continue;
962  TotalNoOfSstHits++;
963  if (hit->flag() >3) noBadSstHits++;
964  if (hit->usedInFit()) noSstHitsUsedInFit++;
965  }
966  }
967  }
968  }
969  if (TotalNoOfSstHits) {
970  LOG_QA << "# SST hits: " << TotalNoOfSstHits
971  << ":\tBad ones(flag>3): " << noBadSstHits
972  << ":\tUsed in Fit: " << noSstHitsUsedInFit << endm;
973  }
974  UInt_t TotalNoOfFtpcHits = 0, noBadFtpcHits = 0, noFtpcHitsUsedInFit = 0;
975  StFtpcHitCollection* ftpchits = event->ftpcHitCollection();
976  if (ftpchits) {
977  StFtpcHit* hit;
978  for (UInt_t plane=0; plane<ftpchits->numberOfPlanes(); ++plane) {
979  StFtpcPlaneHitCollection* planehits = ftpchits->plane(plane);
980  if (!planehits) continue;
981  for (UInt_t sector=0; sector<planehits->numberOfSectors(); ++sector) {
982  StFtpcSectorHitCollection* sectorhits = planehits->sector(sector);
983  if (!sectorhits) continue;
984  const StSPtrVecFtpcHit& hits = sectorhits->hits();
985  for (const_StFtpcHitIterator it=hits.begin(); it!=hits.end(); ++it) {
986  hit = static_cast<StFtpcHit*>(*it);
987  if (!hit) continue;
988  TotalNoOfFtpcHits++;
989  /*
990  bit0:unfolded
991  bit1:unfold failed
992  bit2:saturated
993  bit3:bad shape
994  bit4:cut off
995  bit5:tracked
996  bit6:global coords
997  bit7:don't use for tracking
998 
999  I assume good hits have bit 0 and 5 (if included on a track) on
1000 
1001  Joern and Marcus - is this correct?
1002 
1003  Janet
1004  */
1005  if (! ( hit->flag() & 1 || hit->flag() & (1 << 5))) noBadFtpcHits++;
1006  else if (hit->flag() & (1 << 5)) noFtpcHitsUsedInFit++;
1007  }
1008  }
1009  }
1010  }
1011  if (TotalNoOfFtpcHits) {
1012  LOG_QA << "# FTPC hits: " << TotalNoOfFtpcHits
1013  << ":\tBad ones(!bit0): " << noBadFtpcHits
1014  << ":\tUsed in Fit: " << noFtpcHitsUsedInFit << endm;
1015  }
1016 #ifdef StRnDHit_hh
1017  StRnDHitCollection* rndhits = event->rndHitCollection();
1018  if (rndhits) {
1019  StSPtrVecRnDHit& hits = rndhits->hits();
1020  Int_t NoHits = rndhits->numberOfHits();
1021  if (NoHits) {
1022  struct NoHits_t {
1023  StDetectorId kId;
1024  const Char_t *Name;
1025  Int_t TotalNoOfHits;
1026  Int_t noBadHits;
1027  Int_t noHitsUsedInFit;
1028  };
1029  const Int_t NHtypes = 4;
1030  NoHits_t Hits[7] = {
1031  {kPxlId, "Hft", 0, 0, 0},
1032  {kIstId, "Ist", 0, 0, 0},
1033  {kFgtId, "Fgt", 0, 0, 0},
1034  {kUnknownId,"UnKnown", 0, 0, 0}
1035  };
1036  StRnDHit* hit;
1037  for (Int_t i = 0; i < NoHits; i++) {
1038  hit = hits[i];
1039  Int_t j = 0;
1040  for (j = 0; j < NHtypes-1; j++) if ( Hits[j].kId == hit->detector()) break;
1041  Hits[j].TotalNoOfHits++;
1042  if (hit->flag()) Hits[j].noBadHits++;
1043  if (hit->usedInFit()) Hits[j].noHitsUsedInFit++;
1044  }
1045  for (Int_t j = 0; j < NHtypes; j++) {
1046  if (Hits[j].TotalNoOfHits) {
1047  LOG_QA << "# " << Hits[j].Name << " hits: " << Hits[j].TotalNoOfHits
1048  << ":\tBad ones: " << Hits[j].noBadHits
1049  << ":\tUsed in Fit: " << Hits[j].noHitsUsedInFit << endm;
1050  }
1051  }
1052  }
1053  }
1054 #endif /* StRnDHit_hh */
1055  StEmcCollection* emccol = event->emcCollection();
1056  if (emccol) {
1057  const Char_t *Names[2] = {"EMC ","EEMC"};
1058  for (Int_t be = 0; be < 2; be++) {// Barrel and Endcap
1059  Int_t d1 = 0;
1060  Int_t d2 = 3;
1061  if (be) {d1 = 4; d2 =7;}
1062  Int_t Adcs[4] = {0, 0, 0, 0};
1063  Int_t Cls[4] = {0, 0, 0, 0};
1064  for(Int_t d = d1; d <= d2; d++) {
1065  StDetectorId id = static_cast<StDetectorId>(d+kBarrelEmcTowerId);
1066  const StEmcDetector* detector=emccol->detector(id);
1067  if (detector) {
1068  Int_t maxMod = 121;
1069  if (d > 3) maxMod = 14;
1070  for(Int_t j = 1; j < maxMod; j++) {
1071  const StEmcModule* module = detector->module(j);
1072  if(module) {
1073  const StSPtrVecEmcRawHit& rawHit=module->hits();
1074  for(UInt_t k=0;k<rawHit.size();k++) { //loop on hits in modules
1075  if (rawHit[k]->energy() <= 0.1) continue;
1076  Adcs[d-d1]++;
1077  }
1078  }
1079  }
1080  const StEmcClusterCollection *cl = detector->cluster();
1081  if (cl) {
1082  Cls[d-d1] = cl->numberOfClusters();
1083  }
1084  }
1085  }
1086  Int_t np = 0;
1087  if (! be) np = emccol->barrelPoints().size();
1088  else np = emccol->endcapPoints().size();
1089  if (np ||
1090  Adcs[0] || Adcs[1] || Adcs[2] || Adcs[3] ||
1091  Cls[0] || Cls[1] || Cls[2] || Cls[3] ) {
1092  LOG_QA << Form("# %s points:%5i",Names[be],np);
1093  LOG_QA << Form(": Adc(T/p/E/P) %4i/%4i/%5i/%5i",Adcs[0],Adcs[1],Adcs[2],Adcs[3]);
1094  LOG_QA << Form(": Cls(T/p/E/P) %4i/%4i/%5i/%5i",Cls[0],Cls[1],Cls[2],Cls[3]);
1095  LOG_QA << endm;
1096  }
1097  }
1098  }
1099  const StBTofCollection* tof = event->btofCollection();
1100  if (tof) {
1101  const StSPtrVecBTofHit& tofHits = tof->tofHits();
1102  if (tofHits.size()) {
1103  Int_t n = tofHits.size();
1104  Int_t m = 0;
1105  for(Int_t i=0;i<n;i++) { //loop on hits in modules
1106  StBTofHit *aHit = tofHits[i];
1107  if(!aHit) continue;
1108  if (aHit->associatedTrack()) m++;
1109  }
1110  LOG_QA << Form("# BTof hits:%5i: Matched with tracks:%5i",n,m) << endm;
1111  }
1112  }
1113  const StPhmdCollection* pmdcol = event->phmdCollection();
1114  if (pmdcol) {
1115  const StPhmdDetector* pmd_det = pmdcol->detector(StDetectorId(kPhmdId));
1116  const StPhmdDetector* cpv_det = pmdcol->detector(StDetectorId(kPhmdCpvId));
1117  Int_t n = 0;
1118  Int_t m = 0;
1119  if (pmd_det) {
1120  const StPhmdClusterCollection* pmd_clusters = pmd_det->cluster();
1121  if (pmd_clusters) n = pmd_clusters->numberOfclusters();
1122  }
1123  if (cpv_det) {
1124  const StPhmdClusterCollection* cpv_clusters = cpv_det->cluster();
1125  if (cpv_clusters) m = cpv_clusters->numberOfclusters();
1126  }
1127  LOG_QA << Form("# Pmd clusters:%5i: Cpv clusters:%5i",n,m) << endm;
1128  }
1129  if (event->fpdCollection() && event->fpdCollection()->numberOfADC()) {
1130  UShort_t sum = 0;
1131  for (UInt_t i = 0; i < event->fpdCollection()->numberOfADC(); i++) sum += event->fpdCollection()->adc()[i];
1132  if (sum)
1133  LOG_QA << "# FPD ADC sum: " << sum << endm;
1134  }
1135  if (event->fgtCollection()) {
1136  LOG_QA << "# FGT hits: " << event->fgtCollection()->getNumHits() << endm;
1137  }
1138  if (event->richCollection() && event->richCollection()->getRichHits().size()) {
1139  LOG_QA << "# RICH hits: " << event->richCollection()->getRichHits().size() << endm;
1140  }
1141  if (event->numberOfPsds()) {
1142  LOG_QA << "# PSDs: " << event->numberOfPsds() << endm;
1143  }
1144 #ifdef _ST_GMT_HIT_H_
1145  if (event->gmtCollection() && event->gmtCollection()->getNumHits()) {
1146  LOG_QA << "# GMT hits: " << event->gmtCollection()->getNumHits()
1147  << " points: " << event->gmtCollection()->getNumPoints()
1148  << endm;
1149  }
1150 #endif /* _ST_GMT_HIT_H_ */
1151  if (event->rpsCollection()) {
1152  LOG_QA << "# RPS hits: " << event->rpsCollection()->clusters().size() << endm;
1153  }
1154 
1155  LOG_QA << "*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-" << endm;
1156 }
1157 //________________________________________________________________________________
1158 /* -------------------------------------------------------------------------
1159  * $Log: StAnalysisMaker.cxx,v $
1160  * Revision 2.26 2020/01/28 15:05:10 genevb
1161  * end-of-line needed for nShortTrackForETOF
1162  *
1163  * Revision 2.25 2020/01/27 21:28:30 genevb
1164  * Add short tracks toward ETOF
1165  *
1166  * Revision 2.24 2015/07/19 23:02:44 fisyak
1167  * Add print out for Sst, Gmt, pp2pp
1168  *
1169  * Revision 1.2 2014/01/15 23:12:46 fisyak
1170  * Freeze
1171  *
1172  * Revision 1.1.1.1 2013/08/29 14:21:56 fisyak
1173  * Freeze
1174  *
1175  * Revision 2.23 2012/12/18 17:16:26 fisyak
1176  * Add PrintVertex
1177  *
1178  * Revision 2.22 2012/11/25 22:22:45 fisyak
1179  * Add separators for summary
1180  *
1181  * Revision 2.21 2012/11/08 16:57:53 fisyak
1182  * Add pmd to summary
1183  *
1184  * Revision 2.20 2012/11/07 21:35:26 fisyak
1185  * Add to summary print out for EMC and ToF
1186  *
1187  * Revision 2.19 2012/10/23 19:44:18 fisyak
1188  * Add print out for ToF and Emc hits
1189  *
1190  * Revision 2.18 2012/09/16 21:59:14 fisyak
1191  * Compress print out, add PrintEmcHits
1192  *
1193  * Revision 2.17 2012/05/07 13:59:44 fisyak
1194  * enhance print out for primary vertixes
1195  *
1196  * Revision 2.16 2012/03/22 23:45:16 fisyak
1197  * Compress output for Event summary
1198  *
1199  * Revision 2.15 2010/09/01 14:33:57 fisyak
1200  * Clean ups
1201  *
1202  * Revision 2.14 2010/01/26 20:35:51 fisyak
1203  * use dca print out
1204  *
1205  * Revision 2.13 2009/11/23 15:54:28 fisyak
1206  * Clean-up est tracks
1207  *
1208  * Revision 2.12 2009/11/10 20:17:59 fisyak
1209  * Add print out for StEvent track and hits
1210  *
1211  * Revision 2.11 2009/11/03 15:13:22 fisyak
1212  * Comment print out, wait till StEvent will be mofidied
1213  *
1214  * Revision 2.10 2009/11/03 15:03:56 fisyak
1215  * Add static method to print StEvent
1216  *
1217  * Revision 2.9 2008/04/02 23:15:35 fisyak
1218  * Add protection against allGlobals == 0
1219  *
1220  * Revision 2.8 2004/02/04 01:36:40 jeromel
1221  * Minor change for user's education. Use of gMessMgr
1222  *
1223  * Revision 2.7 2004/02/01 18:01:53 jeromel
1224  * A few message addition
1225  *
1226  * Revision 2.6 2003/03/20 00:29:19 jeromel
1227  * Calling Wite() on 0x0 pointer
1228  *
1229  * Revision 2.5 2003/02/27 15:25:36 jeromel
1230  * Missing check on triggerIdCollection() now added
1231  *
1232  * Revision 2.4 2003/02/18 22:19:09 jeromel
1233  * Added dump of Y3 triggers
1234  *
1235  * Revision 2.3 2002/04/28 00:10:27 jeromel
1236  * doxygen basic dox added. GetCVS() had wrong signature : corrected to avoid
1237  * propagation of this typo in new makers.
1238  *
1239  * Revision 2.2 2000/07/12 05:23:28 ullrich
1240  * Updated for better use as template for actual analysis.
1241  *
1242  * Revision 2.1 1999/12/30 01:54:57 ogilvie
1243  * added countPrimaryPions as example how to use PID
1244  *
1245  * Revision 2.0 1999/11/04 16:10:03 ullrich
1246  * Revision for new StEvent
1247  *
1248  * -------------------------------------------------------------------------
1249  */
1250 
Definition: tof.h:15
const int kIstNumLadders
24 IST Ladders
A typical Analysis Class.
Collection of trigger ids as stored in StEvent.
Definition: Stypes.h:40
StAnalysisMaker(const Char_t *name="analysis")
The constructor. Initialize you data members here.
const int kIstNumSensorsPerLadder
6 sensor per one IST Ladder