StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StSvtClusterMaker.cxx
1 // $Id: StSvtClusterMaker.cxx,v 1.13 2007/04/28 17:57:05 perev Exp $
2 // $Log: StSvtClusterMaker.cxx,v $
3 // Revision 1.13 2007/04/28 17:57:05 perev
4 // Redundant StChain.h removed
5 //
6 // Revision 1.12 2007/03/21 17:22:58 fisyak
7 // Ivan Kotov's drift velocities, use TGeoHMatrix for coordinate transformation
8 //
9 // Revision 1.11 2005/08/04 04:06:38 perev
10 // clear of collection added
11 //
12 // Revision 1.10 2004/01/27 02:32:58 perev
13 // LeakOff
14 //
15 // Revision 1.9 2003/01/28 20:28:34 munhoz
16 // new filters for clusters
17 //
18 // Revision 1.8 2002/03/20 00:33:32 munhoz
19 // temporary fix for memory leaks and new vertex finder params for pp
20 //
21 // Revision 1.7 2001/10/06 00:09:00 caines
22 // Fix deleting that was already done by clear so no longer crashes on exit
23 //
24 // Revision 1.6 2001/09/22 01:07:09 caines
25 // Fixes now that AddData() is cleared everyevent
26 //
27 // Revision 1.5 2001/08/07 20:52:15 caines
28 // Implement better packing of svt hardware and charge values
29 //
30 // Revision 1.4 2001/04/29 20:11:58 caines
31 // Added reset command for Online monitor
32 //
33 // Revision 1.3 2000/08/29 22:46:26 caines
34 // Fixed some memory leaks
35 //
36 // Revision 1.2 2000/08/21 13:06:58 caines
37 // Much improved hit finding and fitting
38 //
39 // Revision 1.1 2000/07/06 03:50:34 caines
40 // First version of cluster finder and fitter
41 //
42 //
44 // //
45 // StSvtClusterMaker class //
46 // //
48 #include "TH1.h"
49 #include "TH2.h"
50 #include "St_DataSetIter.h"
51 #include "TObjectSet.h"
52 #include "StSequence.hh"
53 
54 #include "StMessMgr.h"
55 
56 #include "StSvtClassLibrary/StSvtHybridCollection.hh"
57 #include "StSvtClassLibrary/StSvtHybridData.hh"
58 #include "StSvtClassLibrary/StSvtData.hh"
59 #include "StSvtHybridCluster.hh"
60 #include "StSvtClusterMaker.h"
61 #include "StDbUtilities/St_svtRDOstrippedC.h"
62 
63 
64 
65 //_____________________________________________________________________________
66 StSvtClusterMaker::StSvtClusterMaker(const char *name) : StMaker(name)
67 {
68  mSvtEvent = NULL;
69  mHybridData = NULL;
70  mHybridCluster = NULL;
71  mClusterColl = NULL;
72  mClusterSet = NULL;
73  mClusterFinder = NULL;
74 }
75 //_____________________________________________________________________________
76 StSvtClusterMaker::~StSvtClusterMaker()
77 {
78 
79  delete mClusterFinder;
80 }
81 
82 //_____________________________________________________________________________
83 Int_t StSvtClusterMaker::Init(){
84 
85 
86  if (Debug()) gMessMgr->Debug() << "In StSvtClusterMaker::Init() ..." <<
87  GetName() << endm;
88 
89  if( GetSvtRawData()) gMessMgr->Warning() << "No SVT Raw data..." << endm;
90 
91  SetSvtCluster();
92  mClusterFinder = new StSvtClusterFinder();
93 
94  return StMaker::Init();
95 
96 }
97 
98 //_____________________________________________________________________________
99 Int_t StSvtClusterMaker::GetSvtRawData(){
100 
101  St_DataSet *dataSet = GetDataSet("StSvtData");
102  //assert(dataSet);
103  if( !dataSet) return kStWarn;
104 
105  mSvtEvent = (StSvtData*)(dataSet->GetObject());
106  if( !mSvtEvent) return kStWarn;
107 
108  return kStOk;
109 
110 
111 }
112 //_____________________________________________________________________________
113 Int_t StSvtClusterMaker::SetSvtCluster()
114 {
115  mClusterSet = new St_ObjectSet("StSvtCluster");
116  //AddData(mClusterSet);
117  AddConst(mClusterSet);
118  SetOutput(mClusterSet); //Declare for output
119 
120 
121  mClusterColl = new StSvtHybridCollection(mSvtEvent->getConfiguration());
122  mClusterSet->SetObject(mClusterColl);
123 
124  return kStOK;
125 }
126 
127 
128 //_____________________________________________________________________________
130 
131 
132  if (Debug()) gMessMgr->Debug() << "In StSvtClusterMaker::Make() ..." <<
133  GetName() << endm;
134 
135  if( GetSvtRawData()) {
136  gMessMgr->Warning() << " Problem with SVt raw in ClusterMaker" << endm;
137  return kStWarn;
138  }
139  // SetSvtCluster();
140  SetHybridClusters();
141 
142  return kStOK;
143 }
144 
145 //______________________________________________________________________________________
146 Int_t StSvtClusterMaker::SetHybridClusters()
147 {
148  int index =0;
149 
150  for(int barrel = 1;barrel <= mSvtEvent->getNumberOfBarrels();barrel++) {
151  //cout<<mSvtEvent->getNumberOfBarrels()<<endl;
152  for (int ladder = 1;ladder <= mSvtEvent->getNumberOfLadders(barrel);ladder++) {
153  //cout<<mSvtEvent->getNumberOfLadders(barrel)<<endl;
154  for (int wafer = 1;wafer <= mSvtEvent->getNumberOfWafers(barrel);wafer++) {
155  // cout<<mSvtEvent->getNumberOfWafers(barrel)<<endl;
156 
157  if (St_svtRDOstrippedC::instance()->svtRDOstrippedStatus(barrel,ladder,wafer)) continue; // bad or missing wafer
158 
159  for (int hybrid = 1;hybrid <=mSvtEvent->getNumberOfHybrids();hybrid++){
160  //cout<<mSvtEvent->getNumberOfHybrids()<<endl;
161 
162  index = mSvtEvent->getHybridIndex(barrel,ladder,wafer,hybrid);
163  if(index < 0) continue;
164  // if (index != 69) continue;
165 
166  mHybridData = (StSvtHybridData *)mSvtEvent->at(index);
167  if( !mHybridData) continue;
168 
169  mClusterFinder->setHybridPointer(mHybridData);
170  mClusterFinder->SetHybIndex(index);
171  mClusterFinder->ClusterFinder();
172 
173 
174  mHybridCluster = (StSvtHybridCluster* )mClusterColl->at(index);
175  if( mHybridCluster){
176  delete mHybridCluster;
177  mClusterColl->at(index)=0;
178  }
179 
180  mHybridCluster = new StSvtHybridCluster();
181 
182  mHybridCluster->setCluster(mClusterFinder);
183  mClusterColl->put_at(mHybridCluster,index);
184  mClusterFinder->ResetContainers();
185 
186  }
187  }
188  }
189  }
190 
191  return kStOK;
192 }
193 
194 
195 //_____________________________________________________________________________
196 
198 
199  if (Debug()) gMessMgr->Debug() << "In StSvtClusterMaker::Finish() ..." <<
200  GetName() << endm;
201 
202  return kStOK;
203 }
204 
205 //_____________________________________________________________________________
206 Int_t StSvtClusterMaker::Reset(){
207 
208  if (Debug()) gMessMgr->Debug() << "In StSvtClusterMaker::reset() ..." <<
209  GetName() << endm;
210 
211 //VP delete mHybridCluster;
212  mClusterSet->SetObject(0); //it will be deleted inside
213  delete mClusterFinder;
214 
215  mHybridCluster = NULL;
216  mClusterColl = NULL;
217  mClusterFinder = NULL;
218 
219  return kStOK;
220 }
221 //_____________________________________________________________________________
222 void StSvtClusterMaker::Clear(const char *)
223 {
224  if (mClusterColl) {
225  int n = mClusterColl->size();
226  mClusterColl->clear();
227  mClusterColl->resize(n);
228  }
229  StMaker::Clear();
230 }
231 //_____________________________________________________________________________
232 ClassImp(StSvtClusterMaker)
233 
234 
235 
236 
237 
238 
239 
240 
241 
242 
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: TObjectSet.h:59
virtual TObject * GetObject() const
The depricated method (left here for the sake of the backward compatibility)
Definition: TDataSet.cxx:428
virtual const char * GetName() const
special overload
Definition: StMaker.cxx:237
Definition: Stypes.h:42
Definition: Stypes.h:40
Definition: Stypes.h:41