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