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