]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/EveDet/AliEveEMCALData.cxx
Coding convention fixes from Magali
[u/mrichter/AliRoot.git] / EVE / EveDet / AliEveEMCALData.cxx
1 //
2 // Fill containers for visualisation of EMCAL data structures
3 //    - read and store MC Hits    - read and store digits from esds or runloader
4 //    - read and store clusters from esds or runloader 
5 //
6 // Author: Magali Estienne (magali.estienne@cern.ch)
7 // June 30 2008
8 //
9
10 //#include <Riostream.h>
11 //#include <vector>
12
13 #include <TTree.h>
14 #include <TBranch.h>
15 #include <TObjArray.h>
16 #include <TRefArray.h>
17 #include <TClonesArray.h>
18 #include <TMath.h>
19 #include <TLorentzVector.h>
20
21 #include "AliRunLoader.h"
22 #include "AliEMCALLoader.h"
23 #include "AliESDVertex.h"
24 #include "AliEMCALHit.h"
25 #include "AliEMCALDigit.h"
26
27 #include "AliEMCALRecPoint.h"
28 #include "AliESDCaloCells.h"
29 #include "AliESDCaloCluster.h"
30
31 #include "AliEveEMCALData.h"
32 #include "AliEveEMCALSModuleData.h"
33
34 class Riostream;
35 class TObject;
36 class TEveUtil;
37 class TEvePointSet;
38 class AliRun;
39 class AliESDEvent;
40 class AliEMCAL;
41 class AliEMCALGeometry;
42 class AliEveEMCALSModule;
43
44
45 ClassImp(AliEveEMCALData)
46
47 //______________________________________________________________________________
48 AliEveEMCALData::AliEveEMCALData():
49   TObject(),
50   TEveRefCnt(),
51   fEmcal(0x0),
52   fGeom(0x0),
53   fNode(0x0),
54   fHMatrix(0),
55   fTree(0x0),
56   fESD(0x0),
57   fNsm(12),
58   fNsmfull(10),
59   fNsmhalf(2),
60   fSM(12),
61   fSMfull(10),
62   fSMhalf(2),
63   fRunLoader(0),
64   fDebug(0),
65   fPoint(0)
66 {
67   
68   //
69   // Constructor
70   //
71   //   for(Int_t i=0; i<12; i++)
72   //     fSM[i] = 0x0;
73   //   for(Int_t i=0; i<10; i++)
74   //     fSMfull[i] = 0x0;
75   //   fSMhalf[0] = fSMhalf[1] = 0x0;
76   //  InitEMCALGeom();
77   CreateAllSModules();
78
79
80 }
81
82 //______________________________________________________________________________
83 AliEveEMCALData::AliEveEMCALData(AliRunLoader* rl, TGeoNode* node, TGeoHMatrix* m):
84   TObject(),
85   TEveRefCnt(),
86   fEmcal(0x0),
87   fGeom(0x0),
88   fNode(node),
89   fHMatrix(m),
90   fTree(0x0),
91   fESD(0x0),
92   fNsm(12),
93   fNsmfull(10),
94   fNsmhalf(2),
95   fSM(12),
96   fSMfull(10),
97   fSMhalf(2),
98   fRunLoader(rl),
99   fDebug(0),
100   fPoint(0)
101 {
102   
103   //
104   // Constructor
105   //
106 //   for(Int_t i=0; i<12; i++)
107 //     fSM[i] = 0x0;
108 //   for(Int_t i=0; i<10; i++)
109 //     fSMfull[i] = 0x0;
110 //   fSMhalf[0] = fSMhalf[1] = 0x0;
111   InitEMCALGeom(rl);
112   CreateAllSModules();
113
114
115 }
116
117 //______________________________________________________________________________
118 AliEveEMCALData::~AliEveEMCALData()
119 {
120   //
121   // Destructor
122   //
123
124   DeleteSuperModules();
125   delete fTree;
126   // delete fEmcal; // deleted by run-loader
127   delete fGeom;
128   delete fNode;
129   delete fHMatrix;
130   delete fPoint;
131 }
132
133 //______________________________________________________________________________
134 AliEveEMCALData::AliEveEMCALData(const AliEveEMCALData &edata) :
135   TObject(edata),
136   TEveRefCnt(edata),
137   fEmcal(edata.fEmcal),
138   fGeom(edata.fGeom),
139   fNode(edata.fNode),
140   fHMatrix(edata.fHMatrix),
141   fTree(edata.fTree),
142   fESD(edata.fESD),
143   fNsm(edata.fNsm),
144   fNsmfull(edata.fNsmfull),
145   fNsmhalf(edata.fNsmhalf),
146   fSM(edata.fSM),
147   fSMfull(edata.fSMfull),
148   fSMhalf(edata.fSMhalf),
149   fRunLoader(edata.fRunLoader),
150   fDebug(edata.fDebug),
151   fPoint(edata.fPoint)
152 {
153   //
154   // Copy constructor
155   //
156   InitEMCALGeom(edata.fRunLoader);
157   CreateAllSModules();
158 }
159
160 //______________________________________________________________________________
161 AliEveEMCALData& AliEveEMCALData::operator=(const AliEveEMCALData &edata)
162 {
163   //
164   // Assignment operator
165   //
166
167   if (this != &edata) {
168
169   }
170
171   return *this;
172
173 }
174
175 //______________________________________________________________________________
176 void AliEveEMCALData::SetTree(TTree* const tree)
177 {
178   //
179   // Set digit-tree to be used for digit retrieval. 
180   // Data is loaded on demand.
181   // 
182
183   fTree = tree;
184
185 }
186
187 //______________________________________________________________________________
188 void AliEveEMCALData::SetESD(AliESDEvent* const esd)
189 {
190   //
191   // Set esd
192   //
193
194   fESD = esd;
195 }
196
197 //______________________________________________________________________________
198 void AliEveEMCALData::SetNode(TGeoNode* const node)
199 {
200   //
201   // Set node
202   //
203
204   fNode = node;
205 }
206
207 //______________________________________________________________________________
208 void AliEveEMCALData::InitEMCALGeom(AliRunLoader* const rl)
209 {
210   //
211   // Set data members for EMCAL geometry
212   //
213
214   fEmcal = (AliEMCAL*) rl->GetAliRun()->GetDetector("EMCAL");
215   fGeom  = (AliEMCALGeometry*) fEmcal->GetGeometry();
216
217 }
218
219 //______________________________________________________________________________
220 void AliEveEMCALData::GetGeomInfo(Int_t id, Int_t &iSupMod, Double_t& x, Double_t& y, Double_t& z)
221 {
222   //
223   // Get geometrical information from hit/digit/cluster absolute id
224   //
225
226   Int_t iTower  =  0 ;
227   Int_t iIphi   =  0 ;
228   Int_t iIeta   =  0 ;
229
230   //Geometry methods
231   fGeom->GetCellIndex(id,iSupMod,iTower,iIphi,iIeta);
232   //Gives SuperModule and Tower numbers
233   fGeom->RelPosCellInSModule(id, x, y, z);
234
235 }
236
237 //______________________________________________________________________________
238 void  AliEveEMCALData::CreateAllSModules()
239 {
240   //
241   // Create all fNsm super modules
242   //
243
244   for(Int_t sm = 0; sm < fNsm; sm++)
245     CreateSModule(sm);
246
247 }
248
249 //______________________________________________________________________________
250 void  AliEveEMCALData::CreateSModule(Int_t sm)
251 {
252   //
253   // Create super-module-data for SM if it does not exist already.
254   //
255
256   if(fSM[sm] == 0) fSM[sm] = new AliEveEMCALSModuleData(sm,fGeom,fNode,fHMatrix);
257   if(fSMfull[sm] == 0 && sm < 10) fSMfull[sm] = new AliEveEMCALSModuleData(sm,fGeom,fNode,fHMatrix);
258   if(fSMhalf[sm-10] == 0 && sm > 10) fSMhalf[sm-10] = new AliEveEMCALSModuleData(sm,fGeom,fNode,fHMatrix);
259 }
260
261 //______________________________________________________________________________
262 void AliEveEMCALData::DropAllSModules()
263 {
264   //
265   // Drop data of all existing sectors.
266   //
267
268   for (Int_t sm = 0; sm < fNsm; sm++) {
269     if (fSM[sm] != 0)
270       fSM[sm]->DropData();
271   }
272 }
273
274 //______________________________________________________________________________
275 void AliEveEMCALData::DeleteSuperModules()
276 {
277   //
278   // Delete all super module data
279   //
280
281   for (Int_t sm = 0; sm < fNsm; sm++)
282     {
283       fSM[sm] = 0;
284       delete fSM[sm];
285     }
286
287   for(Int_t smf = 0; smf < fNsmfull; smf++) 
288     {
289       fSMfull[smf] = 0;
290       delete fSMfull[smf];
291     }
292
293   for(Int_t smh = 0; smh < fNsmhalf; smh++)
294     {
295       fSMhalf[smh] = 0;
296       delete fSMhalf[smh];
297     }
298   
299 }
300
301 //______________________________________________________________________________
302 void AliEveEMCALData::LoadHits(TTree* const t)
303 {
304   //
305   // Get hit information from RunLoader
306   //
307
308   /*
309   // These are global coordinates !
310   char form[1000];
311   const char *selection = "";
312   const char *varexp = "fX:fY:fZ";
313   sprintf(form,"EMCAL Hits '%s'", selection);
314   fPoint = new TEvePointSet(form);
315
316   TEvePointSelector ps(t, fPoint, varexp, selection);
317   ps.Select();
318
319   if (fPoint->Size() == 0) {
320     Warning("emcal_hits", Form("No hits match '%s'", selection));
321     delete fPoint;
322     //    return 0;
323   }
324   */
325
326   TObjArray *harr=NULL;
327   TBranch *hbranch=t->GetBranch("EMCAL");
328   hbranch->SetAddress(&harr);
329   
330   if(hbranch->GetEvent(0)) {
331     for(Int_t ih = 0; ih < harr->GetEntriesFast(); ih++) {
332       AliEMCALHit* hit =(AliEMCALHit*)harr->UncheckedAt(ih);
333       if(hit != 0){
334         if(fDebug>1) cout << "Hit info " << hit->GetId() << " " << hit->GetEnergy() << endl;
335         Int_t id = hit->GetId();
336         // These are local coordinates
337         Double_t xl = 0.; Double_t yl = 0.; Double_t zl = 0.;
338         // Get global coordinates
339         Double_t x = hit->X();
340         Double_t y = hit->Y();
341         Double_t z = hit->Z();
342         Double_t amp = hit->GetEnergy();
343         Int_t iSupMod = 0;
344         // Get SM Id
345         GetGeomInfo(id,iSupMod,xl,yl,zl);
346         fSM[iSupMod]->RegisterHit(id,iSupMod,amp,x,y,z);
347       }
348     }
349   }
350 }
351
352 //______________________________________________________________________________
353 void AliEveEMCALData::LoadHitsFromEMCALLoader(AliEMCALLoader* const emcl)
354 {
355   //
356   // Get hit information from EMCAL Loader
357   //
358
359   AliEMCALHit* hit;
360   
361   //Fill array of hits                                                                        
362   TClonesArray *hits = (TClonesArray*)emcl->Hits();
363
364   //Get hits from the list                                                                    
365   for(Int_t ihit = 0; ihit< hits->GetEntries();ihit++){
366
367     hit = static_cast<AliEMCALHit *>(hits->At(ihit)) ;
368     
369     if(hit != 0){
370       if(fDebug>1) cout << "Hit info " << hit->GetId() << " " << hit->GetEnergy() << endl;
371
372       Int_t id = hit->GetId();
373       // These are local coordinates
374       Double_t xl = 0.; Double_t yl = 0.; Double_t zl = 0.;
375       // Get global coordinates
376       Double_t x = hit->X();
377       Double_t y = hit->Y();
378       Double_t z = hit->Z();
379       Double_t amp = hit->GetEnergy();
380       Int_t iSupMod = 0;
381       // Get SM Id
382       GetGeomInfo(id,iSupMod,xl,yl,zl);
383       fSM[iSupMod]->RegisterHit(id,iSupMod,amp,x,y,z);
384     }
385   }
386   
387 }
388
389 //______________________________________________________________________________
390 void AliEveEMCALData::LoadDigits(TTree *t)
391 {
392   //
393   // Get digit information from RunLoader
394   //
395
396   TClonesArray *digits = 0;
397   t->SetBranchAddress("EMCAL", &digits);
398   t->GetEntry(0);
399
400   Int_t nEnt = digits->GetEntriesFast();
401   cout << "nEnt: " << nEnt << endl;
402   AliEMCALDigit * dig;
403
404   //  Double_t amp   = -1 ;
405   Double_t ampInt   = -1 ;
406   Int_t id      = -1 ;
407   Int_t iSupMod =  0 ;
408   Double_t x, y, z;
409
410   for (Int_t idig = 0; idig < nEnt; idig++)
411     {
412       dig = static_cast<AliEMCALDigit *>(digits->At(idig));
413       
414       if(dig != 0) {
415         id   = dig->GetId() ; //cell (digit) label
416         // adc
417         ampInt  = dig->GetAmp(); //amplitude in cell (digit)
418         // GeV
419         //      amp = ampInt*0.0153; // To be modified with correct OCDB conversion     
420
421         GetGeomInfo(id,iSupMod,x,y,z);
422
423 //      // GeV
424 //      fSM[iSupMod]->RegisterDigit(id,iSupMod,amp,x,y,z);
425 // //   fSM[iSupMod]->SaveDigit(dig);
426 // //   if(iSupMod<fNsmfull) fSMfull[iSupMod]->RegisterDigit(id,iSupMod,amp,x,y,z);
427 // //   if(iSupMod>fNsmfull) fSMhalf[iSupMod-10]->RegisterDigit(id,iSupMod,amp,x,y,z);
428         fSM[iSupMod]->RegisterDigit(id,iSupMod,ampInt,x,y,z);
429 //      fSM[iSupMod]->SaveDigit(dig);
430 //      if(iSupMod<fNsmfull) fSMfull[iSupMod]->RegisterDigit(id,iSupMod,amp,x,y,z);
431 //      if(iSupMod>fNsmfull) fSMhalf[iSupMod-10]->RegisterDigit(id,iSupMod,amp,x,y,z);
432       }
433       else {
434         cout << "Digit object empty" << endl;
435         return;
436       }
437     } // end loop digits
438   cout << "after loop on digits !" << endl;
439 }
440
441 //______________________________________________________________________________
442 void AliEveEMCALData::LoadDigitsFromEMCALLoader(AliEMCALLoader* const emcl)
443 {
444
445   //
446   // Get digit information from EMCAL Loader
447   //
448
449   AliEMCALDigit* dig;
450   
451   //Fill array of digits                                                                        
452   TClonesArray *digits = (TClonesArray*)emcl->Digits();
453   
454   //Get digits from the list  
455   
456   //  Double_t amp   = -1 ;
457   Double_t ampInt   = -1 ;
458   Int_t id      = -1 ;
459   Int_t iSupMod =  0 ;
460   Double_t x, y, z;
461                                                
462    for(Int_t idig = 0; idig< digits->GetEntries();idig++){
463
464      dig = static_cast<AliEMCALDigit *>(digits->At(idig)) ;
465
466      if(dig != 0){
467        if(fDebug>1) cout << "Digit info " << dig->GetId() << " " << dig->GetAmp() << endl;
468        id   = dig->GetId() ; //cell (digit) label
469        // adc
470        ampInt  = dig->GetAmp(); //amplitude in cell (digit)
471        // GeV
472        //       amp = ampInt*0.0153.; // To be modified with correct OCDB conversion
473
474        GetGeomInfo(id,iSupMod,x,y,z);
475        
476        //       // GeV
477        //       fSM[iSupMod]->RegisterDigit(id,iSupMod,amp,x,y,z);
478        // adc
479        fSM[iSupMod]->RegisterDigit(id,iSupMod,ampInt,x,y,z);
480      }
481       else {
482         cout << "Digit object empty" << endl;
483         return;
484       }
485    } // end loop on digits
486    
487 }
488
489 //______________________________________________________________________________
490 void AliEveEMCALData::LoadDigitsFromESD()
491 {
492   //
493   // Get digit information from esd
494   //
495   
496   AliESDCaloCells &cells= *(fESD->GetEMCALCells());
497   Int_t ncell = cells.GetNumberOfCells() ;  
498   Int_t iSupMod =  0 ;
499   Double_t x, y, z;
500
501   // Extract digit information from the ESDs
502   for (Int_t icell=  0; icell <  ncell; icell++) 
503     {
504       Int_t id   = cells.GetCellNumber(icell);
505       // adc
506       Double_t ampInt  = cells.GetAmplitude(icell);
507       // GeV
508       //      Double_t amp = ampInt*0.0153; // To be modified with correct OCDB conversion
509
510       GetGeomInfo(id,iSupMod,x,y,z);
511
512 //       // GeV
513 //       fSM[iSupMod]->RegisterDigit(id,iSupMod,amp,x,y,z);
514 //       if(iSupMod<fNsmfull) fSMfull[iSupMod]->RegisterDigit(id,iSupMod,amp,x,y,z);
515 //       if(iSupMod>fNsmfull) fSMhalf[iSupMod-10]->RegisterDigit(id,iSupMod,amp,x,y,z);
516       // adc
517       fSM[iSupMod]->RegisterDigit(id,iSupMod,ampInt,x,y,z);
518       if(iSupMod<fNsmfull) fSMfull[iSupMod]->RegisterDigit(id,iSupMod,ampInt,x,y,z);
519       if(iSupMod>fNsmfull) fSMhalf[iSupMod-10]->RegisterDigit(id,iSupMod,ampInt,x,y,z);
520
521     } // end loop cells
522 }
523
524 //______________________________________________________________________________
525 void AliEveEMCALData::LoadRecPoints(TTree* const t)
526 {
527   //
528   // Get rec point information from RunLoader
529   //
530
531   //*************************************************
532   // To be improved !!!!!
533   // Size and shape of cluster to be implemented
534   // 
535   //*************************************************
536         
537   // From TTreeR
538   TObjArray *carr=NULL;
539   TBranch *cbranch=t->GetBranch("EMCALECARP");
540   cbranch->SetAddress(&carr);
541   
542   if(cbranch->GetEvent(0)) {
543     for(Int_t ic = 0; ic < carr->GetEntriesFast(); ic++) {
544       AliEMCALRecPoint* rp =(AliEMCALRecPoint*)carr->UncheckedAt(ic);
545       if(rp){
546         if(fDebug>1) cout << "RecPoint info " << rp->GetAbsId() << " " << rp->GetEnergy() << endl;
547         Int_t iSupMod = rp->GetSuperModuleNumber();
548         // GeV
549         Double_t amp = (Double_t)rp->GetEnergy();
550         // adc
551         Double_t ampInt = amp/0.0153; // To be modified with correct OCDB conversion
552         TVector3 lpos;
553         rp->GetLocalPosition(lpos);
554
555 //      // GeV
556 //      fSM[iSupMod]->RegisterCluster(iSupMod,amp,lpos[0],lpos[1],lpos[2]);
557         // adc
558         fSM[iSupMod]->RegisterCluster(iSupMod,ampInt,lpos[0],lpos[1],lpos[2]);
559       }
560     }
561   }
562   
563 }
564
565 //______________________________________________________________________________
566 void AliEveEMCALData::LoadRecPointsFromEMCALLoader(AliEMCALLoader* const emcl)
567 {
568   //
569   // Get rec point information from EMCAL Loader
570   //
571
572   //*************************************************
573   // To be improved !!!!!
574   // Size and shape of cluster to be implemented
575   // 
576   //*************************************************
577
578   // From EMCALLoader
579   AliEMCALRecPoint* rp;
580   
581   //Fill array of clusters                                                                        
582   TClonesArray *clusters = (TClonesArray*)emcl->RecPoints();
583   
584   //Get clusters from the list                                                                    
585   for(Int_t iclu = 0; iclu< clusters->GetEntries();iclu++){
586
587     rp = static_cast<AliEMCALRecPoint *>(clusters->At(iclu)) ;
588     
589     if(rp){
590        if(fDebug>1) cout << "RecPoint info " << rp->GetAbsId() << " " << rp->GetEnergy() << endl;
591        Int_t iSupMod = rp->GetSuperModuleNumber();
592        Double_t amp = (Double_t)rp->GetEnergy();
593        Double_t ampInt = amp/0.0153; // To be modified with correct OCDB conversion
594        TVector3 lpos;
595        rp->GetLocalPosition(lpos);
596        
597 //        // GeV
598 //        fSM[iSupMod]->RegisterCluster(iSupMod,amp,lpos[0],lpos[1],lpos[2]);
599        // adc
600        fSM[iSupMod]->RegisterCluster(iSupMod,ampInt,lpos[0],lpos[1],lpos[2]);
601     }
602   }
603   
604 }
605
606 //______________________________________________________________________________
607 void AliEveEMCALData::LoadRecPointsFromESD()
608 {
609   //
610   // Get cluster information from esd
611   //
612
613  Int_t iSupMod =  0 ;
614   Double_t x, y, z;
615   Int_t iSM =  0 ;
616   Int_t iT  =  0 ;
617   Int_t iIp   =  0 ;
618   Int_t iIe   =  0 ;
619   Double_t xd, yd, zd;
620   Float_t pos[3] ; 
621   
622   // Get reconstructed vertex position
623   AliESDVertex* primVertex =(AliESDVertex*) fESD->GetVertex();
624   Double_t vertexPosition[3] ; 
625   primVertex->GetXYZ(vertexPosition) ; 
626
627   //Get the CaloClusters
628   //select EMCAL clusters only 
629   TRefArray * caloClusters  = new TRefArray();
630   fESD->GetEMCALClusters(caloClusters);
631   Int_t nclus = caloClusters->GetEntries();
632   cout << "nclus: " << nclus << endl; 
633   
634   if(!caloClusters) return;
635
636   for (Int_t iclus =  0; iclus <  nclus; iclus++) 
637     {
638       AliESDCaloCluster *clus = (AliESDCaloCluster *) caloClusters->At(iclus) ; 
639       //Get the cluster info
640       
641       Double_t energy = clus->E() ;
642       // adc
643       //      Int_t   eneInt = (Int_t)energy*500+0.5;
644       Double_t eneInt = energy/0.0153; // To be modified with correct OCDB conversion
645       Double_t disp   = clus->GetClusterDisp() ;
646       
647       clus->GetPosition(pos) ; // Global position
648       TVector3 vpos(pos[0],pos[1],pos[2]) ;
649       TLorentzVector p4 ;
650       TVector3 p3;
651       clus->GetMomentum(p4,vertexPosition);
652       p3.SetXYZ(p4[0],p4[1],p4[2]);
653       Double_t eta = p3.Eta();
654       Double_t phi = ( (p3.Phi()) < 0) ? (p3.Phi()) + 2. * TMath::Pi() : (p3.Phi());
655       
656       Int_t mult = clus->GetNCells() ;
657       if(fDebug>2) {
658         cout << "In cluster: " << iclus << ", ncells: " << mult << ", energy : " << energy << 
659           ", disp: " << disp << endl;
660         cout << "Cluster " << iclus << ", eta: " << eta << ", phi: " << phi << endl;
661       }
662
663       Int_t clusId = 0;
664       fGeom->GetAbsCellIdFromEtaPhi(eta,phi,clusId);
665       if(fDebug>2) {
666         cout << "Abs Cluster Id: " << clusId << ", xc: " << pos[0] << 
667           ", yc: " << pos[1] << ", zc: " << pos[2] << endl;
668       }
669
670       GetGeomInfo(clusId,iSupMod,x,y,z);
671       
672       //******** Not used yet but will come  ********
673       // AliESDCaloCells &cells= *(fESD->GetEMCALCells());
674       Int_t     digMult = clus->GetNCells() ;
675       UShort_t *digID   = clus->GetCellsAbsId() ;
676       for(Int_t i=0; i<digMult; i++){
677         // Float_t  digitAmp     = cells.GetCellAmplitude(digID[i]) ;
678         fGeom->RelPosCellInSModule(digID[i], xd, yd, zd);
679         //Geometry methods
680         fGeom->GetCellIndex(digID[i],iSM,iT,iIp,iIe);
681         //Gives SuperModule and Tower numbers
682
683       } // end digit loop
684       //*********************************************
685       //      // GeV
686       //      fSM[iSupMod]->RegisterCluster(iSM,energy,x,y,z);
687       // adc
688       fSM[iSupMod]->RegisterCluster(iSM,eneInt,x,y,z);
689
690     } // end cluster loop
691 }
692
693 //______________________________________________________________________________
694 AliEveEMCALSModuleData* AliEveEMCALData::GetSModuleData(Int_t sm)
695 {
696   //
697   // Return super module data
698   //
699
700   if (sm < 0 || sm > fNsm) 
701     {
702       printf("The number of super modules must be lower or equal to %d",fNsm);
703       return 0;
704     }
705   
706   return fSM[sm];
707 }
708
709 //______________________________________________________________________________
710 void AliEveEMCALData::LoadRaw() const
711 {
712   //
713   // Get raw information
714   //
715
716   // To be implemented !
717 }
718