StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
RawPixels.cxx
1 #include <stdio.h>
2 #include <assert.h>
3 #include <string.h>
4 
5 
6 #include <TObjArray.h>
7 #include <TClonesArray.h>
8 
9 #include <TH2.h>
10 #include <TF2.h>
11 #include <TFile.h>
12 
13 #include "StMuDSTMaker/EZTREE/EztEmcRawData.h"
14 
15 #include "StEEmcUtil/EEfeeRaw/EEname2Index.h"
16 #include "StEEmcUtil/EEfeeRaw/EEdims.h"
17 
18 #include "StEEmcUtil/database/EEmcDbItem.h"
19 #include "StEEmcUtil/database/StEEmcDb.h"
20 
21 #include "RawPixels.h"
22 
23 //-------------------------------------------
24 //-------------------------------------------
25 RawPixels::RawPixels(TObjArray*L,StEEmcDb*dbx) {
26 
27  // clear histo pointers and all other stuff
28  HList=L; assert(HList);
29  hPix=new TH1F * [EEindexMax];
30  memset(hPix,0,sizeof(void*)*EEindexMax);
31 
32  hSmd=new TH2F * [MaxSmdPlains*MaxSectors];
33  memset(hSmd,0,sizeof(void*)*MaxSmdPlains*MaxSectors);
34 
35  eeDb=dbx;
36  setLimits(0,0);
37 }
38 
39 
40 //-------------------------------------------
41 //-------------------------------------------
42 void RawPixels::initHisto() {
43  assert(eeDb);
44  // eeDb->print();
45  assert(c_x1>=0);
46  assert(c_x2>c_x1);
47 
48  TString mode="";
49  switch(c_convMode) {
50  case kRawAdc: mode="rawAdc"; break;
51  case kPedSub: mode="pedSub"; break;
52  case kPedAndGain: mode="pedAndGain"; break;
53  default:
54  assert(2==3);
55  }
56  // ................... histo for individual pixels ...
57  int i,k=0;
58  for(i=0;i<EEindexMax;i++) {
59  const EEmcDbItem *x=eeDb->getByIndex(i);
60  if(x==0) continue;
61  // initialize histos only for pixels acquired from DB
62  k++;
63  char tt1[100],tt2[200];
64  sprintf(tt1,"a%s",x->name);
65  sprintf(tt2,"ADC for %s, cr/chan=%3.3d/%3.3d, tube=%s; ADC (mode=%s)",x->name,
66  x->crate,x->chan,x->tube,mode.Data());
67  TH1F* h=new TH1F(tt1,tt2,c_x2-c_x1+1,c_x1-0.5,c_x2+0.5);
68  hPix[i]=h;
69  HList->Add(h);
70  // printf("k=%d -->%s\n",k,tt1);
71  }
72 
73  // ..... 2D histos for SMD
74  int sec1=eeDb->getFirstSecID();
75  int sec2=eeDb->getLastSecID();
76  int secID;
77  for(secID=sec1; secID<=sec2; secID++) {
78  int iuv;
79  for(iuv=0;iuv<MaxSmdPlains;iuv++) {
80  char tt1[100], tt2[100];
81  sprintf(tt1,"a%2.2d%c",secID,iuv+'U');
82  sprintf(tt2,"%2.2d%c-SMD strip vs. raw ADC",secID,iuv+'U');
83  TH2F* h2=new TH2F(tt1,tt2,200,0,4000,MaxSmdStrips,0.5,MaxSmdStrips+0.5);
84  int key1=secID-1+iuv*MaxSectors;
85  hSmd[key1]=h2;
86  HList->Add(h2);
87  //printf("init 2D %d '%s'\n",key1,tt1);
88  }
89  }
90 
91 
92  { // additional info
93  char tt2[100];
94  sprintf(tt2," stat info ");
95 
96  hInfo=new TH1F("info",tt2,20,-0.5,19.5);
97  HList->Add(hInfo);
98  }
99 
100  printf("RawPixels: Initialized histos for %d pixels\n",k);
101 }
102 
103 //-------------------------------------------
104 //-------------------------------------------
105 void RawPixels::sort(EztEmcRawData *eRaw){
106 
107  assert(eeDb);
108 
109  hInfo->Fill(0); // tot eve
110  if(eRaw==0) return;
111 
112  int icr;
113  for(icr=0;icr<eRaw->getNBlocks();icr++) {
114  if(eRaw->isCrateVoid(icr)) continue;
115  const UShort_t* data=eRaw->data(icr);
116  int crateID=eRaw->getCrateID(icr);
117  hInfo->Fill(crateID); // this crate per eve
118  int nd=eRaw->sizeData(icr);
119  // printf("ic =%d cr=%d \n",icr, crateID);
120 
121  int chan;
122  for(chan=0;chan<nd;chan++) {
123  const EEmcDbItem *x=eeDb->getByCrate(crateID,chan);
124  if(x==0) continue;
125  // process only pixels acquired from DB
126 
127  float value=data[chan]; // raw ADC
128  switch(c_convMode) {
129  case kRawAdc: break;
130  case kPedSub:
131  value-=x->ped; break;
132  case kPedAndGain:
133  value-=x->ped;
134  value/=x->gain;
135  //printf("kPedAndGain %s %s %f %f\n",x->name,x->tube,x->gain,value);
136  break;
137  default:
138  assert(2==3);
139  }
140 
141  hPix[x->key]->Fill(value);
142  // printf("cr=%2d ch=%3d val=%.1f ped=%4.1f '%s' \n",crateID,chan,value,x->ped,x->name);
143 
144  if(x->isSMD()) {
145  int key=(x->plane-'U')*MaxSectors + x->sec-1;
146  //x->print();
147  //printf("key=%d\n",key);
148 #if 0
149  if(key<=0 || key>=MaxSmdPlains*MaxSectors) {
150  printf("cr=%2d ch=%3d val=%.1f ped=%4.1f '%s' \n",crateID,chan,value,x->ped,x->name);
151  x->print();
152  printf("key=%d\n",key);
153  }
154 #endif
155 
156  assert(key>=0 && key<MaxSmdPlains*MaxSectors);
157  assert(hSmd[key]);
158  hSmd[key]->Fill(value,x->strip);
159  }
160 
161 
162  } // channels
163  } // blocks
164 
165 }
const EEmcDbItem * getByIndex(int ikey) const
returns full DB info for one pixel
Definition: StEEmcDb.cxx:584
char name[StEEmcNameLen]
ASCII name of the channel, see Readme.
Definition: EEmcDbItem.h:20
char tube[StEEmcNameLen]
name of PMT or MAPMT pixel
Definition: EEmcDbItem.h:21
int chan
hardware channel
Definition: EEmcDbItem.h:28