StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StSsdPointMaker.cxx
1 // $Id: StSsdPointMaker.cxx,v 1.70 2017/04/26 20:13:03 perev Exp $
2 //
3 // $Log: StSsdPointMaker.cxx,v $
4 // Revision 1.70 2017/04/26 20:13:03 perev
5 // Hide m_DataSet
6 //
7 // Revision 1.69 2015/08/06 17:46:53 smirnovd
8 // Removed unused local variables
9 //
10 // Revision 1.68 2015/05/21 18:49:50 bouchet
11 // arrayIndex out of bounds fixed
12 //
13 // Revision 1.67 2015/05/21 02:54:00 perev
14 // bug #3106 By mistake, fixing memory leak was removed by me line mySsd->Reset()
15 // Now I put it back
16 // Victor
17 //
18 // Revision 1.66 2015/05/19 19:11:18 perev
19 // Avoid leak by useng standard array instead of new and delete
20 //
21 // Revision 1.65 2015/05/19 13:44:17 bouchet
22 // ticket 3105 + warning fixed when printing the size of gain table
23 //
24 // Revision 1.64 2015/05/15 18:31:10 bouchet
25 // possible infinite loop fixed
26 //
27 // Revision 1.63 2012/11/07 21:54:37 fisyak
28 // Remove check for .histos
29 //
30 // Revision 1.62 2009/02/18 21:31:22 bouchet
31 // fix bug for printing the number of hits/ladder and filling histograms
32 //
33 // Revision 1.61 2008/10/20 19:32:39 bouchet
34 // use of new writePointToContainer method for the hit quality calculation
35 //
36 // Revision 1.60 2008/07/16 21:01:57 bouchet
37 // calculation of hits quality removed : call of default writePointToContainer
38 //
39 // Revision 1.59 2008/05/20 03:05:54 bouchet
40 // fix improper STAR logger(#1185) ; thanks to Valeri
41 //
42 // Revision 1.58 2008/05/07 22:45:24 bouchet
43 // add mcEvent dependence for embedding
44 //
45 // Revision 1.57 2008/04/15 21:05:22 bouchet
46 // remove latest change
47 //
48 // Revision 1.56 2008/04/12 14:20:38 bouchet
49 // Add a switch to use constant noise and pedestal ; remove some printing
50 //
51 // Revision 1.55 2008/01/15 13:48:58 bouchet
52 // Set a default value for uninitialized variable
53 //
54 // Revision 1.54 2008/01/11 10:39:39 bouchet
55 // add method to read the Wafer configuration table
56 //
57 // Revision 1.53 2007/09/25 13:40:46 bouchet
58 // Use m_Mode to switch between pedestals used in real data/simulation ; move some message to DEBUG
59 //
60 // Revision 1.52 2007/08/20 06:47:37 bouchet
61 // ssdStripCalib table taken for simulation
62 //
63 // Revision 1.51 2007/07/14 14:29:44 bouchet
64 // forget the Debug condition for the declaration of the tuples
65 //
66 // Revision 1.50 2007/07/14 13:52:16 bouchet
67 // add method to fill with default pedestal/noise values if no table is found
68 //
69 // Revision 1.49 2007/07/13 06:19:43 bouchet
70 // display of number of reconstructed hits corrected
71 //
72 // Revision 1.48 2007/07/12 17:07:18 bouchet
73 // add switch to read old ssdStripCalib Table and new ssdNoise Table
74 //
75 // Revision 1.47 2007/07/02 20:01:03 bouchet
76 // bug fixed for the normalization of reconstruction efficiency histos
77 //
78 // Revision 1.46 2007/07/01 16:18:41 bouchet
79 // add a normalization for the reconstruction efficiency histograms
80 //
81 // Revision 1.45 2007/06/23 04:53:50 bouchet
82 // add 0's to Timestamp which size is less than 6 digits
83 //
84 // Revision 1.44 2007/06/19 18:30:24 bouchet
85 // Add a method to evaluate the reconstruction efficiency (defined as the ratio of the number of matched clusters with all reconstructed clusters) ; some clean-up
86 //
87 // Revision 1.43 2007/06/19 01:19:15 bouchet
88 // cosmetic changes
89 //
90 // Revision 1.42 2007/04/28 17:56:58 perev
91 // Redundant StChain.h removed
92 //
93 // Revision 1.41 2007/04/17 05:09:25 perev
94 // GetTFile()==>StMaker. Jerome request
95 //
96 // Revision 1.40 2007/04/11 22:45:22 perev
97 // 1/0 avoided
98 //
99 // Revision 1.39 2007/03/27 23:15:09 bouchet
100 // Add a switch to use the gain calibration
101 //
102 // Revision 1.38 2007/03/27 18:30:04 fisyak
103 // recover lost access to ssdStripCalib table
104 //
105 // Revision 1.37 2007/03/21 17:19:12 fisyak
106 // use TGeoHMatrix for coordinate transformation, eliminate ssdWafersPostion, ake NTuples only for Debug()>1
107 //
108 // Revision 1.36 2007/03/08 23:04:42 bouchet
109 // add WriteMatchedStrips() method : fill the characteristics of the strips from matched clusters ; Small change for the writing of tuples
110 //
111 // Revision 1.35 2007/03/01 22:19:21 bouchet
112 // add a protection when ssdStripCalib is filled with empty values
113 //
114 // Revision 1.34 2007/02/21 20:36:17 bouchet
115 // add a method WriteMatchedClusters :\ instead of WriteScfTuple() method that fill all the reconstructed clusters,\ this one store the clusters associated to the hits
116 //
117 // Revision 1.33 2007/02/15 14:40:27 bouchet
118 // bug fixed for filling makeScmCtrlHistograms() method
119 //
120 // Revision 1.32 2007/02/14 11:49:12 bouchet
121 // Added control histograms and updated the Cluster and Point Tuple
122 //
123 // Revision 1.31 2007/02/02 20:24:15 bouchet
124 // WriteStripTuple method added, WriteScmTuple method updated
125 //
126 // Revision 1.30 2007/02/02 17:46:58 bouchet
127 // Few changes for the new Logger
128 //
129 // Revision 1.29 2007/01/16 18:01:52 bouchet
130 // Replace printf,cout,gMessMgr with LOG statements
131 //
132 // Revision 1.28 2006/10/16 16:27:49 bouchet
133 // Unify classes ; Methods for all classes (StSsdStrip, StSsdCluster, StSsdPoint) are now in StSsdUtil
134 //
135 // Revision 1.27 2006/09/15 21:03:14 bouchet
136 // id_mctrack is using for setIdTruth and propagated to the hit
137 //
138 // Revision 1.26 2005/12/31 01:43:22 perev
139 // Mack/Upack simplified
140 //
141 // Revision 1.25 2005/12/20 10:35:51 lmartin
142 // ReadStrip method updated and some cosmetic changes
143 //
144 // Revision 1.24 2005/12/20 09:23:35 lmartin
145 // ssdStripCalib table read from the mysql db
146 //
147 // Revision 1.23 2005/09/30 14:28:30 lmartin
148 // add a 0 to myTime if GetTime()<100000
149 //
150 // Revision 1.22 2005/09/26 15:49:54 bouchet
151 // adding a method to the poInt_t maker to check which ssdStripCalib is picked
152 //
153 // Revision 1.21 2005/08/11 13:51:39 lmartin
154 // PrintStripDetails, PrintPackageDetails and PrintPointDetails methods added
155 //
156 // Revision 1.20 2005/06/24 10:19:46 lmartin
157 // preventing crashes if ssdStripCalib is missing
158 //
159 // Revision 1.19 2005/06/16 14:29:22 bouchet
160 // no more makeSsdPedestalHistograms() method
161 //
162 // Revision 1.18 2005/06/14 12:09:15 bouchet
163 // add a histo for the pedestal and new name of the class : SsdPoint
164 //
165 // Revision 1.14 2005/06/07 16:24:47 lmartin
166 // InitRun returns kStOk
167 //
168 // Revision 1.13 2005/06/07 12:04:46 reinnart
169 // Make Stuff moved to Initrun
170 //
171 // Revision 1.12 2005/06/07 11:55:08 reinnart
172 // Initrun and good database connection
173 //
174 // Revision 1.11 2005/04/25 14:13:23 bouchet
175 // new method makeScfCtrlHistograms and makeScmCtrlHistograms and Clusternoise is coded as a float
176 //
177 // Revision 1.10 2005/04/23 08:56:20 lmartin
178 // physics and pedestal data processing separated
179 //
180 // Revision 1.9 2005/03/23 16:07:26 lmartin
181 // PrintClusterSummary and PrintPointSummary methods added
182 //
183 // Revision 1.8 2005/03/22 13:46:43 lmartin
184 // PrintStripSummary method added
185 //
186 // Revision 1.7 2005/03/22 10:38:51 lmartin
187 // HighCut value taken from the db
188 //
189 // Revision 1.6 2005/03/18 13:35:40 lmartin
190 // Partly missing cvs header added
191 //
192 // Revision 1.5 2005/03/18 10:16:34 lmartin
193 // positionSize argument added to the initLadders method
194 //
195 // Revision 1.4 2004/11/04 15:10:19 croy
196 // use the IAttr(".histos") to control histogramming and modification of the SsdHitCollection creation
197 //
198 // Revision 1.3 2004/08/13 07:07:23 croy
199 // Updates to read SSD databases
200 //
201 // Revision 1.2 2004/07/20 14:04:02 croy
202 // Use of new database structure definitions related to SSD config
203 //
204 // Revision 1.1 2004/03/12 06:12:37 jeromel
205 // Peer review closed. Yuri/Frank.
206 //
207 // Revision 1.3 2002/03/25 20:13:05 hippolyt
208 // Merged the two former makers
209 //
210 //
211 #include "StSsdPointMaker.h"
212 
213 #include "TFile.h"
214 #include "TH1.h"
215 #include "TH2.h"
216 #include "TDataSetIter.h"
217 #include "StMessMgr.h"
218 #include "TNtuple.h"
219 #include "StSsdUtil/StSsdPoint.hh"
220 #include "StSsdUtil/StSsdPackage.hh"
221 #include "StSsdUtil/StSsdCluster.hh"
222 #include "StSsdUtil/StSsdStripList.hh"
223 #include "StSsdUtil/StSsdClusterList.hh"
224 #include "StSsdUtil/StSsdPointList.hh"
225 #include "StSsdUtil/StSsdPackageList.hh"
226 #include "StSsdUtil/StSsdWafer.hh"
227 #include "StSsdUtil/StSsdLadder.hh"
228 #include "StSsdUtil/StSsdBarrel.hh"
229 #include "StSsdUtil/StSsdStrip.hh"
230 #include "tables/St_spa_strip_Table.h"
231 #include "tables/St_ssdPedStrip_Table.h"
232 #include "tables/St_scf_cluster_Table.h"
233 #include "tables/St_scm_spt_Table.h"
234 #include "tables/St_slsCtrl_Table.h"
235 #include "tables/St_clusterControl_Table.h"
236 #include "tables/St_ssdDimensions_Table.h"
237 #include "tables/St_ssdConfiguration_Table.h"
238 #include "tables/St_ssdWafersPosition_Table.h"
239 #include "tables/St_ssdLaddersPosition_Table.h"
240 #include "tables/St_ssdSectorsPosition_Table.h"
241 #include "tables/St_ssdBarrelPosition_Table.h"
242 #include "tables/St_ssdStripCalib_Table.h"
243 #include "tables/St_ssdGainCalibWafer_Table.h"
244 #include "tables/St_ssdNoise_Table.h"
245 #include "tables/St_ssdWaferConfiguration_Table.h"
246 #include "StEvent.h"
247 #include "StSsdHitCollection.h"
248 #include "StSsdDbMaker/StSsdDbMaker.h"
249 #include "TMath.h"
250 ClassImp(StSsdPointMaker);
251 
252 //_____________________________________________________________________________
253 Int_t StSsdPointMaker::Init(){
254  LOG_INFO << "Init() : Defining the histograms" << endm;
255  noisDisP = new TH1F("Noise_p","Noise Distribution",250,0,25);
256  snRatioP = new TH1F("SN_p","Signal/Noise (p)",200,0,200);
257  stpClusP = new TH1F("NumberOfStrips_p","Strips per Cluster",8,0,8);
258  totChrgP = new TH1F("ChargeElectron_p","Total Cluster Charge",100,0,300000);
259  ClusNvsClusP = new TH2S("ClusNvsClusP","Number of clusters on the n-side vs Number of clusters on the p-side",200,0,200,200,0,200);
260  ClusNvsClusP->SetXTitle("Number of p-Side Clusters");
261  ClusNvsClusP->SetYTitle("Number of n-Side Clusters");
262  noisDisN = new TH1F("Noise_n","Noise Distribution",250,0,25);
263  snRatioN = new TH1F("SN_n","Signal/Noise",200,0,200);
264  stpClusN = new TH1F("NumberOfStrips_n","Strips per Cluster",8,0,8);
265  totChrgN = new TH1F("ChargeElectron_n","Total Cluster Charge",100,0,300000);
266  ClustMapP = new TH2S("ClustMapP","Number of clusters on the p-side per wafer and ladder",20,0,20,16,0,16);
267  ClustMapP->SetXTitle("Ladder id");
268  ClustMapP->SetYTitle("Wafer id");
269  ClustMapN = new TH2S("ClustMapN","Number of clusters on the n-side per wafer and ladder",20,0,20,16,0,16);
270  ClustMapN->SetXTitle("Ladder id");
271  ClustMapN->SetYTitle("Wafer id");
272  MatchedClusterP = new TH2F("MatchedClusterP","#frac{# clusters matched}{# clusters reconstructed} , wafers on p-side",20,1,21,16,1,17);
273  MatchedClusterP->SetXTitle("Ladder id");
274  MatchedClusterP->SetYTitle("Wafer id");
275  MatchedClusterN = new TH2F("MatchedClusterN","#frac{# clusters matched}{# clusters reconstructed} , wafers on n-side",20,1,21,16,1,17);
276  MatchedClusterN->SetXTitle("Ladder id");
277  MatchedClusterN->SetYTitle("Wafer id");
278  // Create SCM histograms
279  matchisto = new TH2S("matchingHisto","Matching Adc (1p-1n)",500,0,1000,500,0,1000);
280  matchisto->SetXTitle("PSide ADC count");
281  matchisto->SetYTitle("NSide ADC count");
282  matchisto->SetZTitle("(1p-1n) hits");
283 
284  matchisto->SetTitleOffset(2,"X");
285  matchisto->SetTitleOffset(2,"Y");
286  // matchisto->SetTitleOffset(-1,"Z");
287 
288  matchisto->SetLabelSize(0.03,"X");
289  matchisto->SetLabelSize(0.03,"Y");
290  matchisto->SetLabelSize(0.03,"Z");
291 
292  matchisto->SetNdivisions(5,"X");
293  matchisto->SetNdivisions(5,"Y");
294  matchisto->SetNdivisions(10,"Z");
295 
296  orthoproj = new TH1S("ProjectionOrtho","Perfect Matching Deviation",320,-80,80);
297 
298  kind = new TH1S("kind","Kind of hits",11,0,11);
299  kind->SetXTitle("kind");
300  kind->SetYTitle("entries");
301  kind->SetTitleOffset(2,"X");
302  kind->SetTitleOffset(2,"Y");
303 
304  TString Title;
305  Char_t *Name = new Char_t[20];
306  Title ="Matching Adc (1p-1n) for ladder";
307  for(Int_t ii=0;ii<20;ii++) {
308  Title = Form("Matching Adc (1p-1n) for ladder = %i",ii+1);
309  sprintf(Name,"%s%d","matchingHisto_",ii);
310  matchisto_[ii] = new TH2S(Name,Title,500,0,1000, 500, 0, 1000);
311  matchisto_[ii]->SetXTitle("PSide ADC count");
312  matchisto_[ii]->SetYTitle("NSide ADC count");
313  matchisto_[ii]->SetZTitle("(1p-1n) hits");
314  }
315  if (Debug()>1) DeclareNtuple();
316  return StMaker::Init();
317 }
318 //_____________________________________________________________________________
319 Int_t StSsdPointMaker::InitRun(Int_t runumber) {
320  // mDbMgr = StDbManager::Instance();
321  // mDbMgr->setVerbose(false);
322 
323  // maccess = mDbMgr->initConfig(dbGeometry,dbSsd);
324  mode= gStSsdDbMaker->GetMode();
325  LOG_INFO <<"m_Mode = " << mode << endm;
326  NEvent = 0;
327  UseCalibration = 1;
328  UseWaferConfig = 1;
329  LOG_INFO<<Form("UseCalibration =%d UseWaferTable = %d",UseCalibration,UseWaferConfig)<<endm;
330  St_slsCtrl* slsCtrlTable = (St_slsCtrl*) GetDataBase("Geometry/ssd/slsCtrl");
331  if(! slsCtrlTable){LOG_ERROR << "InitRun : No access to slsCtrl table" << endm;}
332  else {
333  mDynamicControl = new StSsdDynamicControl();
334  slsCtrl_st* control = (slsCtrl_st*) slsCtrlTable->GetTable();
335  mDynamicControl -> setnElectronInAMip(control->nElectronInAMip);
336  mDynamicControl -> setadcDynamic(control->adcDynamic);
337  mDynamicControl -> seta128Dynamic(control->a128Dynamic);
338  mDynamicControl -> setnbitEncoding(control->nbitEncoding);
339  mDynamicControl -> setnstripInACluster(control->nstripInACluster);
340  mDynamicControl -> setpairCreationEnergy(control->pairCreationEnergy);
341  mDynamicControl -> setparDiffP(control->parDiffP);
342  mDynamicControl -> setparDiffN(control->parDiffN);
343  mDynamicControl -> setparIndRightP(control->parIndRightP);
344  mDynamicControl -> setparIndRightN(control->parIndRightN);
345  mDynamicControl -> setparIndLeftP(control->parIndLeftP);
346  mDynamicControl -> setparIndLeftN(control->parIndLeftN);
347  mDynamicControl -> setdaqCutValue(control->daqCutValue);
348  mDynamicControl -> printParameters();
349  }
350  St_clusterControl* clusterCtrlTable = (St_clusterControl*) GetDataBase("Geometry/ssd/clusterControl");
351  if (!clusterCtrlTable) {LOG_ERROR << "InitRun : No access to clusterControl table" << endm;}
352  else {
353  mClusterControl = new StSsdClusterControl();
354  clusterControl_st *clusterCtrl = (clusterControl_st*) clusterCtrlTable->GetTable() ;
355  mClusterControl -> setHighCut(clusterCtrl->highCut);
356  mClusterControl -> setTestTolerance(clusterCtrl->testTolerance);
357  mClusterControl -> setClusterTreat(clusterCtrl->clusterTreat);
358  mClusterControl -> setAdcTolerance(clusterCtrl->adcTolerance);
359  mClusterControl -> setMatchMean(clusterCtrl->matchMean);
360  mClusterControl -> setMatchSigma(clusterCtrl->matchSigma);
361  mClusterControl -> printParameters();
362  }
363  year = (GetDate()/10000)-2000;
364  LOG_DEBUG <<Form("TimeStamp is %d Year is =%d\n",GetDate(),year)<<endm;
365  switch(mode)
366  {
367  case 1: {
368  m_noise2 = (St_ssdStripCalib*) GetDataBase("Calibrations/ssd/ssdStripCalib");
369  if (!m_noise2) {LOG_ERROR << "InitRun : No access to ssdStripCalib - will use the default noise and pedestal values" << endm;}
370  else {
371  LOG_INFO<<"InitRun for simu : old Table (ssdStripCalib) is used"<<endm;
372  }
373  break;
374  }
375  case 0 :{
376  if(year<7){
377  m_noise2 = (St_ssdStripCalib*) GetDataBase("Calibrations/ssd/ssdStripCalib");
378  if (!m_noise2) {LOG_ERROR << "InitRun : No access to ssdStripCalib - will use the default noise and pedestal values" << endm;}
379  else {
380  LOG_INFO<<"InitRun for real data : old Table(ssdStripCalib) is used"<<endm;
381  }
382  }
383  else {
384  m_noise3 = (St_ssdNoise*)GetDataBase("Calibrations/ssd/ssdNoise");
385  if (!m_noise3) {LOG_ERROR << "InitRun : No access to ssdNoise - will use the default noise and pedestal values" << endm;}
386  else{
387  LOG_INFO << "InitRun for real data : new Table(ssdNoise) is used" << endm;}
388  }
389  break;
390  }
391  default : {printf("no real data nor simu");}
392  }
393  (UseCalibration==1)?FillCalibTable():FillDefaultCalibTable();
394  (UseWaferConfig==1)?FillWaferTable():FillDefaultWaferTable();
395  /*
396  Init arrays for the reconstruction efficiency
397  */
398  for(Int_t ii=0 ;ii<20;ii++)
399  {
400  for(Int_t jj=0;jj<16;jj++)
401  {
402  ratioP[ii][jj] = 0;
403  ratioN[ii][jj] = 0;
404  }
405  }
406  return kStOk;
407 }
408 //_____________________________________________________________________________
409 void StSsdPointMaker::DeclareNtuple(){
410  TFile *f = GetTFile();
411  if (f){
412  f->cd();
413  string varlist2 = "pulseP:pulseN:ladder:wafer:case:xg:yg:zg:flag:idClusP:idClusN:position_0:position_1:xl:yl";
414  mHitNtuple = new TNtuple("PhysNTuple","Physics Ntuple",varlist2.c_str());
415  string varlist3 = "side:ladder:wafer:nstrip:snratio:noise:first_strip:TotAdc:FirstAdc:LastAdc:TotNoise";
416  nHitNtuple = new TNtuple("ClusTuple","All Clusters stored",varlist3.c_str());
417  string varlist4 = "side:ladder:wafer:nstrip:pedestal:signal:noise:snratio";
418  qHitNtuple = new TNtuple("Strips","All Strips stored",varlist4.c_str());
419  pHitNtuple = new TNtuple("ClustupleIn","Clusters in hits",varlist3.c_str());
420  rHitNtuple = new TNtuple("StripsIn","Strips in hits",varlist4.c_str());
421  }
422 }
423 //_____________________________________________________________________________
425 {
426  LOG_DEBUG << Form("Make : begin")<< endm;
427  // Create output tables
428  Int_t res = 0;
429  char myLabel[100];
430  char myTime[100];
431  char myDate[100];
432  if (GetTime()<999)
433  sprintf(myTime,"000%d",GetTime());
434  else
435  if ((GetTime()<9999)&&(GetTime()>999))
436  sprintf(myTime,"00%d",GetTime());
437  else
438  if ((GetTime()<99999)&&(GetTime()>9999))
439  sprintf(myTime,"0%d",GetTime());
440  else
441  sprintf(myTime,"%d",GetTime());
442  sprintf(myDate,"%d%s",GetDate(),".");
443  sprintf(myLabel,"%s%s",myDate,myTime);
444  // two different tables can exist (physics data or pedestal data)
445 
446  TDataSet *SpaStrip = GetDataSet("SpaStrip");
447  if (! SpaStrip) {
448  LOG_ERROR << "Make : no input data set, wrong chain option" << endm;
449  return kStErr;
450  }
451  St_spa_strip *spa_strip = dynamic_cast<St_spa_strip *> (SpaStrip->Find("spa_strip"));
452  St_ssdPedStrip *spa_ped_strip = dynamic_cast<St_ssdPedStrip *> (SpaStrip->Find("ssdPedStrip"));
453 
454  if (!spa_strip || spa_strip->GetNRows()==0){
455  {
456  LOG_WARN << "Make : no input (fired strip for the SSD)"<<endm;
457  LOG_WARN <<"Make : looking for a pedestal/noise tables"<<endm;
458  }
459  if (!spa_ped_strip || spa_ped_strip->GetNRows()==0) {
460  LOG_WARN<<"Make : no pedestal/noise data..."<<endm;
461  return kStWarn;
462  }
463  else
464  { LOG_WARN<<"Make : pedestal/noise data found : "<<spa_ped_strip->GetNRows()<<endm;}
465  }
466 
467  St_scm_spt *scm_spt = new St_scm_spt("scm_spt",5000);
468  AddData(scm_spt);
469 
470  St_scf_cluster *scf_cluster = new St_scf_cluster("scf_cluster",5000);//09/13
471  AddData(scf_cluster);
472 
473  mCurrentEvent = (StEvent*) GetInputDS("StEvent");
474  if(mCurrentEvent)
475  {
476  mSsdHitColl = mCurrentEvent->ssdHitCollection();
477  if (!mSsdHitColl) {
478  LOG_WARN << "Make : The SSD hit collection does not exist - creating a new one" << endm;
479  mSsdHitColl = new StSsdHitCollection;
480  mCurrentEvent->setSsdHitCollection(mSsdHitColl);
481  }
482  }
483  else
484  mSsdHitColl = 0;
485 
486  LOG_INFO<<"#################################################"<<endm;
487  LOG_INFO<<"#### START OF NEW SSD POINT MAKER ####"<<endm;
488  LOG_INFO<<"#### SSD BARREL INITIALIZATION ####"<<endm;
489  LOG_INFO<<"#### BEGIN INITIALIZATION ####"<<endm;
490  StSsdBarrel *mySsd =gStSsdDbMaker->GetSsd();
491  mySsd->setClusterControl(mClusterControl);
492  //The full SSD object is built only if we are processing physics data
493  if((! spa_ped_strip || spa_ped_strip->GetNRows()==0) && (spa_strip->GetNRows()!=0))
494  {
495  Int_t stripTableSize = mySsd->readStripFromTable(spa_strip);
496  LOG_INFO<<"#### NUMBER OF SPA STRIPS "<<stripTableSize<<" ####"<<endm;
497  mySsd->sortListStrip();
498  PrintStripSummary(mySsd);
499  noiseTableSize = 0;
500  noiseTableSize = ReadNoiseTable(mySsd,year);
501  LOG_INFO<<"#### NUMBER OF DB ENTRIES "<<noiseTableSize<<" ####"<<endm;
502  Int_t nClusterPerSide[2];
503  nClusterPerSide[0] = 0;
504  nClusterPerSide[1] = 0;
505  mySsd->doSideClusterisation(nClusterPerSide,WafStatus);
506  LOG_INFO<<"#### NUMBER OF CLUSTER P SIDE "<<nClusterPerSide[0]<<" ####"<<endm;
507  LOG_INFO<<"#### NUMBER OF CLUSTER N SIDE "<<nClusterPerSide[1]<<" ####"<<endm;
508  mySsd->sortListCluster();
509  Int_t nClusterWritten = mySsd->writeClusterToTable(scf_cluster,spa_strip);
510  LOG_INFO<<"#### -> "<<nClusterWritten<<" CLUSTERS WRITTEN INTO TABLE ####"<<endm;
511  PrintClusterSummary(mySsd);
512  //PrintStripDetails(mySsd,8310);
513  //PrintClusterDetails(mySsd,8310);
514  makeScfCtrlHistograms(mySsd);
515  //debugUnPeu(mySsd);
516  Int_t nPackage = mySsd->doClusterMatching(CalibArray);
517  LOG_INFO<<"#### -> "<<nPackage<<" PACKAGES IN THE SSD ####"<<endm;
518  mySsd->convertDigitToAnalog(mDynamicControl);
519  mySsd->convertUFrameToOther();
520  PrintPointSummary(mySsd);
521  //Int_t nSptWritten = mySsd->writePointToContainer(scm_spt,mSsdHitColl);
522  if(Debug()){
523  for(Int_t i=1;i<=20;i++)
524  {
525  for(Int_t j=1;j<=16;j++)
526  {
527  //PrintStripDetails(mySsd,7000+(100*j)+i);
528  //PrintClusterDetails(mySsd,7000+(100*j)+i);
529  //PrintPointDetails(mySsd,7000+(100*j)+i);
530  //PrintPackageDetails(mySsd,7000+(100*j)+i);
531  }
532  }
533  }
534  //get McEvent here
535  Int_t nSptWritten = 0;
536  StMcEvent* mcEvent = 0;
537  mcEvent = (StMcEvent*) GetDataSet("StMcEvent");
538  if(mcEvent)
539  {
540  LOG_DEBUG << " mcEvent exists " << endm;
541  nSptWritten = mySsd->writePointToContainer(scm_spt,mSsdHitColl,scf_cluster,spa_strip,mDynamicControl,mcEvent);
542  }
543  else{
544  nSptWritten = mySsd->writePointToContainer(scm_spt,mSsdHitColl,scf_cluster);
545  }
546  LOG_INFO<<"#### -> "<<nSptWritten<<" HITS WRITTEN INTO TABLE ####"<<endm;
547  if(mSsdHitColl){
548  if (mSsdHitColl->numberOfHits()>0) {
549  NEvent++;
550  LOG_INFO<<"#### -> "<<mSsdHitColl->numberOfHits()<<" HITS WRITTEN INTO CONTAINER ####"<<endm;
551  makeScmCtrlHistograms(mySsd);
552  EvaluateEfficiency(mySsd);
553  NormalizeEfficiency();
554  scm_spt->Purge();
555  }
556  else {
557  LOG_INFO<<" ######### NO SSD HITS WRITTEN INTO CONTAINER ####"<<endm;
558  }
559  }
560  LOG_INFO<<"#### END OF SSD NEW POINT MAKER ####"<<endm;
561  LOG_INFO<<"#################################################"<<endm;
562  if(Debug() >1){
563  if (qHitNtuple) WriteStripTuple(mySsd);
564  if (nHitNtuple) WriteScfTuple(mySsd);
565  if (mHitNtuple) WriteScmTuple(mySsd);
566  if (rHitNtuple) WriteMatchedStrips(mySsd);
567  if (pHitNtuple) WriteMatchedClusters(mySsd);
568  }
569  if (nSptWritten) res = kStOK;
570  }
571  else
572  {
573  if((spa_strip->GetNRows()==0)&&(spa_ped_strip && spa_ped_strip->GetNRows()!=0))
574  {
575  LOG_INFO <<"###### WRITING SSD PEDESTAL HISTOGRAMS##########"<<endm;
576  if(year<7){
577  mySsd->writeNoiseToFile(spa_ped_strip,myLabel);}
578  else{mySsd->writeNewNoiseToFile3(spa_ped_strip,myLabel);//new method
579  printf("done\n");
580  }
581  }
582  }
583  //clear stuff
584  mySsd->Reset();
585  if(res!=kStOK){
586  LOG_WARN <<"Make : no output" << endm;;
587  return kStWarn;
588  }
589  return kStOK;
590 }
591 //_____________________________________________________________________________
592 void StSsdPointMaker::makeScfCtrlHistograms(StSsdBarrel *mySsd)
593 {
594  //Int_t LadderIsActive[20]={1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1};
595  Float_t convAdcToE = (mDynamicControl->getadcDynamic()*mDynamicControl->getnElectronInAMip())/(pow(2.0,mDynamicControl->getnbitEncoding()));
596  Int_t ClustersP_tot = 0;
597  Int_t ClustersN_tot = 0;
598  Int_t pSize = 0;
599  Int_t nSize = 0;
600  for (Int_t i=0;i<20;i++)
601  if (mySsd->isActiveLadder(i)>0) {
602  for (Int_t j=0; j<mySsd->mLadders[i]->getWaferPerLadder();j++) {
603  StSsdClusterList *pList = mySsd->mLadders[i]->mWafers[j]->getClusterP();
604  pSize = pList->getSize();
605  ClustersP_tot+= pSize;
606  StSsdClusterList *nList = mySsd->mLadders[i]->mWafers[j]->getClusterN();
607  nSize = nList->getSize();
608  ClustersN_tot+= nSize;
609  ClusNvsClusP->Fill(pSize,nSize);
610  ClustMapP->Fill(i,j,pSize);
611  ClustMapN->Fill(i,j,nSize);
612  pSize = 0;
613  nSize = 0;
614  StSsdCluster *pClusterP = mySsd->mLadders[i]->mWafers[j]->getClusterP()->first();
615  while (pClusterP)
616  {
617  stpClusP->Fill(pClusterP->getClusterSize());
618  totChrgP->Fill(convAdcToE*pClusterP->getTotAdc());
619  if (pClusterP->getClusterSize()>0)
620  noisDisP->Fill(pClusterP->getTotNoise()/pClusterP->getClusterSize());
621  if (pClusterP->getTotNoise()>0)
622  snRatioP->Fill((pClusterP->getTotAdc()*pClusterP->getClusterSize())/pClusterP->getTotNoise());
623  pClusterP = mySsd->mLadders[i]->mWafers[j]->getClusterP()->next(pClusterP);
624  }
625  StSsdCluster *pClusterN = mySsd->mLadders[i]->mWafers[j]->getClusterN()->first();
626  while (pClusterN)
627  {
628  stpClusN->Fill(pClusterN->getClusterSize());
629  totChrgN->Fill(convAdcToE*pClusterN->getTotAdc());
630  noisDisN->Fill(pClusterN->getTotNoise()/(3e-33+pClusterN->getClusterSize()));
631  snRatioN->Fill((pClusterN->getTotAdc()*pClusterN->getClusterSize())/(3e-33+pClusterN->getTotNoise()));
632  pClusterN = mySsd->mLadders[i]->mWafers[j]->getClusterN()->next(pClusterN);
633  }
634  }
635  }
636  LOG_DEBUG <<"totclusters P="<<ClustersP_tot<<" totclusters N="<<ClustersN_tot<<endm;
637 }
638 //_____________________________________________________________________________
639 
640 void StSsdPointMaker::makeScmCtrlHistograms(StSsdBarrel *mySsd)
641 {
642  //Int_t LadderIsActive[20]={1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1};
643  Int_t conversion[11]={11,12,21,13,31,221,222,223,23,32,33};
644  Float_t convMeVToAdc = (int)pow(2.0,mDynamicControl->getnbitEncoding())/(mDynamicControl->getpairCreationEnergy()*mDynamicControl->getadcDynamic()*mDynamicControl->getnElectronInAMip());
645  for (Int_t i=0;i<20;i++)
646  if (mySsd->isActiveLadder(i)>0) {
647  for (Int_t j=0; j<mySsd->mLadders[i]->getWaferPerLadder();j++) {
648  if (mySsd->mLadders[i]->mWafers[j]->getPoint()->getSize()==0) {
649  //LOG_INFO <<"PrintPointDetails() - No hit in this wafer "<< endm;
650  }
651  else {
652  StSsdPoint *pSpt = mySsd->mLadders[i]->mWafers[j]->getPoint()->first();
653  while (pSpt){
654  if (pSpt->getNMatched() == 11)// case 11
655  {
656  Float_t a = 0, b = 0;
657  a = convMeVToAdc*(pSpt->getDe(0)+pSpt->getDe(1));
658  b = convMeVToAdc*(pSpt->getDe(0)-pSpt->getDe(1));
659  matchisto->Fill(a,b);
660  orthoproj->Fill((b-a)/TMath::Sqrt(2.));
661  matchisto_[i]->Fill(a,b);
662  }
663 
664  for(Int_t k=0;k<11;k++)
665  {
666  if(pSpt->getNMatched()==conversion[k])
667  {
668  kind->Fill(k);
669  }
670  }
671  pSpt = mySsd->mLadders[i]->mWafers[j]->getPoint()->next(pSpt);
672  }
673  }
674  }
675  }
676 }
677 //_____________________________________________________________________________
678 
679 void StSsdPointMaker::PrintStripSummary(StSsdBarrel *mySsd)
680 {
681  Int_t ladderCountN[20]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0} ;
682  Int_t ladderCountP[20]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0} ;
683  for (Int_t i=0;i<20;i++)
684  if (mySsd->isActiveLadder(i)>0) {
685  for (Int_t j=0; j<mySsd->mLadders[i]->getWaferPerLadder();j++) {
686  ladderCountP[i]=ladderCountP[i]+mySsd->mLadders[i]->mWafers[j]->getStripP()->getSize();
687  ladderCountN[i]=ladderCountN[i]+mySsd->mLadders[i]->mWafers[j]->getStripN()->getSize();
688  }
689  }
690 
691  LOG_INFO <<"PrintStripSummary : Number of raw data in the SSD" << endm;
692  LOG_INFO << "PrintStripSummary : Active Ladders : ";
693  for (Int_t i=0;i<20;i++)
694  if (mySsd->isActiveLadder(i)>0) {
695  LOG_DEBUG.width(5);
696  LOG_DEBUG<<i+1;
697  }
698 
699  LOG_DEBUG<<endm;
700  LOG_INFO << "PrintStripSummary : Counts (p-side): ";
701  for (Int_t i=0;i<20;i++)
702  if (mySsd->isActiveLadder(i)>0) {
703  LOG_DEBUG.width(5);
704  LOG_DEBUG <<ladderCountP[i];
705  }
706  LOG_DEBUG<<endm;
707  LOG_INFO << "PrintStripSummary : Counts (n-side): ";
708  for (Int_t i=0;i<20;i++)
709  if (mySsd->isActiveLadder(i)>0) {
710  LOG_DEBUG.width(5);
711  LOG_DEBUG <<ladderCountN[i];
712  }
713  LOG_DEBUG<<endm;
714 }
715 
716 //_____________________________________________________________________________
717 void StSsdPointMaker::debugUnPeu(StSsdBarrel *mySsd)
718 {
719  Int_t monladder,monwafer;
720  monladder=7;
721  monwafer=6;
722  mySsd->debugUnPeu(monladder,monwafer);
723 }
724 
725 //_____________________________________________________________________________
726 void StSsdPointMaker::PrintClusterSummary(StSsdBarrel *mySsd)
727 {
728  Int_t ladderCountN[20]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0} ;
729  Int_t ladderCountP[20]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0} ;
730  for (Int_t i=0;i<20;i++)
731  if (mySsd->isActiveLadder(i)>0) {
732  for (Int_t j=0; j<mySsd->mLadders[i]->getWaferPerLadder();j++) {
733  ladderCountP[i]=ladderCountP[i]+mySsd->mLadders[i]->mWafers[j]->getClusterP()->getSize();
734  ladderCountN[i]=ladderCountN[i]+mySsd->mLadders[i]->mWafers[j]->getClusterN()->getSize();
735  }
736  }
737 
738  LOG_INFO <<"PrintClusterSummary : Number of clusters in the SSD" << endm;
739  LOG_INFO << "PrintClusterSummary : Active Ladders : ";
740  for (Int_t i=0;i<20;i++)
741  if (mySsd->isActiveLadder(i)>0) {
742  LOG_DEBUG.width(5);
743  LOG_DEBUG<<i+1;
744  }
745 
746  LOG_DEBUG<<endm;
747  LOG_INFO << "PrintClusterSummary : Counts (p-side): ";
748  for (Int_t i=0;i<20;i++)
749  if (mySsd->isActiveLadder(i)>0) {
750  LOG_DEBUG.width(5);
751  LOG_DEBUG <<ladderCountP[i];
752  }
753  LOG_DEBUG<<endm;
754  LOG_INFO << "PrintClusterSummary : Counts (n-side): ";
755  for (Int_t i=0;i<20;i++)
756  if (mySsd->isActiveLadder(i)>0) {
757  LOG_DEBUG.width(5);
758  LOG_DEBUG <<ladderCountN[i];
759  }
760  LOG_DEBUG<<endm;
761 }
762 //_____________________________________________________________________________
763 void StSsdPointMaker::PrintPointSummary(StSsdBarrel *mySsd)
764 {
765  Int_t ladderCount[20]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0} ;
766  Int_t ladderCount11[20]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0} ;
767  for (Int_t i=0;i<20;i++)
768  if (mySsd->isActiveLadder(i)>0) {
769  for (Int_t j=0; j<mySsd->mLadders[i]->getWaferPerLadder();j++) {
770  ladderCount[i]=ladderCount[i]+mySsd->mLadders[i]->mWafers[j]->getPoint()->getSize();
771  StSsdPoint *pSpt = mySsd->mLadders[i]->mWafers[j]->getPoint()->first();
772  while (pSpt){
773  if (pSpt->getNMatched()==11) ladderCount11[i]++;
774  pSpt = mySsd->mLadders[i]->mWafers[j]->getPoint()->next(pSpt);
775  }
776  }
777  }
778 
779  LOG_INFO<<"PrintPointSummary : Number of hits in the SSD" << endm;
780  LOG_INFO<< "PrintPointSummary : Active Ladders : ";
781  for (Int_t i=0;i<20;i++)
782  if (mySsd->isActiveLadder(i)>0) {
783  LOG_DEBUG.width(5);
784  LOG_DEBUG<<i+1;
785  }
786 
787  LOG_DEBUG<<endm;
788  LOG_INFO << "PrintPointSummary : Counts : ";
789  for (Int_t i=0;i<20;i++)
790  if (mySsd->isActiveLadder(i)>0) {
791  LOG_DEBUG.width(5);
792  LOG_DEBUG <<ladderCount[i];
793  }
794  LOG_DEBUG<<endm;
795  LOG_INFO << "PrintPointSummary : Counts (11) : ";
796  for (Int_t i=0;i<20;i++)
797  if (mySsd->isActiveLadder(i)>0) {
798  LOG_DEBUG.width(5);
799  LOG_DEBUG <<ladderCount11[i];
800  }
801  LOG_DEBUG<<endm;
802 }
803 //_____________________________________________________________________________
804 void StSsdPointMaker::WriteStripTuple(StSsdBarrel *mySsd)
805 {
806  //Int_t LadderIsActive[20]={1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1} ;
807  for (Int_t i=0;i<20;i++)
808  if (mySsd->isActiveLadder(i)>0) {
809  for (Int_t j=0; j<mySsd->mLadders[i]->getWaferPerLadder();j++) {
810  StSsdStrip *pStripP = mySsd->mLadders[i]->mWafers[j]->getStripP()->first();
811  while (pStripP){
812  Strips_hits[0] = 0;
813  Strips_hits[1] = i+1;
814  Strips_hits[2] = j+1;
815  Strips_hits[3] = pStripP->getNStrip();
816  Strips_hits[4] = pStripP->getPedestal();
817  Strips_hits[5] = pStripP->getDigitSig();
818  Strips_hits[6] = pStripP->getSigma();
819  Strips_hits[7] = (float)(pStripP->getDigitSig()/pStripP->getSigma());
820  qHitNtuple->Fill(Strips_hits);
821  pStripP = mySsd->mLadders[i]->mWafers[j]->getStripP()->next(pStripP);
822  }
823  StSsdStrip *pStripN = mySsd->mLadders[i]->mWafers[j]->getStripN()->first();
824  while (pStripN){
825  Strips_hits[0] = 1;
826  Strips_hits[1] = i+1;
827  Strips_hits[2] = j+1;
828  Strips_hits[3] = pStripN->getNStrip();
829  Strips_hits[4] = pStripN->getPedestal();
830  Strips_hits[5] = pStripN->getDigitSig();
831  Strips_hits[6] = pStripN->getSigma();
832  Strips_hits[7] = (float)(pStripN->getDigitSig()/pStripN->getSigma());
833  qHitNtuple->Fill(Strips_hits);
834  pStripN = mySsd->mLadders[i]->mWafers[j]->getStripN()->next(pStripN);
835  }
836  }
837  }
838 }
839 
840 //_____________________________________________________________________________
841 void StSsdPointMaker::WriteScfTuple(StSsdBarrel *mySsd)
842 {
843  //Int_t LadderIsActive[20]={1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1};
844  for (Int_t i=0;i<20;i++)
845  if (mySsd->isActiveLadder(i)>0) {
846  for (Int_t j=0; j<mySsd->mLadders[i]->getWaferPerLadder();j++) {
847  StSsdCluster *pClusterP = mySsd->mLadders[i]->mWafers[j]->getClusterP()->first();
848  while (pClusterP)
849  {
850  ClusterNtuple[0] = 0;
851  ClusterNtuple[1] = i+1;
852  ClusterNtuple[2] = j+1;
853  ClusterNtuple[3] = pClusterP->getClusterSize();
854  ClusterNtuple[4] = ((pClusterP->getTotAdc()*pClusterP->getClusterSize())/pClusterP->getTotNoise());
855  ClusterNtuple[5] = pClusterP->getTotNoise()/pClusterP->getClusterSize();
856  ClusterNtuple[6] = pClusterP->getFirstStrip();
857  ClusterNtuple[7] = pClusterP->getTotAdc();
858  ClusterNtuple[8] = pClusterP->getFirstAdc();
859  ClusterNtuple[9] = pClusterP->getLastAdc();
860  ClusterNtuple[10]= pClusterP->getTotNoise();
861  pClusterP = mySsd->mLadders[i]->mWafers[j]->getClusterP()->next(pClusterP);
862  nHitNtuple->Fill(ClusterNtuple);
863  }
864  StSsdCluster *pClusterN = mySsd->mLadders[i]->mWafers[j]->getClusterN()->first();
865  while (pClusterN)
866  {
867  ClusterNtuple[0] = 1;
868  ClusterNtuple[1] = i+1;
869  ClusterNtuple[2] = j+1;
870  ClusterNtuple[3] = pClusterN->getClusterSize();
871  ClusterNtuple[4] = ((pClusterN->getTotAdc()*pClusterN->getClusterSize())/pClusterN->getTotNoise());
872  ClusterNtuple[5] = pClusterN->getTotNoise()/pClusterN->getClusterSize();
873  ClusterNtuple[6] = pClusterN->getFirstStrip();
874  ClusterNtuple[7] = pClusterN->getTotAdc();
875  ClusterNtuple[8] = pClusterN->getFirstAdc();
876  ClusterNtuple[9] = pClusterN->getLastAdc();
877  ClusterNtuple[10]= pClusterN->getTotNoise();
878  pClusterN = mySsd->mLadders[i]->mWafers[j]->getClusterN()->next(pClusterN);
879  nHitNtuple->Fill(ClusterNtuple);
880  }
881  }
882  }
883 }
884 //_____________________________________________________________________________
885 void StSsdPointMaker::WriteScmTuple(StSsdBarrel *mySsd)
886 {
887  //Int_t LadderIsActive[20]={1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1};
888  Int_t conversion[11]={11,12,21,13,31,221,222,223,23,32,33};
889  Float_t convMeVToAdc = (int)pow(2.0,mDynamicControl->getnbitEncoding())/(mDynamicControl->getpairCreationEnergy()*mDynamicControl->getadcDynamic()*mDynamicControl->getnElectronInAMip());
890  for (Int_t i=0;i<20;i++)
891  if (mySsd->isActiveLadder(i)>0) {
892  for (Int_t j=0; j<mySsd->mLadders[i]->getWaferPerLadder();j++) {
893  if (mySsd->mLadders[i]->mWafers[j]->getPoint()->getSize()==0) {
894  }
895  else {
896  StSsdPoint *pSpt = mySsd->mLadders[i]->mWafers[j]->getPoint()->first();
897  while (pSpt){
898  Float_t a = 0, b = 0;
899  a = convMeVToAdc*(pSpt->getDe(0)+pSpt->getDe(1));
900  b = convMeVToAdc*(pSpt->getDe(0)-pSpt->getDe(1));
901  hitNtuple[0] = a;
902  hitNtuple[1] = b;
903  hitNtuple[2] = i+1;
904  hitNtuple[3] = j+1;
905  for(Int_t k=0;k<11;k++)
906  {
907  if(pSpt->getNMatched()==conversion[k])
908  {
909  hitNtuple[4]=k;
910  }
911  }
912  hitNtuple[5] = pSpt->getXg(0);
913  hitNtuple[6] = pSpt->getXg(1);
914  hitNtuple[7] = pSpt->getXg(2);
915  hitNtuple[8] = pSpt->getFlag();
916 
917  Int_t IdP = pSpt->getIdClusterP();
918  Int_t IdN = pSpt->getIdClusterN();
919 
920  StSsdClusterList *currentListP_j = mySsd->mLadders[i]->mWafers[j]->getClusterP();
921  StSsdCluster *cluster_P_j = currentListP_j->first();
922  while(cluster_P_j)
923  {
924  if(cluster_P_j->getNCluster()==IdP)
925  break;
926  cluster_P_j = currentListP_j->next(cluster_P_j);
927  }
928 
929  StSsdClusterList *currentListN_j = mySsd->mLadders[i]->mWafers[j]->getClusterN();
930  StSsdCluster *cluster_N_j = currentListN_j->first();
931  while(cluster_N_j)
932  {
933  if(cluster_N_j->getNCluster()==IdN)
934  break;
935  cluster_N_j = currentListN_j->next(cluster_N_j);
936  }
937 
938  hitNtuple[9] = cluster_P_j->getStripMean();
939  hitNtuple[10]= cluster_N_j->getStripMean();
940  hitNtuple[11] = pSpt->getPositionU(0);
941  hitNtuple[12] = pSpt->getPositionU(1);
942  hitNtuple[13] = pSpt->getXl(0);
943  hitNtuple[14] = pSpt->getXl(1);
944  mHitNtuple->Fill(hitNtuple);
945  pSpt = mySsd->mLadders[i]->mWafers[j]->getPoint()->next(pSpt);
946  }
947  }
948  }
949  }
950 }
951 //_____________________________________________________________________________
952 void StSsdPointMaker::PrintStripDetails(StSsdBarrel *mySsd, Int_t mywafer)
953 {
954  //Int_t LadderIsActive[20]={1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1} ;
955  Int_t found = 0 ;
956  LOG_DEBUG <<"PrintStripDetails() - Wafer "<<mywafer<< endm;
957  for (Int_t i=0;i<20;i++)
958  if (mySsd->isActiveLadder(i)>0) {
959  for (Int_t j=0; j<mySsd->mLadders[i]->getWaferPerLadder();j++) {
960  if (mySsd->mLadders[i]->mWafers[j]->getIdWafer()==mywafer) {
961  found=1;
962  //Looking for the P-side strip informations
963  if (mySsd->mLadders[i]->mWafers[j]->getStripP()->getSize()==0) {
964  LOG_DEBUG <<"PrintStripDetails() - No strip on the P-side of this wafer "<< endm;
965  }
966  else {
967  LOG_DEBUG<<"PrintStripDetails() - "
968  <<mySsd->mLadders[i]->mWafers[j]->getStripP()->getSize()<<" strip(s) on the P-side of this wafer "<< endm;
969  LOG_DEBUG<<"PrintStripDetails() - Strip/Adc/Ped/Noise/Analog"<< endm;
970  StSsdStrip *pStripP = mySsd->mLadders[i]->mWafers[j]->getStripP()->first();
971  while (pStripP){
972  LOG_DEBUG<<"PrintStripDetails() - "
973  <<pStripP->getNStrip()<<" "
974  <<pStripP->getDigitSig()<<" "
975  <<pStripP->getPedestal()<<" "
976  <<pStripP->getSigma()<<" "
977  <<pStripP->getAnalogSig()<<" "
978  <<endm;
979  for(Int_t e=0;e<5;e++){printf("e=%d idMcHit=%d idMcTrack=%d\n",e,pStripP->getIdMcHit(e),pStripP->getIdMcTrack(e));}
980  pStripP = mySsd->mLadders[i]->mWafers[j]->getStripP()->next(pStripP);
981  }
982  }
983  //Looking for the N-side strip informations
984  if (mySsd->mLadders[i]->mWafers[j]->getStripN()->getSize()==0) {
985  LOG_DEBUG <<"PrintStripDetails() - No strip on the N-side of this wafer "<< endm;
986  }
987  else {
988  LOG_DEBUG<<"PrintStripDetails() - "
989  <<mySsd->mLadders[i]->mWafers[j]->getStripN()->getSize()<<" strip(s) on the N-side of this wafer "<< endm;
990  LOG_DEBUG <<"StSsdPointMaker::PrintStripDetails() - Strip/Adc/Ped/Noise/Analog"<< endm;
991  StSsdStrip *pStripN = mySsd->mLadders[i]->mWafers[j]->getStripN()->first();
992  while (pStripN){
993  LOG_DEBUG<<"PrintStripDetails() - "
994  <<pStripN->getNStrip()<<" "
995  <<pStripN->getDigitSig()<<" "
996  <<pStripN->getPedestal()<<" "
997  <<pStripN->getSigma()<<" "
998  <<pStripN->getAnalogSig()<<" "
999  <<endm;
1000  for(Int_t e=0;e<5;e++){printf("e=%d idMcHit=%d idMcTrack=%d\n",e,pStripN->getIdMcHit(e),pStripN->getIdMcTrack(e));}
1001  pStripN = mySsd->mLadders[i]->mWafers[j]->getStripN()->next(pStripN);
1002  }
1003  }
1004  }
1005  }
1006  }
1007  if (found==0) {LOG_DEBUG <<"PrintStripDetails() - Wafer not found !!!"<<endm;}
1008 }
1009 //_____________________________________________________________________________
1010 void StSsdPointMaker::PrintClusterDetails(StSsdBarrel *mySsd, Int_t mywafer)
1011 {
1012  //Int_t LadderIsActive[20]={1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1} ;
1013  Int_t found = 0;
1014  LOG_INFO <<"PrintClusterDetails() - Wafer "<<mywafer<< endm;
1015  for (Int_t i=0;i<20;i++)
1016  if (mySsd->isActiveLadder(i)>0) {
1017  for (Int_t j=0; j<mySsd->mLadders[i]->getWaferPerLadder();j++) {
1018  if (mySsd->mLadders[i]->mWafers[j]->getIdWafer()==mywafer) {
1019  found=1;
1020  //Looking for the P-side cluster informations
1021  if (mySsd->mLadders[i]->mWafers[j]->getClusterP()->getSize()==0) {
1022  LOG_INFO <<"PrintClusterDetails() - No cluster on the P-side of this wafer "<< endm;
1023  }
1024  else {
1025  LOG_INFO<<"PrintClusterDetails() - "
1026  <<mySsd->mLadders[i]->mWafers[j]->getClusterP()->getSize()<<" cluster(s) on the P-side of this wafer "<< endm;
1027  LOG_INFO<<"PrintClusterDetails() - Cluster/Flag/Size/1st Strip/Strip Mean/TotAdc/1st Adc/Last Adc/TotNoise"<< endm;
1028  StSsdCluster *pClusterP = mySsd->mLadders[i]->mWafers[j]->getClusterP()->first();
1029  while (pClusterP){
1030  LOG_INFO<<"PrintClusterDetails() - "
1031  <<pClusterP->getNCluster()<<" "
1032  <<pClusterP->getFlag()<<" "
1033  <<pClusterP->getClusterSize()<<" "
1034  <<pClusterP->getFirstStrip()<<" "
1035  <<pClusterP->getStripMean()<<" "
1036  <<pClusterP->getTotAdc()<<" "
1037  <<pClusterP->getFirstAdc()<<" "
1038  <<pClusterP->getLastAdc()<<" "
1039  <<pClusterP->getTotNoise()<<" "
1040  <<endm;
1041  for(Int_t e=0;e<5;e++){printf("e=%d idMcHit=%d \n",e,pClusterP->getIdMcHit(e));}
1042  pClusterP = mySsd->mLadders[i]->mWafers[j]->getClusterP()->next(pClusterP);
1043  }
1044  }
1045  //Looking for the N-side cluster informations
1046  if (mySsd->mLadders[i]->mWafers[j]->getClusterN()->getSize()==0) {
1047  LOG_INFO <<"PrintClusterDetails() - No cluster on the N-side of this wafer "<< endm;
1048  }
1049  else {
1050  LOG_INFO<<"PrintClusterDetails() - "
1051  <<mySsd->mLadders[i]->mWafers[j]->getClusterN()->getSize()<<" cluster(s) on the N-side of this wafer "<< endm;
1052  LOG_INFO<<"PrintClusterDetails() - Cluster/Flag/Size/1st Strip/Strip Mean/TotAdc/1st Adc/Last Adc/TotNoise"<< endm;
1053  StSsdCluster *pClusterN = mySsd->mLadders[i]->mWafers[j]->getClusterN()->first();
1054  while (pClusterN){
1055  LOG_INFO<<"PrintClusterDetails() - "
1056  <<pClusterN->getNCluster()<<" "
1057  <<pClusterN->getFlag()<<" "
1058  <<pClusterN->getClusterSize()<<" "
1059  <<pClusterN->getFirstStrip()<<" "
1060  <<pClusterN->getStripMean()<<" "
1061  <<pClusterN->getTotAdc()<<" "
1062  <<pClusterN->getFirstAdc()<<" "
1063  <<pClusterN->getLastAdc()<<" "
1064  <<pClusterN->getTotNoise()<<" "
1065  <<endm;
1066  for(Int_t e=0;e<5;e++){printf("e=%d idMcHit=%d \n",e,pClusterN->getIdMcHit(e));}
1067  pClusterN = mySsd->mLadders[i]->mWafers[j]->getClusterN()->next(pClusterN);
1068  }
1069  }
1070  }
1071  }
1072  }
1073  if (found==0){ LOG_INFO <<"PrintClusterDetails() - Wafer not found !!!"<<endm; }
1074 }
1075 
1076 //_____________________________________________________________________________
1077 void StSsdPointMaker::PrintPackageDetails(StSsdBarrel *mySsd, Int_t mywafer)
1078 {
1079  //Int_t LadderIsActive[20]={1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1} ;
1080  Int_t found;
1081  found=0;
1082  LOG_INFO <<"PrintPackageDetails() - Wafer "<<mywafer<< endm;
1083  for (Int_t i=0;i<20;i++)
1084  if (mySsd->isActiveLadder(i)>0) {
1085  for (Int_t j=0; j<mySsd->mLadders[i]->getWaferPerLadder();j++) {
1086  if (mySsd->mLadders[i]->mWafers[j]->getIdWafer()==mywafer) {
1087  found=1;
1088  if (mySsd->mLadders[i]->mWafers[j]->getPackage()->getSize()==0) {
1089  LOG_INFO <<"PrintPackageDetails() - No package in this wafer "<< endm;
1090  }
1091  else {
1092  LOG_INFO <<"PrintPackageDetails() - "<<mySsd->mLadders[i]->mWafers[j]->getPackage()->getSize()<<" package(s) in this wafer "<< endm;
1093  LOG_INFO <<"PrintPackageDetails() - Package/Kind/Size"<< endm;
1094  StSsdPackage *pPack = mySsd->mLadders[i]->mWafers[j]->getPackage()->first();
1095  while (pPack){
1096  LOG_INFO<<"PrintPackageDetails() - "<<pPack->getNPackage()<<" "
1097  <<pPack->getKind()<<" "
1098  <<pPack->getSize()<<" "<<endm;
1099  for (Int_t k=0;k<pPack->getSize();k++) {
1100  LOG_INFO<<"PrintPackageDetails() - "<<k<<" "<<pPack->getMatched(k)<<" "<<pPack->getMatched(k)->getNCluster()<<endm;
1101  }
1102  pPack = mySsd->mLadders[i]->mWafers[j]->getPackage()->next(pPack);
1103  }
1104  }
1105  }
1106  }
1107  }
1108  if (found==0){ LOG_INFO <<"PrintPackageDetails() - Wafer not found !!!"<<endm;}
1109 }
1110 
1111 //_____________________________________________________________________________
1112 void StSsdPointMaker::PrintPointDetails(StSsdBarrel *mySsd, Int_t mywafer)
1113 {
1114  //Int_t LadderIsActive[20]={1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1} ;
1115  Int_t found;
1116  found=0;
1117  Float_t convMeVToAdc = (int)pow(2.0,mDynamicControl->getnbitEncoding())/(mDynamicControl->getpairCreationEnergy()*mDynamicControl->getadcDynamic()*mDynamicControl->getnElectronInAMip());
1118  LOG_INFO <<"PrintPointDetails() - Wafer "<<mywafer<< endm;
1119  for (Int_t i=0;i<20;i++)
1120  if (mySsd->isActiveLadder(i)>0) {
1121  for (Int_t j=0; j<mySsd->mLadders[i]->getWaferPerLadder();j++) {
1122  if (mySsd->mLadders[i]->mWafers[j]->getIdWafer()==mywafer) {
1123  found=1;
1124  if (mySsd->mLadders[i]->mWafers[j]->getPoint()->getSize()==0) {
1125  LOG_INFO <<"PrintPointDetails() - No hit in this wafer "<< endm;
1126  }
1127  else {
1128  LOG_INFO<<"PrintPointDetails() - "<<mySsd->mLadders[i]->mWafers[j]->getPoint()->getSize()<<" hit(s) in this wafer "<< endm;
1129 
1130  LOG_INFO<<"PrintPointDetails() - Hit/Flag/NMatched/IdClusP/IdClusN/idMcHit[0]/idMcHit[1]/idMcHit[2]/idMcHit[3]/idMcHit[4]/Xg[0]/Xg[1]/Xg[2]/Xl[0]/Xl[1]/Xl[2]/a/b"<<endm;
1131  StSsdPoint *pSpt = mySsd->mLadders[i]->mWafers[j]->getPoint()->first();
1132  while (pSpt){
1133  Float_t a = 0, b = 0;
1134  a = convMeVToAdc*(pSpt->getDe(0)+pSpt->getDe(1));
1135  b = convMeVToAdc*(pSpt->getDe(0)-pSpt->getDe(1));
1136  LOG_INFO<<"PrintPointDetails() - "
1137  <<pSpt->getNPoint() <<" "
1138  <<pSpt->getFlag() <<" "
1139  <<pSpt->getNMatched() <<" "
1140  <<pSpt->getIdClusterP()<<" "
1141  <<pSpt->getIdClusterN()<<" "
1142  <<pSpt->getNMchit(0) <<" "
1143  <<pSpt->getNMchit(1) <<" "
1144  <<pSpt->getNMchit(2) <<" "
1145  <<pSpt->getNMchit(3) <<" "
1146  <<pSpt->getNMchit(4) <<" "
1147  <<pSpt->getXg(0) <<" "
1148  <<pSpt->getXg(1) <<" "
1149  <<pSpt->getXg(2) <<" "
1150  <<pSpt->getXl(0) <<" "
1151  <<pSpt->getXl(1) <<" "
1152  <<pSpt->getXl(2) <<" "
1153  <<a <<" "
1154  <<b <<" "
1155  <<endm;
1156  printf("pulseP =%f pulseN = %f\n",a,b);
1157  pSpt = mySsd->mLadders[i]->mWafers[j]->getPoint()->next(pSpt);
1158  }
1159  }
1160  }
1161  }
1162  }
1163  if (found==0) {LOG_INFO <<"PrintPointDetails() - Wafer not found !!!"<<endm; }
1164 }
1165 
1166 //_____________________________________________________________________________
1167 void StSsdPointMaker::PrintInfo()
1168 {
1169  if (Debug()) StMaker::PrintInfo();
1170 }
1171 //_____________________________________________________________________________
1172 void StSsdPointMaker::Read_Strip(St_ssdStripCalib *strip_calib)
1173 {
1174  ssdStripCalib_st *noise = strip_calib->GetTable();
1175  Int_t mSsdLayer = 7;
1176  LOG_INFO << "Read_Strip : printing few pedestal/noise values " << endm;
1177  Int_t idWaf = 0;
1178  Int_t iWaf = 0;
1179  Int_t iLad = 0;
1180  Int_t iZero = 0;
1181  Int_t nStrip = 0;
1182  Int_t iSide = 0;
1183  Int_t iOver = 0;
1184  Int_t iUnder = 0;
1185  Int_t iGood = 0;
1186  for (Int_t i = 0 ; i < strip_calib->GetNRows(); i++)
1187  {
1188  if (noise[i].id>0 && noise[i].id<=76818620) {
1189  nStrip = (int)(noise[i].id/100000.);
1190  idWaf = noise[i].id-10000*((int)(noise[i].id/10000.));
1191  iWaf = (int)((idWaf - mSsdLayer*1000)/100 - 1);
1192  iLad = (int)(idWaf - mSsdLayer*1000 - (iWaf+1)*100 - 1);
1193  iSide = (noise[i].id - nStrip*100000 - idWaf)/10000;
1194  if (iLad==11 && iWaf==8 && nStrip <10) {
1195  LOG_DEBUG<<"ReadStrip: iLad,idWaf,nStrip,iSide,pedestal,rms = "<<iLad
1196  <<" "<<idWaf
1197  <<" "<<nStrip
1198  <<" "<<iSide
1199  <<" "<<(float)(noise[i].pedestals)
1200  <<" "<<(float)(noise[i].rms)<<endm;
1201  iGood++;
1202  }
1203  }
1204  else {
1205  if (noise[i].id<0) iUnder++;
1206  else {
1207  if (noise[i].id==0) iZero++;
1208  else iOver++;
1209  }
1210  }
1211  }
1212  LOG_INFO<<"ReadStrip: Number of rows in the table : "<<strip_calib->GetNRows()<<endm;
1213  LOG_INFO<<"ReadStrip: Number of good id : "<<iGood<<endm;
1214  LOG_INFO<<"ReadStrip: Number of id = 0 : "<<iZero<<endm;
1215  if (iUnder>0){
1216  LOG_WARN <<"ReadStrip: Number of underf : "<<iUnder<<endm;}
1217  if (iOver>0){
1218  LOG_WARN <<"ReadStrip: Number of overf : "<<iOver<<endm;}
1219  }
1220 
1221 //_____________________________________________________________________________
1222 void StSsdPointMaker::Read_Strip(St_ssdNoise *strip)
1223 {
1224  ssdNoise_st *noise = strip->GetTable();
1225  LOG_INFO << " ReadStrip : printing few pedestal/noise values "<< endm;
1226  Int_t iLad = 0;
1227  for (Int_t i = 0 ; i < strip->GetNRows(); i++)
1228  {
1229  iLad = noise[i].id/16;
1230  Int_t idWaf = noise[i].id;
1231  if (idWaf==1)
1232  {
1233  for(Int_t nStrip=0;nStrip<10;nStrip++)
1234  {
1235  LOG_DEBUG<<"ReadStrip: iLad,idWaf,nStrip,rmsP,rmsN = "<<iLad+1
1236  <<" "<<idWaf+1
1237  <<" "<<nStrip
1238  <<" "<<(int)(noise[i].rmsp[nStrip])
1239  <<" "<<(int)(noise[i].rmsn[nStrip])<<endm;
1240  }
1241  }
1242  }
1243  LOG_INFO<<"ReadStrip: Number of rows in the table : "<<strip->GetNRows()<<endm;
1244  }
1245 //_____________________________________________________________________________
1246 
1247 void StSsdPointMaker::WriteMatchedStrips(StSsdBarrel *mySsd)
1248 {
1249  //Int_t LadderIsActive[20]={1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1};
1250  for (Int_t i=0;i<20;i++)
1251  if (mySsd->isActiveLadder(i)>0) {
1252  for (Int_t j=0; j<mySsd->mLadders[i]->getWaferPerLadder();j++) {
1253  LOG_DEBUG << " in ladder= "<<i << " wafer = "<<j<<endm;
1254 
1255  if (mySsd->mLadders[i]->mWafers[j]->getClusterP()->getSize()==0) {
1256  }
1257  else{
1258  LOG_DEBUG << " Size of the cluster P list = " << mySsd->mLadders[i]->mWafers[j]->getClusterP()->getSize()<<endm;
1259  StSsdCluster *pclusterP = mySsd->mLadders[i]->mWafers[j]->getClusterP()->first();
1260  while (pclusterP)
1261  {
1262  Int_t IdP = pclusterP->getNCluster();
1263  //LOG_DEBUG << " we are looking for clusterId= " << pclusterP->getNCluster()<< endm;
1264  if (mySsd->mLadders[i]->mWafers[j]->getPoint()->getSize()==0) {
1265  if (pclusterP!=mySsd->mLadders[i]->mWafers[j]->getClusterP()->last())
1266  pclusterP = mySsd->mLadders[i]->mWafers[j]->getClusterP()->next(pclusterP);
1267  }
1268  else {
1269  StSsdPoint *pSpt = mySsd->mLadders[i]->mWafers[j]->getPoint()->first();
1270  Int_t stop = 0;
1271  while ((pSpt)&&(stop==0)){
1272  if(IdP==pSpt->getIdClusterP())
1273  {
1274  StSsdStrip *pStripP = mySsd->mLadders[i]->mWafers[j]->getStripP()->first();
1275  StSsdStrip *pStripLastP = mySsd->mLadders[i]->mWafers[j]->getStripP()->last();
1276  while ((pStripP->getNStrip()!=pStripLastP->getNStrip())&&(pStripP->getNStrip()!=pclusterP->getFirstStrip())){
1277  pStripP = mySsd->mLadders[i]->mWafers[j]->getStripP()->next(pStripP);
1278  }
1279  Int_t stripInCluster = pclusterP->getClusterSize();
1280  Int_t k = 0;
1281  while(k<stripInCluster){
1282  LOG_DEBUG<<"PrintStripDetails() - "
1283  << k << " "
1284  <<pStripP->getNStrip()<<" "
1285  <<pStripP->getDigitSig()<<" "
1286  <<pStripP->getPedestal()<<" "
1287  <<pStripP->getSigma()<<" "
1288  <<endm;
1289  StripsIn[0] = 0;
1290  StripsIn[1] = i+1;
1291  StripsIn[2] = j+1;
1292  StripsIn[3] = pStripP->getNStrip();
1293  StripsIn[4] = pStripP->getPedestal();
1294  StripsIn[5] = pStripP->getDigitSig();
1295  StripsIn[6] = pStripP->getSigma();
1296  StripsIn[7] = (float)(pStripP->getDigitSig()/pStripP->getSigma());
1297  rHitNtuple->Fill(StripsIn);
1298  k++;
1299  stop=1;// 1 cluster matched to 1 hit
1300  pStripP = mySsd->mLadders[i]->mWafers[j]->getStripP()->next(pStripP);
1301  }
1302  }
1303  pSpt = mySsd->mLadders[i]->mWafers[j]->getPoint()->next(pSpt);
1304  }
1305  }
1306  pclusterP = mySsd->mLadders[i]->mWafers[j]->getClusterP()->next(pclusterP);
1307  }
1308  }
1309  if (mySsd->mLadders[i]->mWafers[j]->getClusterP()->getSize()==0) {
1310  }
1311  else{
1312  LOG_DEBUG << " Size of the cluster N list = " << mySsd->mLadders[i]->mWafers[j]->getClusterN()->getSize()<<endm;
1313  StSsdCluster *pclusterN = mySsd->mLadders[i]->mWafers[j]->getClusterN()->first();
1314  while (pclusterN)
1315  {
1316  Int_t IdN = pclusterN->getNCluster();
1317  if (mySsd->mLadders[i]->mWafers[j]->getPoint()->getSize()==0) {
1318  if (pclusterN!=mySsd->mLadders[i]->mWafers[j]->getClusterN()->last())
1319  pclusterN = mySsd->mLadders[i]->mWafers[j]->getClusterN()->next(pclusterN);
1320  }
1321  else {
1322  StSsdPoint *pSpt = mySsd->mLadders[i]->mWafers[j]->getPoint()->first();
1323  Int_t stop = 0;
1324  while ((pSpt)&&(stop==0)){
1325  if(IdN==pSpt->getIdClusterN())
1326  {
1327  StSsdStrip *pStripN = mySsd->mLadders[i]->mWafers[j]->getStripN()->first();
1328  StSsdStrip *pStripLastN = mySsd->mLadders[i]->mWafers[j]->getStripN()->last();
1329  while ((pStripN->getNStrip()!=pStripLastN->getNStrip())&&(pStripN->getNStrip()!=pclusterN->getFirstStrip())){
1330  pStripN = mySsd->mLadders[i]->mWafers[j]->getStripN()->next(pStripN);
1331  }
1332  Int_t stripInClusterN = pclusterN->getClusterSize();
1333  Int_t kk = 0;
1334  while(kk<stripInClusterN){
1335  LOG_DEBUG<<"PrintStripDetails() - "
1336  << kk << " "
1337  <<pStripN->getNStrip()<<" "
1338  <<pStripN->getDigitSig()<<" "
1339  <<pStripN->getPedestal()<<" "
1340  <<pStripN->getSigma()<<" "
1341  <<endm;
1342  StripsIn[0] = 1;
1343  StripsIn[1] = i+1;
1344  StripsIn[2] = j+1;
1345  StripsIn[3] = pStripN->getNStrip();
1346  StripsIn[4] = pStripN->getPedestal();
1347  StripsIn[5] = pStripN->getDigitSig();
1348  StripsIn[6] = pStripN->getSigma();
1349  StripsIn[7] = (float)(pStripN->getDigitSig()/pStripN->getSigma());
1350  rHitNtuple->Fill(StripsIn);
1351  kk++;
1352  stop=1;// 1 cluster matched to 1 hit
1353  pStripN = mySsd->mLadders[i]->mWafers[j]->getStripN()->next(pStripN);
1354  }
1355  }
1356  pSpt = mySsd->mLadders[i]->mWafers[j]->getPoint()->next(pSpt);
1357  }
1358  }
1359  pclusterN = mySsd->mLadders[i]->mWafers[j]->getClusterN()->next(pclusterN);
1360  }
1361  }
1362  }
1363  }
1364 }
1365 //_____________________________________________________________________________
1366 void StSsdPointMaker::WriteMatchedClusters(StSsdBarrel *mySsd)
1367 {
1368  //Int_t LadderIsActive[20]={1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1};
1369  Int_t clusP = 0;
1370  Int_t clusN = 0;
1371  for (Int_t i=0;i<20;i++)
1372  if (mySsd->isActiveLadder(i)>0) {
1373  for (Int_t j=0; j<mySsd->mLadders[i]->getWaferPerLadder();j++) {
1374  LOG_DEBUG << " in ladder= "<<i << " wafer = "<<j<<endm;
1375  if (mySsd->mLadders[i]->mWafers[j]->getClusterP()->getSize()==0) {
1376  }
1377  else{
1378  LOG_DEBUG << " Size of the cluster P list = " << mySsd->mLadders[i]->mWafers[j]->getClusterP()->getSize()<<endm;
1379  StSsdCluster *pclusterP = mySsd->mLadders[i]->mWafers[j]->getClusterP()->first();
1380  while (pclusterP)
1381  {
1382  LOG_DEBUG << " clusterId= " << pclusterP->getNCluster()<< endm;
1383  Int_t IdP = pclusterP->getNCluster();
1384  LOG_DEBUG << " we are looking for clusterId= " << pclusterP->getNCluster()<< endm;
1385  if (mySsd->mLadders[i]->mWafers[j]->getPoint()->getSize()==0) {
1386  pclusterP = mySsd->mLadders[i]->mWafers[j]->getClusterP()->next(pclusterP);
1387  }
1388  else {
1389  StSsdPoint *pSpt = mySsd->mLadders[i]->mWafers[j]->getPoint()->first();
1390  while (pSpt){
1391  if(IdP==pSpt->getIdClusterP())
1392  {
1393  clusP++;
1394  LOG_DEBUG << "ok found the corresponding hit to this cluster id = "<<pSpt->getIdClusterP()<<endm;
1395  ClustupleIn[0] = 0;
1396  ClustupleIn[1] = i+1;
1397  ClustupleIn[2] = j+1;
1398  ClustupleIn[3] = pclusterP->getClusterSize();
1399  ClustupleIn[4] = ((pclusterP->getTotAdc()*pclusterP->getClusterSize())/pclusterP->getTotNoise());
1400  ClustupleIn[5] = pclusterP->getTotNoise()/pclusterP->getClusterSize();
1401  ClustupleIn[6] = pclusterP->getFirstStrip();
1402  ClustupleIn[7] = pclusterP->getTotAdc();
1403  ClustupleIn[8] = pclusterP->getFirstAdc();
1404  ClustupleIn[9] = pclusterP->getLastAdc();
1405  ClustupleIn[10]= pclusterP->getTotNoise();
1406  pHitNtuple->Fill(ClustupleIn);
1407  //break;
1408  }
1409  pSpt = mySsd->mLadders[i]->mWafers[j]->getPoint()->next(pSpt);
1410  }
1411  pclusterP = mySsd->mLadders[i]->mWafers[j]->getClusterP()->next(pclusterP);
1412  }
1413  }
1414  clusP = 0;
1415  }
1416  if (mySsd->mLadders[i]->mWafers[j]->getClusterN()->getSize()==0) {
1417  }
1418  else{
1419  LOG_DEBUG << " Size of the cluster N list = " << mySsd->mLadders[i]->mWafers[j]->getClusterN()->getSize()<<endm;
1420  StSsdCluster *pclusterN = mySsd->mLadders[i]->mWafers[j]->getClusterN()->first();
1421  while (pclusterN)
1422  {
1423  //LOG_DEBUG << " clusterId= " << pclusterN->getNCluster()<< endm;
1424  Int_t IdN = pclusterN->getNCluster();
1425  //LOG_DEBUG << " we are looking for clusterId= " << pclusterN->getNCluster()<< endm;
1426  if (mySsd->mLadders[i]->mWafers[j]->getPoint()->getSize()==0) {
1427  pclusterN = mySsd->mLadders[i]->mWafers[j]->getClusterN()->next(pclusterN);
1428  }
1429  else {
1430  StSsdPoint *pSpt = mySsd->mLadders[i]->mWafers[j]->getPoint()->first();
1431  while (pSpt){
1432  if(IdN==pSpt->getIdClusterN())
1433  {
1434  clusN++;
1435  LOG_DEBUG << "ok found the corresponding hit to this cluster id = "<<pSpt->getIdClusterN()<<endm;
1436  ClustupleIn[0] = 1;
1437  ClustupleIn[1] = i+1;
1438  ClustupleIn[2] = j+1;
1439  ClustupleIn[3] = pclusterN->getClusterSize();
1440  ClustupleIn[4] = ((pclusterN->getTotAdc()*pclusterN->getClusterSize())/pclusterN->getTotNoise());
1441  ClustupleIn[5] = pclusterN->getTotNoise()/pclusterN->getClusterSize();
1442  ClustupleIn[6] = pclusterN->getFirstStrip();
1443  ClustupleIn[7] = pclusterN->getTotAdc();
1444  ClustupleIn[8] = pclusterN->getFirstAdc();
1445  ClustupleIn[9] = pclusterN->getLastAdc();
1446  ClustupleIn[10]= pclusterN->getTotNoise();
1447  pHitNtuple->Fill(ClustupleIn);
1448  //break;
1449  }
1450  pSpt = mySsd->mLadders[i]->mWafers[j]->getPoint()->next(pSpt);
1451  }
1452  pclusterN = mySsd->mLadders[i]->mWafers[j]->getClusterN()->next(pclusterN);
1453  }
1454  }
1455  clusN = 0;
1456  }
1457  }
1458  }
1459  LOG_DEBUG << "Number of p-side clusters in hits = " << clusP << endm ;
1460  LOG_DEBUG << "Number of n-side clusters in hits = " << clusN << endm ;
1461 }
1462 //_____________________________________________________________________________
1463 void StSsdPointMaker::EvaluateEfficiency(StSsdBarrel *mySsd)
1464 {
1465  //Int_t LadderIsActive[20]={1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1};
1466  Int_t clusP = 0;
1467  Int_t clusN = 0;
1468  for (Int_t i=0;i<20;i++)
1469  if (mySsd->isActiveLadder(i)>0) {
1470  for (Int_t j=0; j<mySsd->mLadders[i]->getWaferPerLadder();j++) {
1471  if (mySsd->mLadders[i]->mWafers[j]->getClusterP()->getSize()==0) {
1472  }
1473  else{
1474  StSsdCluster *pclusterP = mySsd->mLadders[i]->mWafers[j]->getClusterP()->first();
1475  while (pclusterP)
1476  {
1477  Int_t IdP = pclusterP->getNCluster();
1478  if (mySsd->mLadders[i]->mWafers[j]->getPoint()->getSize()==0) {
1479  pclusterP = mySsd->mLadders[i]->mWafers[j]->getClusterP()->next(pclusterP);
1480  }
1481  else {
1482  StSsdPoint *pSpt = mySsd->mLadders[i]->mWafers[j]->getPoint()->first();
1483  while (pSpt){
1484  if(IdP==pSpt->getIdClusterP())
1485  {
1486  clusP++;
1487  //break;
1488  }
1489  pSpt = mySsd->mLadders[i]->mWafers[j]->getPoint()->next(pSpt);
1490  }
1491  pclusterP = mySsd->mLadders[i]->mWafers[j]->getClusterP()->next(pclusterP);
1492  }
1493  }
1494  ratioP[i][j]+=(float)clusP/mySsd->mLadders[i]->mWafers[j]->getClusterP()->getSize();
1495  clusP = 0;
1496  }
1497  if (mySsd->mLadders[i]->mWafers[j]->getClusterN()->getSize()==0) {
1498  }
1499  else{
1500  StSsdCluster *pclusterN = mySsd->mLadders[i]->mWafers[j]->getClusterN()->first();
1501  while (pclusterN)
1502  {
1503  Int_t IdN = pclusterN->getNCluster();
1504  if (mySsd->mLadders[i]->mWafers[j]->getPoint()->getSize()==0) {
1505  pclusterN = mySsd->mLadders[i]->mWafers[j]->getClusterN()->next(pclusterN);
1506  }
1507  else {
1508  StSsdPoint *pSpt = mySsd->mLadders[i]->mWafers[j]->getPoint()->first();
1509  while (pSpt){
1510  if(IdN==pSpt->getIdClusterN())
1511  {
1512  clusN++;
1513  //break;
1514  }
1515  pSpt = mySsd->mLadders[i]->mWafers[j]->getPoint()->next(pSpt);
1516  }
1517  pclusterN = mySsd->mLadders[i]->mWafers[j]->getClusterN()->next(pclusterN);
1518  }
1519  }
1520  ratioN[i][j]+=(float)clusN/mySsd->mLadders[i]->mWafers[j]->getClusterN()->getSize();
1521  clusN = 0;
1522  }
1523  }
1524  }
1525 }
1526 //_____________________________________________________________________________
1527 void StSsdPointMaker::NormalizeEfficiency(){
1528  Double_t NClusP = 0;
1529  Double_t NClusN = 0;
1530  for(Int_t ii=0 ;ii< 20 ;ii++)
1531  {
1532  for(Int_t jj=0 ;jj<16;jj++)
1533  {
1534  NClusP = ratioP[ii][jj];
1535  //LOG_INFO<<Form("NEvent=%d side P ladder=%d wafer=%d bin content=%f\n",NEvent,ii,jj,NClusP)<<endm;
1536  if(NEvent!=0){
1537  NClusP = NClusP/NEvent;
1538  MatchedClusterP->SetBinContent(ii+1,jj+1,NClusP);}
1539  NClusN = ratioN[ii][jj];
1540  if(NEvent!=0){
1541  //LOG_INFO<<Form("NEvent=%d side N ladder=%d wafer=%d bin content=%f\n",NEvent,ii,jj,NClusN)<<endm;
1542  NClusN = NClusN/NEvent;
1543  MatchedClusterN->SetBinContent(ii+1,jj+1,NClusN);}
1544  }
1545  }
1546 }
1547 //_____________________________________________________________________________
1548 void StSsdPointMaker::FillCalibTable(){
1549  mGain = (St_ssdGainCalibWafer*)GetDataBase("Calibrations/ssd/ssdGainCalibWafer");
1550  if(mGain){
1551  ssdGainCalibWafer_st *g = mGain->GetTable() ;
1552  Int_t size = mGain->GetNRows();
1553  LOG_INFO<<Form("Size of gain table = %d",(int)mGain->GetNRows())<<endm;
1554  for(Int_t i=0; i<size;i++){
1555  LOG_DEBUG<<Form(" Print entry %d : ladder=%d gain =%lf wafer=%d",i,g[i].nLadder,g[i].nGain,g[i].nWafer)<<endm;
1556  CalibArray[i] = g[i].nGain;
1557  }
1558  }
1559  else {
1560  LOG_WARN << "InitRun : No access to Gain Calib - will use the default gain" << endm;
1561  LOG_WARN << "We will use the default table" <<endm;
1562  for(Int_t i=0; i<320;i++){
1563  CalibArray[i] = 1;
1564  }
1565  }
1566 }
1567 //_____________________________________________________________________________
1568 void StSsdPointMaker::FillDefaultCalibTable(){
1569  LOG_INFO << " The calibration gain will not be used." << endm;
1570  for(Int_t i=0; i<320;i++){
1571  CalibArray[i] = 1;
1572  //LOG_INFO << Form("wafer=%d gain=%f",i,CalibArray[i])<<endm;
1573  }
1574 }
1575 //_____________________________________________________________________________
1576 void StSsdPointMaker::FillWaferTable(){
1577  mWafConfig = (St_ssdWaferConfiguration*) GetDataBase("Calibrations/ssd/ssdWaferConfiguration");
1578  if(mWafConfig){
1579  ssdWaferConfiguration_st *g = mWafConfig->GetTable() ;
1580  Int_t size = mWafConfig->GetNRows();
1581  for(Int_t i=0; i<size;i++){
1582  LOG_DEBUG<<Form(" Print entry=%d : ladder=%d wafer=%d status=%d",i,g[i].nLadder,g[i].nWafer,g[i].nStatus)<<endm;
1583  WafStatus[g[i].nLadder][g[i].nWafer] = g[i].nStatus;
1584  }
1585  }
1586  else {
1587  LOG_WARN << "InitRun : No access to Wafer Config - will use the default wafer config" << endm;
1588  LOG_WARN << "We will use the default table" <<endm;
1589  for(Int_t i=0; i<20;i++){
1590  for(Int_t j=0; j<16;j++){
1591  WafStatus[i][j] = 1;
1592  }
1593  }
1594  }
1595 }
1596 //_____________________________________________________________________________
1597 void StSsdPointMaker::FillDefaultWaferTable(){
1598  LOG_INFO << " The wafer configuration table will not be used." << endm;
1599  for(Int_t i=0; i<20;i++){
1600  for(Int_t j=0; j<16;j++){
1601  WafStatus[i][j] = 1;
1602  LOG_DEBUG << Form("wafer=%d gain=%f",i,CalibArray[i])<<endm;
1603  }
1604  }
1605 }
1606 //_____________________________________________________________________________
1607 
1608 Int_t StSsdPointMaker::ReadNoiseTable(StSsdBarrel *mySsd,Int_t year){
1609  Int_t noiseTableSize = 0;
1610 
1611  //ssdStripCalib is used for year <2007 and for the simulation
1612  if((year<7)||(mode==1))
1613  {
1614  if (!m_noise2)
1615  {
1616  LOG_WARN << "Make : No pedestal and noise values (ssdStripCalib table missing), will use default values" <<endm;
1617  noiseTableSize = mySsd->readNoiseDefault(mDynamicControl);
1618  }
1619  else
1620  {
1621  if(Debug()) {Read_Strip(m_noise2);}
1622  noiseTableSize = mySsd->readNoiseFromTable(m_noise2,mDynamicControl);
1623  }
1624  }
1625  else if (year>=7){
1626  if(!m_noise3)
1627  {
1628  LOG_WARN << "Make : No pedestal and noise values (ssdNoise table missing), will use default values" << endm;
1629  noiseTableSize = mySsd->readNoiseDefault(mDynamicControl);
1630  }
1631  else
1632  if(m_noise3)
1633  {
1634  if (Debug()){Read_Strip(m_noise3);}
1635  noiseTableSize = mySsd->readNoiseFromTable(m_noise3,mDynamicControl);
1636  }
1637  }
1638  return noiseTableSize;
1639 }
1640 //____________________________________________________________________________
1642  LOG_INFO << Form("Finish()") << endm;
1643  return kStOK;
1644 }
1645 
StSsdPointList * getPoint()
Returns the point list attached to this wafer.
Definition: StSsdWafer.hh:106
TH2S * ClustMapN
Map of number of clusters on the p-side ladders.
StSsdStripList * getStripP()
Returns the P-side strip list attached to this wafer.
Definition: StSsdWafer.hh:103
TH1S * orthoproj
(1p-1n) packages control matching.
TH1F * noisDisN
p-side distribution of cluster total charge.
virtual void AddData(TDataSet *data, const char *dir=".data")
User methods.
Definition: StMaker.cxx:332
TH1F * snRatioN
n-side distribution of noise.
StSsdStripList * getStripN()
Returns the N-side strip list attached to this wafer.
Definition: StSsdWafer.hh:104
TH1F * stpClusP
p-side distribution of signal to noise ratio.
TH1F * totChrgN
n-side distribution of strips per cluster.
TH1S * kind
orthonormal projection and perfect matching deviation.
TH2S * ClustMapP
p-side clusters entries vs n-side clusters entries
TH1F * totChrgP
p-side distribution of strips per cluster.
TH1F * snRatioP
p-side distribution of noise.
Int_t readNoiseFromTable(St_sdm_calib_db *spa_noise, StSsdDynamicControl *dynamicControl)
Definition: StSsdBarrel.cc:345
Definition: Stypes.h:42
Definition: Stypes.h:40
TH2S * ClusNvsClusP
n-side distribution of cluster total charge.
virtual Int_t Finish()
virtual Int_t Make()
StSsdPackageList * getPackage()
Returns the package list attached to this wafer.
Definition: StSsdWafer.hh:105
StSsdClusterList * getClusterN()
Returns the N-side cluster list attached to this wafer.
Definition: StSsdWafer.hh:99
StSsdClusterList * getClusterP()
Returns the P-side cluster list attached to this wafer.
Definition: StSsdWafer.hh:98
TH2S * matchisto_[20]
kind of hits –&gt;see StSsdWafer for definition
Definition: Stypes.h:44
Event data structure to hold all information from a Monte Carlo simulation. This class is the interfa...
Definition: StMcEvent.hh:169
Definition: Stypes.h:41
virtual TDataSet * Find(const char *path) const
Definition: TDataSet.cxx:362
TH1F * stpClusN
n-side distribution of signal to noise ratio.