StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StEEmcMixerMaker.cxx
1 #include <stdlib.h>
2 #include <string.h>
3 #include <TFile.h>
4 
5 #include <Stiostream.h>
6 #include <StMessMgr.h>
7 #include <StEventTypes.h>
8 #include <StEvent.h>
9 
10 #include "StEEmcSimulatorMaker/StEEmcFastMaker.h"
11 #include "StEEmcSimulatorMaker/StEEmcSlowMaker.h"
12 #include "StEEmcMixerMaker.h"
13 #include "StEEmcUtil/database/EEmcDbItem.h"
14 #include "StEEmcUtil/database/StEEmcDb.h"
15 #include "StEEmcUtil/database/cstructs/eemcConstDB.hh" // def of status bits
16 #include "StEEmcUtil/EEmcGeom/EEmcGeomDefs.h"
17 
18 #include "StDAQMaker/StDAQReader.h"
19 #include "StEmcRawMaker/StEemcRaw.h"
20 
21 #include "StEmcRawMaker/StEEmcPrint.h"
22 
23 
24 ClassImp(StEEmcMixerMaker)
25 
26 //-------------------------------------------------------------------
27 //-------------------------------------------------------------------
28 StEEmcMixerMaker::StEEmcMixerMaker(const char *name):StMaker(name){
29  panicOff=false; // once activated disables Endcap embedding
30  mEEDb = 0;
31 }
32 
33 //-------------------------------------------------------------------
34 //-------------------------------------------------------------------
35 StEEmcMixerMaker::~StEEmcMixerMaker() {
36 
37 }
38 
39 //-------------------------------------------------------------------
40 //-------------------------------------------------------------------
42  return StMaker::Finish();
43 }
44 
45 //-------------------------------------------------------------------
46 //-------------------------------------------------------------------
47 Int_t
48 StEEmcMixerMaker::Init(){
49  mEEDb=(StEEmcDb*)GetDataSet("StEEmcDb");
50  if( mEEDb==0){
51  panicOff=true;
52  LOG_FATAL<< "::Init()\n\n Fatal Error - Eemc_DbMaker is not in the chain,\n panicOff="<<panicOff<<endm;
53  return kStErr;
54  }
55  assert(!strcmp(mEEDb->ClassName(),"StEEmcDb"));
56  return StMaker::Init();
57 }
58 
59 //-------------------------------------------------------------------
60 //-------------------------------------------------------------------
61 Int_t
63 /* Embedding is performed in the level of StEmcCollection.
64  Data input is daq file and MC is from EEmc fast Simulator StEEmcFastMaker.
65 */
66 
67  LOG_INFO <<"::Make() start...."<<endm;
68  StEvent* mEvent = (StEvent*)GetInputDS("StEvent");
69  if(mEvent==0) panicOff=true;
70  StEmcCollection* ecolA =(StEmcCollection*)mEvent->emcCollection();
71  if(ecolA==0) {
72  LOG_WARN<<"::Make(), raw2pixels() no emc collection, skip"<<endm;
73  return kStErr;
74  }
75 
76  if(panicOff) {
77  LOG_FATAL<< "::Make()\n\n Fatal Error was encounter earlier, Endcap embedding disabled panicOff="<<panicOff<<endm;
78  return kStErr;
79  }
80 
81 /*
82  Fill data from StEEMCReader to StEmcRawData first, then from raw blocks
83  to StEmcCollection of StEvent using the mehtods copyRawData() and
84  raw2pixels() of StEemcRaw. They are changed from private to pubilc in
85  StEemcRaw for being used legally by StEEmcMixerMaker. Also, raw2pixels
86  is modified to accept three parameters: raw2pixels(bool mixer,
87  StEmcCollection*, StEvent*), allowing input of embedded StEmcCollection.
88  If mixer=false (default) , raw2pixels builds a StEmcCollection of a StEvent
89  from DAQ as event/ecolA is built here. If mix=true, the second input
90  StEEmcCollection* is used to build an embedded StEvent.
91 
92 */
93 
94  StEEmcPrint eemcPrint;
95  eemcPrint.setMode(15);//bits= 1:Tower, 2:Pre, 4:SmdU, 8:SmdV, 15:all
96 
97  LOG_DEBUG<<"::Make() -------------- print data: Ecoll-A ---- real backg eve ---------"<<endm;
98  if(Debug()) eemcPrint.print(ecolA);
99 
100 
101  /* If the second source is the EEMC simulator,
102  it owns the StEmcCollection from simulator.
103  */
104  StEmcCollection *ecolB = 0;
105  StEEmcFastMaker *sim = (StEEmcFastMaker*)GetMakerInheritsFrom("StEEmcFastMaker");
106  // one can use fast simu as the source even if slow simu is used, since slow simu is set up in the overwrite mode.
107  if(!sim) {
108  LOG_WARN<<"::Make() No fast EEmcSimulator found, nothing to embed"<<endm; return kStWarn; }
109  ecolB = sim->GetLocalEmcCollection();
110 
111  if(!ecolB) {
112  LOG_WARN<<"::Make() No second EmcCollection to embed"<<endm; return kStWarn; }
113 
114 
115  LOG_DEBUG <<"::Make() -------------- print data: Ecoll-B ----- M-C physics probe eve -----------"<<endm;
116  if(Debug()) eemcPrint.print(ecolB);
117 
118 
119  LOG_DEBUG<<GetName() <<"::Make() -------------- print data: Ecoll-A+B ----- before mrging -----"<<endm;
120  if(Debug()) eemcPrint.printChange(ecolA,ecolB,Form("before merging"));
121 
122  if(!mergeADCs(ecolA,ecolB)) return kStErr;
123 
124  mMixerEmcCollection = ecolA;
125 
126  LOG_DEBUG <<"::Make() -------------- print Ecoll-A after mixing ----------------"<<endm;
127  if(Debug()) eemcPrint.print(ecolA);
128 
129  LOG_DEBUG <<"::Make() -------------- print data: Ecoll-A+B ----- after mrging -----"<<endm;
130  if(Debug()) eemcPrint.printChange(ecolA,ecolB,Form("after merging"));
131  return kStOK;
132 }
133 
134 
135 //--------------------------------------------------------------------
136 //--------------------------------------------------------------------
137 bool
138 StEEmcMixerMaker::mergeADCs(StEmcCollection*emccolA,StEmcCollection*emccolB){
143  LOG_DEBUG <<"::Make() merging ADCs ..."<<endm;
144 
145 /*
146  For EEMC, module = sector(1-12).
147  Tower (detector-name eemc): sub(1-5), eta(1-12)
148  Pre/Post (detector-name eprs): sub(1-15), eta(1-12)
149  Pre1 sub(1-5), Pre2: sub(6-10), Post: sub(11-15)
150  Smd U & V (detector-name: esmdu and esmdv): sub=1, eta=strip(1-288)
151 */
152 
153  // output of merging: collA=collA+collB
154  for(int det = kEndcapEmcTowerId; det<= kEndcapSmdVStripId; det++){
155  // Loop over all four sub detectors
156  StDetectorId id = StDetectorId(det);
157  StEmcDetector* detectorA=emccolA->detector(id);
158  StEmcDetector* detectorB=emccolB->detector(id);
159  if(!detectorA)
160  LOG_WARN<<"detectorA not loaded"<<endm;
161  if(!detectorB)
162  LOG_WARN<<"detectorB not loaded"<<endm;
163 
164  if(!detectorA || !detectorB) continue;// nothing to mix for such layer
165  for(int secID=1; secID<=kEEmcNumSectors; secID++){
166  // if(secID!=6) continue;//tmp
167  StEmcModule* sectorA = detectorA->module(secID);
168  StEmcModule* sectorB = detectorB->module(secID);
169  StSPtrVecEmcRawHit& rawHitA=sectorA->hits();
170  StSPtrVecEmcRawHit& rawHitB=sectorB->hits();
171 
172  // clone pointers to hits collection B
173  vector<StEmcRawHit*> myHitB;
174  for(UInt_t k2=0;k2<rawHitB.size();k2++)
175  myHitB.push_back(rawHitB[k2]);
176  vector<StEmcRawHit*>::iterator hitB;
177 
178  // printf("\nmixIn idet=%d sect=%d Nhit A=%d B=%d\n",det,secID,rawHitA.size(),myHitB.size());
179 
180  /*....................
181  1) merge hits for pixels present in both collections */
182 
183  for(UInt_t k1=0;k1<rawHitA.size();k1++) {
184  /* do not bother with differences between BEMC & EEMC
185  numbering scheme, pretend it is barrel */
186  uint Bmod=rawHitA[k1]->module();
187 
188  if((int)Bmod==!secID) {
189  LOG_FATAL << "::Make()\n\n Fatal Error "<<Bmod<<"= Bmod==!secID ="<<secID<<" StEvent internal consistency failed - corrupted event file, abort"<<endm;
190  panicOff=true;
191  return false;
192 
193  }
194 
195  // StEvent internal consistency check
196  uint Beta=rawHitA[k1]->eta();
197  int Bsub=rawHitA[k1]->sub();
198 
199  // erase old energy for every hit, use ADC2E-maker to get energy back
200  rawHitA[k1]->setEnergy(-654.3210); // just in case, make more non-physical
201 
202  for( hitB=myHitB.begin() ; hitB< myHitB.end(); hitB++) {
203  if( (*hitB)->module()!=Bmod) continue;
204  if( (*hitB)->eta() !=Beta) continue;
205  if( (*hitB)->sub() !=Bsub) continue;
206  // ......... found match to the same element
207  int adc=(*hitB)->adc();
208 
209  myHitB.erase(hitB); // to shrink the list for next search
210  /* Note, *hitB is now a stray pointer, do NOT use it
211  for anything. The loop must terminate with 'break',
212  since hitB++ makes no sense any more.
213 
214  It is assumed there is only one hit per pixel in
215  collectionB and it was just found.
216  */
217  hitB=(vector<StEmcRawHit*>::iterator)NULL; // siher ist siher
218 
219  const EEmcDbItem *x=mEEDb->StBarrelIndex2Item(det,Bmod,Beta,Bsub);
220  if(!x) break;// DB info not avaliable, drop this secondary hit
221 
222  /* it is assumed GEANT energy deposit is stored
223  in secondary collection
224  */
225  int adcAdd=adc;
226  if(adcAdd<=0) break;
227  int adcSum=rawHitA[k1]->adc() +adcAdd;
228  //printf("D: %d %g %d %d %s\n",rawHitA[k1]->adc(),deB,adcAdd,adcSum,x->name);
229  rawHitA[k1]->setAdc(adcSum);
230  break;
231  } // loop overB
232 
233  } // loop overA
234  // printf("mixMid idet=%d sect=%d Nhit A=%d B=%d\n",det,secID,rawHitA.size(),myHitB.size());
235 
236 
237  /* 2) verify no hits remain in collection B.
238  some hits may be left for real events:
239  * in 2004 or earlier (due to the code implementation)
240  * in 2005+ if a crate is masked out due to any form of corruption or being OFF
241  */
242 
243  if(myHitB.size())
244  LOG_WARN<<Form("%s::Make() merging: %d hits in collB for sect=%d are dropped\n since those channels are not avaliable in collA, \n a create is probably masked out\n",GetName(), myHitB.size(),secID)<<endm;
245  // printf("mixEnd idet=%d sect=%d Nhit A=%d B=%d\n",det,secID,rawHitA.size(),myHitB.size());
246 
247  } // loop over sector
248  } // loop over detectors
249 
250  return true; // all finished w/o problems
251 }
252 
253 
254 
256 //
257 // $Id: StEEmcMixerMaker.cxx,v 1.10 2010/01/28 14:16:46 jwebb Exp $
258 // $Log: StEEmcMixerMaker.cxx,v $
259 // Revision 1.10 2010/01/28 14:16:46 jwebb
260 // (1) Use "Form" to get around deprecated conversion from string constant to
261 // char *
262 // (2) eemcPrint called on debug only.
263 //
264 // Revision 1.9 2009/11/23 23:44:32 ogrebeny
265 // At Pibero's request, for the embedding infrastructure.
266 //
267 // Revision 1.8 2009/02/05 20:06:52 ogrebeny
268 // Changed StEEmcDbMaker -> StEEmcDb
269 //
270 // Revision 1.7 2007/04/28 17:56:02 perev
271 // Redundant StChain.h removed
272 //
273 // Revision 1.6 2007/03/23 03:26:23 balewski
274 // Corretions from Victor
275 //
276 // Revision 1.5 2007/02/07 02:24:34 balewski
277 // fix logic error found by Wei-Ming
278 //
279 // Revision 1.4 2007/01/24 21:07:03 balewski
280 // 1) no cout or printf, only new Logger
281 // 2) EndcapMixer:
282 // - no assert()
283 // - locks out on first fatal error til the end of the job
284 //
285 // Revision 1.3 2006/12/22 15:20:45 balewski
286 // ignore M-C hits for crates masked for real events, as suggested by Wei-Ming
287 //
288 // Revision 1.2 2006/12/13 13:24:38 balewski
289 // fix wrong header path
290 //
291 // Revision 1.1 2006/12/12 20:29:13 balewski
292 // added hooks for Endcap embedding
293 //
294 // Revision 1.1.1.1 2005/05/31 18:53:25 wzhang
295 // First version
296 //
297 //
const EEmcDbItem * StBarrelIndex2Item(int StDetId, int Bmod, int Beta, int Bsub) const
Definition: StEEmcDb.cxx:938
virtual Int_t Make()
virtual Int_t Finish()
virtual const char * GetName() const
special overload
Definition: StMaker.cxx:237
Definition: Stypes.h:42
Definition: Stypes.h:40
virtual Int_t Finish()
Definition: StMaker.cxx:776
Definition: Stypes.h:44