StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StMuPmdUtil.cxx
1 /*****************************************************************
2  * $Id: StMuPmdUtil.cxx,v 1.4 2004/10/19 01:43:37 mvl Exp $
3  *
4  * Class : StMuPmdUtil
5  * Author: Supriya Das
6  * ****************************************************************
7  *
8  * Description: This is the utility class for PMD to convert
9  * StEvent to StMuDst and vice versa
10  * ****************************************************************
11  * $Log: StMuPmdUtil.cxx,v $
12  * Revision 1.4 2004/10/19 01:43:37 mvl
13  * Added code to copy hits
14  *
15  * Revision 1.3 2004/04/26 00:12:33 perev
16  * include StMuPmdCluster.h added
17  *
18  * Revision 1.2 2004/04/14 17:15:57 subhasis
19  * Xin's TOF reinclusion
20  *
21  * Revision 1.1 2004/04/02 03:36:21 jeromel
22  * New files for PMD
23  *
24  * ****************************************************************/
25 
26 #include "StMuPmdUtil.h"
27 #include "StEvent.h"
28 #include "StEventTypes.h"
29 #include "StMuPmdCluster.h"
30 #include "StMuPmdHit.h"
31 #include "StMuPmdCollection.h"
32 
33 #include "StPhmdDetector.h"
34 #include "StPhmdClusterCollection.h"
35 #include "StPhmdCluster.h"
36 #include "StPhmdModule.h"
37 #include "StPhmdHit.h"
38 
39 ClassImp(StMuPmdUtil)
40 
42 {
43 }
44 
45 StMuPmdUtil::~StMuPmdUtil()
46 {
47 }
48 
49 StMuPmdCollection* StMuPmdUtil::getMuPmd(StPhmdCollection *phmdColl)
50 {
51  if(!phmdColl) return NULL;
52  StMuPmdCollection* muPmd = new StMuPmdCollection();
53  fillMuPmd(phmdColl,muPmd);
54  return muPmd;
55 }
56 
57 StPhmdCollection* StMuPmdUtil::getPmd(StMuPmdCollection* muPmd)
58 {
59  if(!muPmd) return NULL;
60 
61  StPhmdCollection *phmdColl=new StPhmdCollection();
62  fillPmd(muPmd,phmdColl);
63  return phmdColl;
64 }
65 
66 void StMuPmdUtil::fillMuPmd(StPhmdCollection *phmdColl, StMuPmdCollection *muPmd)
67 {
68  if(!muPmd)return;
69  if(!phmdColl)return;
70 
71  for(Int_t d=0; d<2; d++){
72  Int_t PmdDet=d+1;
73 
74  Int_t totalHit=0;
75  StDetectorId pdet = static_cast<StDetectorId>(kPhmdId-d);
76  StPhmdDetector* detector = (StPhmdDetector*)phmdColl->detector(StDetectorId(pdet));
77 
78  if(detector) {
79  // Added to accommodate PmdHits : Supriya 11/10/04
80  for(UInt_t j=0;j<12;j++) {
81  StPhmdModule* module = detector->module(j);
82  // module=NULL;
83  if(module) {
84  //container for Hit
85  StSPtrVecPhmdHit& rawHit=module->hits();
86  if(rawHit.size()>0) {
87  //loop over Hits
88  totalHit+=rawHit.size();
89 
90  for(Int_t k=0;k<(Int_t)rawHit.size();k++) {
91 
92  muPmd->addHit(PmdDet); //adding to muPmd
93 
94  StPhmdHit *pmdhit = (StPhmdHit*)rawHit[k];
95  Int_t sm=pmdhit->superModule();
96  Int_t det=pmdhit->subDetector();
97  Int_t row=pmdhit->row();
98  Int_t col=pmdhit->column();
99  Float_t energy=pmdhit->energy();
100  Int_t adc=pmdhit->adc();
101  // Int_t vol=det*1000000+sm*100000+row*100+col;
102  //cout<<"sm,det,row,col,adc,energy "<<sm<<" "<<det<<" "<<row<<" "<<col<<" "<<adc<<" "<<energy<<endl;
103  StMuPmdHit *mupmdhit = muPmd->getHit(muPmd->getNHits(PmdDet)-1,PmdDet);
104 
105  mupmdhit->setEnergy(energy);
106  mupmdhit->setADC(adc);
107  mupmdhit->setSuperModule((Short_t)sm);
108  mupmdhit->setSubDetector((Short_t)det);
109  mupmdhit->setRow((Short_t)row);
110  mupmdhit->setColumn((Short_t)col);
111  } // Loop over Hits
112  } // checking rawHit.saize()>0
113  } // checking module existance
114  } //Loop over 12 Modules
115 
116  // Now fill clusters
117  StPhmdClusterCollection* cluscol = detector->cluster();
118 //cout<<"NHITS, nclust******** "<<totalHit<<" "<<cluscol->numberOfclusters()<<endl;
119 // cluscol=NULL;
120  if(cluscol)
121  {
122  Int_t Ncluster0 = cluscol->numberOfclusters();
123  //cout<<"Ncluster0, "<<Ncluster0<<" "<<PmdDet<<endl;
124  if(Ncluster0 > 0)
125  {
126  const StSPtrVecPhmdCluster& pmdclusters = cluscol->clusters();
127  //cout<<"cluster size "<<pmdclusters.size()<<endl;
128  for(UInt_t i=0; i<pmdclusters.size(); i++)
129  {
130  muPmd->addCluster(PmdDet); // adding to muPmd
131  //cout<<"clusterId "<<muPmd->getNClusters(PmdDet)<<endl;
132  StPhmdCluster *pmdcl = (StPhmdCluster*)pmdclusters[i];
133  Int_t sm = pmdcl->module();
134  Int_t ncell = pmdcl->numberOfCells();
135  Float_t eta = pmdcl->eta();
136  Float_t phi = pmdcl->phi();
137  Float_t sigma = pmdcl->sigma();
138  Float_t energy = pmdcl->energy();
139  Int_t energyPID = pmdcl->energyPid(); //PID based on energy cut
140  Int_t PID = pmdcl->pid(); //PID based on other methods like Neural Network
141  Int_t mcPID = pmdcl->mcPid(); //PID got from GEANT, 0 for Data
142 
143  StMuPmdCluster *mupmdcl=muPmd->getCluster(muPmd->getNClusters(PmdDet)-1,PmdDet);
144  // cout<<"size "<<i<<" "<<muPmd->getNClusters(PmdDet)<<" "<<mupmdcl<<endl;
145 
146  mupmdcl->setSuperModule(Int_t(sm));
147  mupmdcl->setNcell(Int_t(ncell));
148  mupmdcl->setEta(Float_t(eta));
149  mupmdcl->setPhi(Float_t(phi));
150  mupmdcl->setSigma(Float_t(sigma));
151  mupmdcl->setEnergy(Float_t(energy));
152  mupmdcl->setEnergyPID(Int_t(energyPID));
153  mupmdcl->setPID(Int_t(PID));
154  mupmdcl->setMcPID(Int_t(mcPID));
155 
156  } //Loop over clusters
157  } //Check for non-zero cluster.size()
158  } //Check for existence of ClusterCollection
159  } //Check for Detector
160  } //Loop over detectors
161 
162  return;
163 }
164 
165 void StMuPmdUtil::fillPmd(StMuPmdCollection* muPmd, StPhmdCollection* phmdColl)
166 {
167  if(!muPmd) return;
168  if(!phmdColl) return;
169 
170  for(Int_t d=0; d<2; d++) {
171  StDetectorId pdet = static_cast<StDetectorId>(kPhmdId-d);
172  StPhmdDetector* detector = (StPhmdDetector*)phmdColl->detector(StDetectorId(pdet));
173  // phmdColl->setDetector(detector);
174 
175  Int_t PmdDet=d+1;
176 
177  // Getting Clusters
178  Int_t nClusters = muPmd->getNClusters(PmdDet);
179  // cout<<"PmdDet, nclusters "<<PmdDet<<" "<<nClusters<<endl;
180  if(nClusters>0)
181  {
182  StPhmdClusterCollection* phmdClusColl = new StPhmdClusterCollection();
183  //cout<<"cluscoll "<<phmdClusColl<<endl;
184  //phmdClusColl->setDetector(pdet);
185  detector->setCluster(phmdClusColl);
186  //cout<<"cluscoll Set "<<phmdClusColl<<endl;
187 
188  for(Int_t n=0; n<nClusters; n++)
189  {
190  StMuPmdCluster *muPmdCl = muPmd->getCluster(n,PmdDet);
191  if(muPmdCl)
192  {
193  Int_t sm = muPmdCl->superModule();
194  Int_t nCell = muPmdCl->ncell();
195  Float_t eta = muPmdCl->eta();
196  Float_t phi = muPmdCl->phi();
197  Float_t energy = muPmdCl->energy();
198  Int_t energyPID = muPmdCl->energyPID();
199  Int_t PID = muPmdCl->PID();
200  Int_t McPID = muPmdCl->mcPID();
201 
202  StPhmdCluster* phmdCl = new StPhmdCluster();
203 
204  phmdCl->setModule(sm);
205  phmdCl->setNumberOfCells(nCell);
206  phmdCl->setEta(eta);
207  phmdCl->setPhi(phi);
208  phmdCl->setEnergy(energy);
209  phmdCl->setEnergyPid(energyPID);
210  phmdCl->setPid(PID);
211  phmdCl->setMcPid(McPID);
212 
213  phmdClusColl->addCluster(phmdCl);
214  }
215  }
216  }
217 
218  //Getting Hits
219 
220  Int_t nHits = muPmd->getNHits(PmdDet);
221 
222  if (nHits>0) {
223  for(Int_t nh=0; nh<nHits; nh++) {
224 
225  StMuPmdHit *muPmdHit = muPmd->getHit(nh,PmdDet);
226 
227  if(muPmdHit) {
228  Float_t energy= muPmdHit->energy();
229  Int_t ADC= muPmdHit->adc();
230  Int_t sm= muPmdHit->superModule();
231  Int_t subDet= muPmdHit->subDetector();
232  Int_t row= muPmdHit->row();
233  Int_t col= muPmdHit->column();
234  StPhmdHit* phmdHit = new StPhmdHit();
235 
236  phmdHit->setSuperModule(sm);
237  phmdHit->setEnergy(energy);
238  phmdHit->setAdc(ADC);
239  phmdHit->setSubDetector(subDet);
240  phmdHit->setRow(row);
241  phmdHit->setColumn(col);
242  detector->addHit(phmdHit);
243  }
244  }
245  }
246  }
247  return;
248 }
249 
250 
251 
252