StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StFgtJanGainMaker.cxx
1 
11 #include <string>
12 #include "StFgtJanGainMaker.h"
13 #include "StRoot/StEvent/StFgtCollection.h"
14 #include "StRoot/StEvent/StFgtStrip.h"
15 #include "StRoot/StEvent/StEvent.h"
16 #include "StRoot/StFgtUtil/geometry/StFgtGeom.h"
17 #include "StRoot/StFgtDbMaker/StFgtDbMaker.h"
18 
19 //========================================================
20 //========================================================
21 //========================================================
22 
23 StFgtJanGainMaker::StFgtJanGainMaker( const Char_t* dbMkrName ,const Char_t* name ) : StMaker( name ), mDbMkrName( dbMkrName ), mFgtDbMkr(0) {
24  HList=0;
25 };
26 
27 
28 //========================================================
29 //========================================================
30 //========================================================
31 Int_t StFgtJanGainMaker::Init(){
32  assert(HList);
33 
34  GetEvtHddr()->SetEventNumber(1);
35 
36  StEvent* eventPtr=0;
37  eventPtr= (StEvent*)GetInputDS("StEvent");
38  assert(eventPtr);
39 
40  mFgtCollectionPtr=eventPtr->fgtCollection();
41  assert(mFgtCollectionPtr);
42  mFgtDbMkr = static_cast< StFgtDbMaker* >( GetMaker( mDbMkrName.data() ));
43  assert(mFgtDbMkr );
44 
45  LOG_INFO << "Using date and time " << mFgtDbMkr->GetDateTime().GetDate() << ", "
46  << mFgtDbMkr->GetDateTime().GetTime() << endm;
47 
48  iEvt=-1;
49  initHistos();
50  return kStOk;
51 };
52 
53 
54 //________________________________________________
55 //________________________________________________
56 void
57 StFgtJanGainMaker::initHistos(){
58  // const float PI=TMath::Pi();
59 
60  //...... data histograms
61  memset(hA,0,sizeof(hA));
62  hA[5]=new TH1F("fgt1","Seen APVs per event; # APVs/event",150,-0.5,149.5);
63 
64 
65  // add histos to the list (if provided)
66  for(int i=0;i<mxHA;i++) {
67  if( hA[i]==0) continue;
68  HList->Add( hA[i]);
69  }
70  // HList->ls();
71 
72 }
73 
74 //========================================================
75 //========================================================
76 //========================================================
78 
79  enum {Ntimebin=7};
80  Int_t adcA[Ntimebin];
81  int nPulse=0;
82  iEvt++;
83  // cout << "iEvt = " << iEvt << endl;
84  assert( mFgtDbMkr );
85  StFgtDb *fgtTables = mFgtDbMkr->getDbTables();
86  assert(fgtTables);
87 
88  StEvent* eventPtr = (StEvent*)GetInputDS("StEvent");
89  assert(eventPtr);
90 
91  mFgtCollectionPtr = mFgtCollectionPtr=eventPtr->fgtCollection();
92  assert( mFgtCollectionPtr);
93 
94  int lastGeo=-999;
95  for( UInt_t discIdx=0; discIdx<mFgtCollectionPtr->getNumDiscs(); ++discIdx ){
96  StFgtStripCollection *stripCollectionPtr = mFgtCollectionPtr->getStripCollection( discIdx );
97  if( stripCollectionPtr ){
98  StSPtrVecFgtStrip& stripVec = stripCollectionPtr->getStripVec();
99  StSPtrVecFgtStripIterator stripIter;
100 
101  for( stripIter = stripVec.begin(); stripIter != stripVec.end(); ++stripIter ){
102  int rdo=0, arm=-1, apv=-1, chn=-1;
103  Short_t disk=-1, quad=-1, strip=-1; char layer='x';
104  double ordinate=-1,lowerSpan=-1,upperSpan=-1.;
105  (*stripIter)->getElecCoords( rdo, arm, apv, chn );
106  int stat=fgtTables->getStatusFromElecCoord(rdo,arm,apv,chn);
107  if(stat) continue; // drop bad strips
108  int geoId=fgtTables->getGeoIdFromElecCoord(rdo, arm, apv, chn);
109 
110  assert (geoId>=0);
111  StFgtGeom::decodeGeoId(geoId,disk,quad,layer,strip);
112  StFgtGeom::getPhysicalCoordinate(geoId,disk,quad,layer,ordinate,lowerSpan,upperSpan);
113 
114  // drop unstable or dead APVs - Jan's private list
115  if(rdo==2 && arm==1 && apv==5) continue;
116  if(rdo==2 && arm==1 ) continue;
117  if(rdo==2 && arm==2 && apv==18) continue;
118  if(rdo==1 && arm==3 && apv==9) continue;
119  if(rdo==2 && arm==3 ) continue;
120  if(rdo==1 && arm==4 && apv==19) continue;
121  if(rdo==1 && arm==4 && apv==21) continue;
122  if(rdo==2 && arm==4 ) continue;
123 
124  double ped=fgtTables->getPedestalFromElecCoord(rdo,arm,apv,chn);
125  double pedSig=fgtTables->getPedestalSigmaFromElecCoord(rdo,arm,apv,chn);
126 
127  memset(adcA,0,sizeof(adcA));
128  float minAdc=9999, maxAdc=-9999, sum=0;
129  int iMin=-1;
130  for(Int_t is=0;is<Ntimebin;is++){
131  adcA[is]=(*stripIter)->getAdc(is)-ped;
132  if(adcA[is]<minAdc) { minAdc=adcA[is]; iMin=is;}
133  if(adcA[is]>maxAdc) maxAdc=adcA[is];
134  sum+=adcA[is];
135  }
136 
137  //printf("mm %f %f\n",ped,maxAdc);
138  if(maxAdc<600) continue;
139  if(maxAdc>3000) continue;
140  if(iMin>1) continue; // require ped is in time bin 0 or 1
141 
142  // save histo with pulse-shape
143  nPulse++;
144  TH1F * hsp=new TH1F(Form("ps%d_%d",iEvt,nPulse), Form("ieve=%d rdo=%d arm=%d apv=%d ch=%d geoId=%d strip=%d%c%c%03d; time bin",iEvt,rdo,arm,apv,chn,geoId,disk+1,quad+'A',layer,strip),Ntimebin+2,0,Ntimebin+2);
145  for(Int_t is=0;is<Ntimebin;is++){
146  hsp->SetBinContent(is+1,adcA[is]-ped);
147  hsp->SetBinError(is+1,pedSig);
148  }
149  HList->Add(hsp);
150  // end of histo-save
151 
152  char star=' ';
153  if(abs(lastGeo-geoId)==1) star='*';
154  lastGeo=geoId;
155  printf("\nieve=%d rdo=%d arm=%d apv=%d ch=%d geoId=%d strip=%d%c%c%03d %csumADC-ped=%.1f\n adc-ped[0...6]=",iEvt,rdo,arm,apv,chn,geoId,disk+1,quad+'A',layer,strip,star,sum);
156  for(Int_t is=0;is<Ntimebin;is++) printf(" %.0f ",adcA[is]-minAdc);
157  printf("\n dbPd=%.1f minAdc=%.1f del=%.1f\n",ped,minAdc,ped-minAdc);
158  }
159  }
160  }
161 
162 
163  return kStOK;
164 };
165 
166 
168 
169  cout << "StFgtJanGainMaker::Finish()" << endl;
170  return StMaker::Finish();
171 };
172 
173 
174 
175 
176 ClassImp( StFgtJanGainMaker );
177 
178 /**************************************************************************
179  *
180  * $Log: StFgtJanGainMaker.cxx,v $
181  * Revision 1.1 2012/02/07 08:25:29 balewski
182  * *** empty log message ***
183  *
184  * Revision 1.4 2012/02/07 06:14:42 balewski
185  * *** empty log message ***
186  *
187  * Revision 1.3 2012/02/07 05:33:30 balewski
188  * *** empty log message ***
189  *
190  * Revision 1.2 2012/02/06 04:17:36 balewski
191  * added 2012 APV exclusions
192  *
193  * Revision 1.1 2012/02/04 22:03:40 balewski
194  * start
195  *
196  *
197  **************************************************************************/
Definition: Stypes.h:40
virtual Int_t Finish()
Definition: StMaker.cxx:776
Definition: Stypes.h:41