StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StEmcVirtualFinder.cxx
1 
2 #include "StEmcVirtualFinder.h"
3 #include "StEvent.h"
4 #include "StEventTypes.h"
5 #include "TString.h"
6 #include "StEmcUtil/others/emcDetectorName.h"
7 
8 ClassImp(StEmcVirtualFinder)
9 
10 //_____________________________________________________________________________
12 {
13  TString name,title;
14  for(Int_t i = 0; i<MAXDETBARREL; i++)
15  {
16  mColl[i] = new StEmcPreClusterCollection(i+1);
17 
18  name = detname[i];
19  name += "_NClusters";
20  title = "Number of clusters for detector ";
21  title+= detname[i];
22  mHist1D[0][i] = new TH1F(name.Data(),title.Data(),500,0,500);
23 
24  name = detname[i];
25  name += "_NHits";
26  title = "Number of hits/cluster for detector ";
27  title+= detname[i];
28  mHist1D[1][i] = new TH1F(name.Data(),title.Data(),20,0,20);
29 
30  name = detname[i];
31  name += "_Energy";
32  title = "Cluster energy for detector ";
33  title+= detname[i];
34  mHist1D[2][i] = new TH1F(name.Data(),title.Data(),500,0,50);
35 
36  name = detname[i];
37  name += "_RMSEta";
38  title = "Cluster RMS/eta for detector ";
39  title+= detname[i];
40  mHist1D[3][i] = new TH1F(name.Data(),title.Data(),20,0,0.1);
41 
42  name = detname[i];
43  name += "_RMSPhi";
44  title = "Cluster RMS/phi for detector ";
45  title+= detname[i];
46  mHist1D[4][i] = new TH1F(name.Data(),title.Data(),20,0,0.1);
47 
48  name = detname[i];
49  name += "_EtaPhi";
50  title = "Cluster (eta,phi) for detector ";
51  title+= detname[i];
52  mHist2D[0][i] = new TH2F(name.Data(),title.Data(),40,-1,1,120,-3.1415,3.1415);
53  }
54 }
55 //_____________________________________________________________________________
56 StEmcVirtualFinder::~StEmcVirtualFinder()
57 {
58  for(Int_t i = 0; i<MAXDETBARREL; i++)
59  if(mColl[i])
60  delete mColl[i];
61 }
63 {
64  return kFALSE;
65 }
67 {
68  StEmcCollection *emc = event->emcCollection();
69  if(!emc)
70  return kFALSE;
71 
72  StSPtrVecEmcPoint& pvec = emc->barrelPoints();
73  if(pvec.size()>0)
74  pvec.clear();
75 
76  for(Int_t i=0; i<MAXDETBARREL; i++)
77  {
78  StDetectorId id = static_cast<StDetectorId>(i+kBarrelEmcTowerId);
79  StEmcDetector* detector=emc->detector(id);
80  if(detector)
81  if(detector->cluster())
82  {
83  StSPtrVecEmcCluster& cluster=detector->cluster()->clusters();
84  if(cluster.size()>0)
85  cluster.clear();
86  }
87  }
88  return kTRUE;
89 }
91 {
92  for(Int_t i = 0;i<MAXDETBARREL; i++)
93  mColl[i]->empty();
94  return kTRUE;
95 }
97 {
98  StEmcCollection *emc = event->emcCollection();
99  if(!emc)
100  return kFALSE;
101 
102  for(Int_t i = 0;i<MAXDETBARREL;i++)
103  if(mColl[i])
104  {
105  StDetectorId id = static_cast<StDetectorId>(i+kBarrelEmcTowerId);
106  StEmcDetector* detector = emc->detector(id);
107  if(detector)
108  {
109  Int_t n = mColl[i]->getNClusters();
110  if(n>0)
111  {
112  StEmcClusterCollection* coll = detector->cluster();
113  if(!coll)
114  {
115  coll = new StEmcClusterCollection();
116  detector->setCluster(coll);
117  }
118  coll->setDetector(id);
119  coll->setClusterFinderId(1);
120  coll->setClusterFinderParamVersion(1);
121 
122  for(Int_t j = 0;j<n;j++)
123  {
124  StEmcPreCluster *pc = mColl[i]->getCluster(j);
125  if(pc)
126  {
127  StEmcCluster *c = pc->makeStCluster();
128  if(c)
129  coll->addCluster(c);
130  }
131  }
132  }
133  else
134  detector->setCluster(NULL);
135  }
136  }
137  return kTRUE;
138 }
140 {
141  StEmcCollection *emc = event->emcCollection();
142  if(!emc)
143  return kFALSE;
144  for(Int_t i = 0;i<MAXDETBARREL;i++)
145  {
146  StDetectorId id = static_cast<StDetectorId>(i+kBarrelEmcTowerId);
147  StEmcDetector* detector = emc->detector(id);
148  if(detector)
149  {
150  StEmcClusterCollection* coll = detector->cluster();
151  if(coll)
152  {
153  StSPtrVecEmcCluster& clusters = coll->clusters();
154  Int_t n = clusters.size();
155  mHist1D[0][i]->Fill(n);
156  for(Int_t j = 0;j<n;j++)
157  {
158  StEmcCluster *c =clusters[j];
159  if(c)
160  {
161  mHist1D[1][i]->Fill(c->nHits());
162  mHist1D[2][i]->Fill(c->energy());
163  mHist1D[3][i]->Fill(c->sigmaEta());
164  mHist1D[4][i]->Fill(c->sigmaPhi());
165  mHist2D[0][i]->Fill(c->eta(),c->phi());
166  }
167  }
168  }
169  }
170  }
171  return kFALSE;
172 }
virtual Bool_t clear()
clear the pre cluster collections
Int_t getNClusters()
gets the number of clusters in the collection
virtual Bool_t fillStEvent(StEvent *)
fills the StEvent object with the StEmcPreCluster objects in the collections
virtual Bool_t fillHistograms(StEvent *)
fills the QA histograms
virtual Bool_t findClusters(StEvent *)
finds clusters in a StEvent object
StEmcCluster * makeStCluster()
creates an StEmcCluster from the information in this pre cluster
StEmcPreCluster * getCluster(Int_t)
gets a cluster in the collection by its index