StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
St_emc_Maker.cxx
1 // $Id: St_emc_Maker.cxx,v 1.15 2017/04/26 20:22:16 perev Exp $
2 // $Log: St_emc_Maker.cxx,v $
3 // Revision 1.15 2017/04/26 20:22:16 perev
4 // Hide m_DataSet
5 //
6 // Revision 1.14 2007/04/28 17:56:05 perev
7 // Redundant StChain.h removed
8 //
9 // Revision 1.13 2003/09/02 17:59:29 perev
10 // gcc 3.2 updates + WarnOff
11 //
12 // Revision 1.12 2003/04/30 20:39:16 perev
13 // Warnings cleanup. Modified lines marked VP
14 //
15 // Revision 1.11 2000/01/29 00:04:58 akio
16 // temprary fix for endcap. need more work, but no more junk messages and crash
17 //
18 // Revision 1.10 1999/07/16 18:24:27 pavlinov
19 // Little correction for StEclMake
20 //
21 // Revision 1.9 1999/07/15 13:58:01 perev
22 // cleanup
23 //
24 // Revision 1.8 1999/07/02 03:01:56 pavlinov
25 // Little corrections for Linux
26 //
27 // Revision 1.7 1999/07/01 16:17:57 pavlinov
28 // class StEmcGeom was created and maker was remade for new maker scheme
29 //
30 // Revision 1.6 1999/03/03 04:12:15 fisyak
31 // replace kStErr to kStWarn
32 //
33 // Revision 1.5 1999/02/26 17:28:50 kathy
34 // fix histograms
35 //
36 // Revision 1.4 1999/02/16 18:15:44 fisyak
37 // Check in the latest updates to fix them
38 //
39 // Revision 1.3 1999/02/12 19:16:30 akio
40 // *** empty log message ***
41 //
42 // Revision 1.2 1998/12/15 22:39:49 akio
43 // Add emc_hit object and adc_to_energy in here.
44 //
46 // //
47 // St_emc_Maker class for <FONT COLOR="RED"> EMC</FONT> dataset //
48 // //
50 #include <Stiostream.h>
51 #include <math.h>
52 #include "St_emc_Maker.h"
53 #include "St_DataSetIter.h"
54 ClassImp(St_emc_Maker)
55 
56 const TString detname[] = {
57  "bemc", "bprs", "bsmde", "bsmdp",
58  "eemc", "eprs", "esmde", "esmdp"
59 };
60 //_____________________________________________________________________________
61 St_emc_Maker::St_emc_Maker(const char *name, const char *title):StMaker(name,title){
62  drawinit=kFALSE;
63 }
64 //_____________________________________________________________________________
65 St_emc_Maker::~St_emc_Maker(){ /* Nobody */ }
66 //_____________________________________________________________________________
67 Int_t St_emc_Maker::Init(){
68  Int_t i;
69  //Making QA histgrams
70  const char *title_h[] = {
71  "Barrel EMC hits","Barrel PRS hits","Barrel SMD-Eta hits","Barrel SMD-Phi hits",
72  "Endcap EMC hits","Endcap PRS hits","Endcap SMD-u hits", "Endcap SMD-v hits"};
73  const char *title_e[] = {
74  "Barrel EMC energy","Barrel PRS energy","Barrel SMD-Eta energy","Barrel SMD-Phi energy",
75  "Endcap EMC energy","Endcap PRS energy","Endcap SMD-u energy", "Endcap SMD-v energy"};
76  const Int_t nx[] = {40,40,300,20,12,12,12,12};
77  const Float_t xl[] = {-1.0,-1.0,-1.0,-1.0, 0.5 , 0.5, 0.5, 0.5};
78  const Float_t xu[] = { 1.0, 1.0, 1.0, 1.0, 12.5,12.5,12.5,12.5};
79  const Int_t ny[] = {120, 120, 60, 900, 60, 60, 60, 60};
80  m_nhit = new TH2F("EmcNHitsVsDet" ,"Number of hit(log) .vs. Detector #",100,0.0,4.5,8,0.5,8.5);
81  m_etot = new TH2F("EmcEtotVsDet" ,"Total energy(log) .vs. Detector #",100,-4.0,4.5,8,0.5,8.5);
82  for (i=0; i<MAXDET; i++){
83  TString name_h = detname[i] + "Hits";
84  TString name_e = detname[i] + "Energy";
85  Float_t rpi = M_PI + 0.00001;
86  m_hits[i] = new TH2F(name_h,title_h[i],nx[i],xl[i],xu[i],ny[i],-rpi, rpi);
87  m_energy[i] = new TH2F(name_e,title_e[i],nx[i],xl[i],xu[i],ny[i],-rpi, rpi);
88  }
89 
90  return StMaker::Init();
91 }
92 //_____________________________________________________________________________
94 
95  mEmcCalib = GetInputDB("calib");
96  if(!mEmcCalib){
97  cout << "Warning in St_emc_Maker: Database not found" << endl;
98  return kStWarn;;
99  }else{
100  assert(mEmcCalib);
101  }
102 
103  if (!GetData()->GetList()){ //if DataSet is empty, create object and fill it
104  //Making Hits
105  St_DataSet *emcRaw = GetDataSet("emc_raw/.data");
106  if(!emcRaw) return kStWarn;
107  St_DataSetIter itr(emcRaw);
108 
109  St_emc_hits *adc = 0;
110  St_emc_hits dummy; TString tit = dummy.GetTitle();
111  while ((adc = (St_emc_hits *)itr())) {
112  if(adc->GetTitle() == tit){ // Get only emc_hits
113  TString name = adc->GetName();
114  name.ReplaceAll("emc_hits_","");
115  for(int i=4; i<8; i++){if(detname[i]==name){goto SKIP;}} //only barrel for now
117  hit->setEmcCalib(mEmcCalib);
118  AddData(hit);
119  if(hit->fill(adc) != kStOK) return kStWarn;
120  }
121  SKIP:;
122  }
123  }
124 
125  if(m_mode){ // make a copy in STAF table 'emc_hits' for STAF module
126  St_DataSetIter itr(GetData());
127  StEmcHitCollection *hit = 0;
128  StEmcHitCollection dummy; TString tit = dummy.GetTitle();
129  while ((hit = (StEmcHitCollection *)itr())) {
130  if(hit->GetTitle() == tit){
131  TString name_hits = "emc_hits_" + TString(hit->GetName());
132  AddData(hit->copyToTable(name_hits));
133  }
134  }
135  }
136 
137  MakeHistograms(); // Fill QA histgrams
138 
139  return kStOK;
140 }
141 //_____________________________________________________________________________
142 void St_emc_Maker::MakeHistograms(){
143  //Filling QA Histograms
144  if (GetData()) {
145  Int_t det, i, n, id;
146  Float_t eta, phi, ene;
147  St_DataSetIter itr(GetData());
148  StEmcHitCollection *hit = 0;
149  StEmcHitCollection dummy;
150  TString tit = dummy.GetTitle();
151  while ((hit = (StEmcHitCollection *)itr())) {
152  if(hit->GetTitle() == tit){
153  n = hit->NHit();
154  if(n>0){
155  det = hit->Detector();
156  m_nhit->Fill(log10((Float_t)hit->NHit()),(Float_t)det);
157  m_etot->Fill(log10(hit->EnergySum()),(Float_t)det);
158  for(i = 0; i < n; i++){
159  ene = hit->HitEnergy(i);
160  if(ene > 0.0){
161  id = hit->HitId(i);
162  if(hit->getEtaPhi(id, eta, phi) == 0){
163  m_hits[det-1]->Fill(eta,phi);
164  m_energy[det-1]->Fill(eta,phi,ene);
165  }
166  else {
167  printf(" <W> St_emc_Maker::MakeHistograms() => det %i bad id %i id energy %f \n"
168  ,det,id,ene);
169  }
170  }
171  }
172  }
173  }
174  }
175  }
176 }
177 //_____________________________________________________________________________
178 //_____________________________________________________________________________
virtual Int_t Make()
virtual void AddData(TDataSet *data, const char *dir=".data")
User methods.
Definition: StMaker.cxx:332
Definition: Stypes.h:42
Definition: Stypes.h:40