]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDAlla.cxx
3c51d7b49e393f40f9aafd8bb5cc69b0fb11eb78
[u/mrichter/AliRoot.git] / FMD / AliFMDAlla.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 /////////////////////////////////////////////////////////////////////
19 //                                                                 //
20 // Forward Multiplicity detector based on Silicon version 0        //
21 //
22 //Begin Html       
23 /*
24   <img src="gif/AliFMDv0Class.gif">
25 */
26 //End Html
27 //                                                                  //
28 //                                                                  //
29 //////////////////////////////////////////////////////////////////////
30
31 #include <Riostream.h>
32 #include <stdlib.h>
33
34 #include <TClonesArray.h>
35 #include <TDirectory.h>
36 #include <TFile.h>
37 #include <TGeometry.h>
38 #include <TLorentzVector.h>
39 #include <TMath.h>
40 #include <TNode.h>
41 #include <TTUBE.h>
42 #include <TTree.h>
43 #include <TVirtualMC.h>
44 #include <TDatabasePDG.h>
45 #include <AliFMDHit.h>
46 #include "AliFMDv0.h"
47 #include "AliFMDAlla.h"
48 #include "AliFMDDigitizerAlla.h"
49 #include "AliMagF.h"
50 #include "AliRun.h"
51 #include "AliMC.h"
52
53 ClassImp(AliFMDAlla)
54 #if 0
55   ;
56 #endif
57
58 //--------------------------------------------------------------------
59 AliFMDAlla::AliFMDAlla(const char *name, const char *title)
60   : AliFMD(name,title)
61 {
62   //
63   // Standart constructor for Forward Multiplicity Detector version 0
64   //
65   fIdSens1=0;
66   fIdSens2=0;
67   fIdSens3=0;
68   fIdSens4=0;
69   fIdSens5=0;
70   //  setBufferSize(128000);
71 }
72 //-------------------------------------------------------------------------
73 void 
74 AliFMDAlla::CreateGeometry()
75 {
76   //
77   // Create the geometry of Forward Multiplicity Detector version 0
78   //
79   //Detector consists of 6 volumes: 
80   // 1st covered pseudorapidity interval from 3.3 to 2.0
81   // and placed on 65cm in Z-direction;
82   // 2nd - from 2.0 to 1.6 and Z=85 cm;
83   // 3d  - the same pseudorapidity interval as the 1st 
84   // but on the other side from the interaction point z=-65cm;
85   // 4th - simmetricaly with the 2nd : 
86   // pseudorapidity from 2.0 to 1.6, Z=-85cm   
87   // 5th - from 3.6 to 4.7, Z=-270cm
88   // 6th - from 4.5 to 5.5 , Z=-630cm.
89   // Each part has 400mkm Si (sensetive area, detector itself),
90   // 0.75cm of plastic simulated electronics material,
91   // Al support ring 2cm thickness and 1cm width placed on 
92   // the outer radius of each Si disk;
93   //    
94   // begin Html
95   /*
96     <img src="gif/AliFMDv0.gif">
97   */
98   //
99
100
101
102   Int_t *idtmed = fIdtmed->GetArray();
103    
104   Int_t ifmd;
105   Int_t idrotm[999];
106   Float_t zFMD;
107   Float_t par[3];
108   Float_t ppcon[15];
109   Float_t z[5]             = {-62.8, -75.2, 83.4, 75.2, 340.};
110   Float_t NylonTube[3]     = {0.2,0.6,0.45};
111   Float_t zPCB             = 0.12; 
112   Float_t zHoneyComb       = 0.5; 
113   Float_t zSi              = 0.03;
114   Float_t rin[5]           = {4.2,15.4,4.2,15.4,4.2};
115   Float_t rout[5]          = {17.4,28.4,17.4,28.4,17.4};
116   Float_t RinHoneyComb[5]  = { 5.15,16.4,  5.15,16.4,  5.15};
117   Float_t RoutHoneyComb[5] = {20.63,34.92,22.3, 32.02,20.63};
118   Float_t zInside;
119   Float_t zCooper          = 0.01; 
120   Float_t zChips           = 0.01;
121   Float_t yNylonTube[5]    = {10,20,10,20,10};
122  
123   char    nameFMD[5];
124   char    nameSi[5];
125   char    nameSector[5];
126   char    nameRing[5];
127   Char_t  nameHoney[5];
128   char    nameHoneyIn[5];
129   char    nameHoneyOut[5];
130   Char_t  namePCB[5];
131   char    nameCopper[5];
132   char    nameChips[5];
133   char    nameG10[5];
134   Char_t  nameLPCB[5];
135   char    nameLCopper[5];
136   char    nameLChips[5];
137   char    nameGL10[5];
138
139   AliMatrix(idrotm[901], 90, 0, 90, 90, 180, 0);
140   
141   
142   // Nylon tubes
143   gMC->Gsvolu("GNYL","TUBE", idtmed[1], NylonTube, 3);  //support nylon tube
144   Float_t wideSupport=zSi+3*zPCB+2*NylonTube[2]+zHoneyComb;
145   AliDebug(1, Form("Support width: %f", wideSupport));
146
147   for (ifmd=0; ifmd<5; ifmd++)  {
148     sprintf(nameFMD,"FMD%d",ifmd+1);
149     ppcon[0]=0;
150     ppcon[1]=360;
151     ppcon[2]=4;
152       
153     ppcon[3]=-wideSupport;
154     ppcon[4]=rin[ifmd]-0.1;
155     ppcon[5]=rout[ifmd]+0.1;
156       
157     ppcon[6]=ppcon[3]+2*zSi+2*zPCB+2*NylonTube[2];
158     ppcon[7]=rin[ifmd]-0.1;
159     ppcon[8]=rout[ifmd]+0.1;
160       
161     ppcon[9]=ppcon[6];
162     ppcon[10]=RinHoneyComb[ifmd]-0.1;
163     ppcon[11]=RoutHoneyComb[ifmd]+0.1;
164
165     ppcon[12]=ppcon[9]+2*zHoneyComb+zPCB;
166     ppcon[13]=RinHoneyComb[ifmd]-0.1;
167     ppcon[14]=RoutHoneyComb[ifmd]+0.1;
168     gMC->Gsvolu(nameFMD,"PCON",idtmed[0],ppcon,15);
169     if (z[ifmd] >0){  
170       zFMD=z[ifmd]+wideSupport;
171       gMC->Gspos(nameFMD,1,"ALIC",0,0,zFMD,0, "ONLY");}
172     else {
173       zFMD=z[ifmd]-wideSupport;
174       gMC->Gspos(nameFMD,1,"ALIC",0,0,zFMD,idrotm[901], "ONLY");}
175
176     //silicon
177     sprintf(nameSi,"GSI%d",ifmd+1);
178     sprintf(nameSector,"GSC%d",ifmd+1);
179     sprintf(nameRing,"GRN%d",ifmd+1);
180       
181     //honeycomb support
182     sprintf(nameHoney,"GSU%d",ifmd+1);
183     gMC->Gsvolu(nameHoney,"TUBE", idtmed[0], par, 0);  //honeycomb 
184     sprintf(nameHoneyIn,"GHI%d",ifmd+1);
185     gMC->Gsvolu(nameHoneyIn,"TUBE", idtmed[7], par, 0);  //honey comb inside 
186     sprintf(nameHoneyOut,"GHO%d",ifmd+1);
187     gMC->Gsvolu(nameHoneyOut,"TUBE", idtmed[6], par, 0);  //honey comb skin
188
189     //PCB
190     sprintf(namePCB,"GPC%d",ifmd+1);
191     gMC->Gsvolu(namePCB,"TUBE", idtmed[0], par, 0); //PCB
192     sprintf(nameCopper,"GCO%d",ifmd+1);
193     gMC->Gsvolu(nameCopper,"TUBE", idtmed[3], par, 0);  // Cooper
194     sprintf(nameChips,"GCH%d",ifmd+1);
195     gMC->Gsvolu(nameChips,"TUBE", idtmed[5], par, 0); // Si chips
196     sprintf(nameG10,"G10%d",ifmd+1);
197     gMC->Gsvolu(nameG10,"TUBE", idtmed[2], par, 0);  //G10 plate
198
199     //last PCB
200     sprintf(nameLPCB,"GPL%d",ifmd+1);
201     gMC->Gsvolu(nameLPCB,"TUBE", idtmed[0], par, 0); //PCB
202     sprintf(nameLCopper,"GCL%d",ifmd+1);
203     gMC->Gsvolu(nameLCopper,"TUBE", idtmed[3], par, 0);  // Cooper
204     sprintf(nameLChips,"GHL%d",ifmd+1);
205     gMC->Gsvolu(nameLChips,"TUBE", idtmed[5], par, 0); // Si chips
206     sprintf(nameGL10,"G1L%d",ifmd+1);
207     gMC->Gsvolu(nameGL10,"TUBE", idtmed[2], par, 0); // Last G10
208     par[0]=rin[ifmd]; // pipe size
209     par[1]=rout[ifmd];
210     par[2]=zSi/2;
211     gMC->Gsvolu(nameSi,"TUBE", idtmed[4], par, 3);
212     zInside=ppcon[3]+par[2];
213     gMC->Gspos(nameSi,ifmd+1,nameFMD,0,0,zInside,0, "ONLY");
214
215     //PCB 1
216     zInside += par[2]+zPCB/2;
217     par[2]=zPCB/2;
218     gMC->Gsposp(namePCB,1,nameFMD,0,0,zInside,0, "ONLY",par,3);
219     zInside += zPCB;
220     gMC->Gsposp(namePCB,2,nameFMD,0,0,zInside,0, "ONLY",par,3);
221     Float_t NulonTubeBegin=zInside+2.5*zPCB;
222     par[2]=zPCB/2-0.02;
223     Float_t zInPCB = -zPCB/2+par[2];
224     gMC->Gsposp(nameG10,1,namePCB,0,0,zInPCB,0, "ONLY",par,3);
225     zInPCB+=par[2]+zCooper/2 ;
226     par[2]=zCooper/2;
227     gMC->Gsposp(nameCopper,1,namePCB,0,0,zInPCB,0, "ONLY",par,3);
228     zInPCB += zCooper/2 + zChips/2;
229     par[2]=zChips/2;
230     gMC->Gsposp(nameChips,1,namePCB,0,0,zInPCB,0, "ONLY",par,3);
231
232     //HoneyComb
233     zHoneyComb=0.8;   
234     par[0] = RinHoneyComb[ifmd];
235     par[1] = RoutHoneyComb[ifmd];
236     par[2] = zHoneyComb/2;
237     zInside += 2*NylonTube[2]+par[2];
238     gMC->Gsposp(nameHoney,1,nameFMD,0,0,zInside,0, "ONLY",par,3);
239     par[2]=0.1/2;
240     Float_t zHoney=-zHoneyComb/2+par[2];
241     gMC->Gsposp(nameHoneyOut,1,nameHoney,0,0,zHoney,0,
242                 "ONLY",par,3); //shkurki
243     zHoney=zHoneyComb/2-par[2];
244     gMC->Gsposp(nameHoneyOut,2,nameHoney,0,0,zHoney,0, "ONLY",par,3);
245     par[2]=(zHoneyComb-2.*0.1)/2; //soty vnutri
246     gMC->Gsposp(nameHoneyIn,1,nameHoney,0,0,0,0, "ONLY",par,3);
247       
248     gMC->Gspos("GNYL",1,nameFMD,0,yNylonTube[ifmd],
249                NulonTubeBegin+NylonTube[2]/2.,0, "ONLY");
250     gMC->Gspos("GNYL",2,nameFMD,0,-yNylonTube[ifmd],
251                NulonTubeBegin+NylonTube[2]/2.,0, "ONLY");
252          
253     //last PCB
254     par[0]=RoutHoneyComb[ifmd]-9;
255     par[1]=RoutHoneyComb[ifmd];
256     par[2]=zPCB/2;
257     zInside += zHoneyComb/2+par[2];
258     gMC->Gsposp(nameLPCB,1,nameFMD,0,0,zInside,0, "ONLY",par,3);
259       
260     par[2]=zPCB/2-0.02;
261     zInPCB = -zPCB/2+par[2];
262     gMC->Gsposp(nameGL10,1,nameLPCB,0,0,zInPCB,0, "ONLY",par,3);
263     zInPCB+=par[2]+zCooper/2 ;
264     par[2]=zCooper/2;
265     gMC->Gsposp(nameLCopper,1,nameLPCB,0,0,zInPCB,0, "ONLY",par,3);
266     zInPCB += zCooper/2 + zChips/2;
267     par[2]=zChips/2;
268     gMC->Gsposp(nameLChips,1,nameLPCB,0,0,zInPCB,0, "ONLY",par,3);
269       
270            
271     //Granularity
272     fSectorsSi1 = 20;
273     fRingsSi1   = 256 * 2;
274     //  fRingsSi1=3; // for drawing only
275     fSectorsSi2 = 40;
276     fRingsSi2   = 128 * 2;
277     //   fRingsSi2=3; //for  drawing onl
278     if(ifmd==1||ifmd==3) { 
279       gMC->Gsdvn(nameSector, nameSi , fSectorsSi2, 2);
280       gMC->Gsdvn(nameRing, nameSector, fRingsSi2, 1);
281     }
282     else {
283       gMC->Gsdvn(nameSector, nameSi , fSectorsSi1, 2);
284       gMC->Gsdvn(nameRing, nameSector , fRingsSi1, 1);
285     }
286   }
287 }    
288
289
290 //------------------------------------------------------------------------
291 void 
292 AliFMDAlla::CreateMaterials() 
293 {
294   Int_t isxfld   = gAlice->Field()->Integ();
295   Float_t sxmgmx = gAlice->Field()->Max();
296
297   // Plastic CH
298   Float_t aPlastic[2]={1.01,12.01};
299   Float_t zPlastic[2]={1,6};
300   Float_t wPlastic[2]={1,1};
301   Float_t denPlastic=1.03;
302   //
303   //     60% SiO2 , 40% G10FR4 
304   // PC board
305   Float_t apcb[3]  = { 28.0855,15.9994,17.749 };
306   Float_t zpcb[3]  = { 14.,8.,8.875 };
307   Float_t wpcb[3]  = { .28,.32,.4 };
308   Float_t denspcb  = 1.8;
309   //
310   // AIR
311   Float_t aAir[4]  = {12.0107,14.0067,15.9994,39.948};
312   Float_t zAir[4]  = {6.,7.,8.,18.};
313   Float_t wAir[4]  = {0.000124,0.755267,0.231781,0.012827};
314   Float_t dAir     = 1.20479E-3;
315   //*** Definition Of avaible FMD materials ***
316   AliMixture(0, "FMD Air$", aAir, zAir, dAir, 4,wAir);
317   AliMixture(1, "Plastic$",aPlastic,zPlastic,denPlastic,-2,wPlastic);
318   AliMixture(2, "SSD PCB$",   apcb, zpcb, denspcb, 3, wpcb);
319   AliMaterial(3, "SSD Copper$", 63.546, 29., 8.96, 1.43, 999.);
320   AliMaterial(4, "SSD Si$",      28.0855, 14., 2.33, 9.36, 999.);
321   AliMaterial(5, "SSD Si chip$", 28.0855, 14., 2.33, 9.36, 999.);
322   AliMaterial(6, "SSD C$",       12.011,   6., 2.265,18.8, 999.);
323   AliMaterial(7, "SSD Kapton$", 12.011, 6., 0.01, 31.27, 999.);//honeycomb
324   AliMaterial(8, "SSD G10FR4$", 17.749, 8.875, 1.8, 21.822, 999.);
325    
326
327   //**
328   AliMedium(0, "FMD air$", 0, 0, isxfld, sxmgmx, 1., .001, 1., .001, .001);
329   AliMedium(1, "Plastic$", 1, 0,isxfld, sxmgmx,  10., .01, 1., .003, .003);
330   AliMedium(2, "SSD PCB$", 2, 0, isxfld, sxmgmx, 1., .001, 1., .001, .001);
331   AliMedium(3, "SSD Copper$", 3, 0,isxfld, sxmgmx,  10., .01, 1., .003, .003);
332   AliMedium(4, "SSD Si$", 4, 1, isxfld, sxmgmx, 1., .001, 1., .001, .001);
333   AliMedium(5, "SSD Si chip$", 5, 0,isxfld, sxmgmx,  10., .01, 1., .003, .003);
334   AliMedium(6, "SSD C$", 6, 0,isxfld, sxmgmx,  10., .01, 1., .003, .003);
335   AliMedium(7, "SSD Kapton$", 7, 0, isxfld, sxmgmx, 1., .001, 1., .001, .001);
336   AliMedium(8, "SSD G10FR4$", 8, 0,isxfld, sxmgmx,  10., .01, 1., .003, .003);
337  
338
339
340 }
341 //---------------------------------------------------------------------
342 void 
343 AliFMDAlla::DrawDetector()
344 {
345   //
346   // Draw a shaded view of the Forward multiplicity detector version 0
347   //
348
349   //Set ALIC mother transparent
350   gMC->Gsatt("ALIC","SEEN",0);
351   //
352   //Set volumes visible
353   gMC->Gsatt("FMD1","SEEN",1);
354   gMC->Gsatt("FMD2","SEEN",1);
355   gMC->Gsatt("FMD3","SEEN",1);
356   gMC->Gsatt("FMD4","SEEN",1);
357   gMC->Gsatt("FMD5","SEEN",1);
358
359   //
360   gMC->Gdopt("hide","on");
361   gMC->Gdopt("shad","on");
362   gMC->SetClipBox(".");
363   gMC->SetClipBox("*",0,1000,-1000,1000,-1000,1000);
364   gMC->DefaultRange();
365   gMC->Gdraw("alic",40,30,0,12,9.5,.2,0.2);
366   gMC->Gdhead(1111,"Forward multiplicity detector");
367   gMC->Gdopt("hide","off");
368 }
369 //-------------------------------------------------------------------
370 void 
371 AliFMDAlla::Init()
372 {
373   // Initialises version 0 of the Forward Multiplicity Detector
374   //
375   AliFMD::Init();
376   fIdSens1=gMC->VolId("GRN1");
377   fIdSens2=gMC->VolId("GRN2");
378   fIdSens3=gMC->VolId("GRN3");
379   fIdSens4=gMC->VolId("GRN4");
380   fIdSens5=gMC->VolId("GRN5");
381   printf("*** FMD version 1 initialized ***\n");
382 }
383
384 //-------------------------------------------------------------------
385 void 
386 AliFMDAlla::StepManager()
387 {
388   //
389   // Called for every step in the Forward Multiplicity Detector
390   //
391   Int_t id,copy,copy1,copy2;
392   static Float_t hits[9];
393   static Int_t vol[3];
394   static Float_t de;
395   static TLorentzVector pos;
396   static TLorentzVector mom;
397   static Int_t iPart;
398   
399   TClonesArray &lhits = *fHits;
400   if(!gMC->IsTrackAlive()) return; // particle has disappeared
401
402   Float_t charge = gMC->TrackCharge();
403   if(TMath::Abs(charge)<=0.) return; //take only charged particles
404
405   //  printf(" in StepManeger \n");
406   id=gMC->CurrentVolID(copy);
407   //((TGeant3*)gMC)->Gpcxyz();
408   
409   // Check the sensetive volume
410   if(id==fIdSens1||id==fIdSens2||id==fIdSens3||id==fIdSens4||id==fIdSens5) {
411
412     if(gMC->IsTrackEntering())  {
413       vol[2]=copy-1;
414       gMC->CurrentVolOffID(1,copy1);
415       vol[1]=copy1;
416       gMC->CurrentVolOffID(2,copy2);
417       vol[0]=copy2;
418
419       gMC->TrackPosition(pos);
420       hits[0]=pos[0];
421       hits[1]=pos[1];
422       hits[2]=pos[2];
423       
424       gMC->TrackMomentum(mom);
425       hits[3]=mom[0];
426       hits[4]=mom[1];
427       hits[5]=mom[2];
428       
429       iPart= gMC->TrackPid();
430       Int_t partId=gMC->IdFromPDG(iPart);
431       hits[7]=partId;
432       hits[8]=1e9*gMC->TrackTime();
433       de=0.;
434     }
435     Float_t edep = 1000 * gMC->Edep();
436     Float_t p    = mom.P();
437     TParticlePDG* pdg = TDatabasePDG::Instance()->GetParticle(iPart);
438     Float_t mass = pdg ? pdg->Mass() : 1;
439     if (edep > 1 && p/mass > 1) {
440       TArrayI procs;
441       gMC->StepProcesses(procs);
442       TString processes;
443       for (Int_t ip = 0; ip < procs.fN; ip++) {
444         if (ip != 0) processes.Append(",");
445         processes.Append(TMCProcessName[procs.fArray[ip]]);
446       }
447       TString what;
448       if (gMC->IsTrackEntering())    what.Append("entering ");
449       if (gMC->IsTrackExiting())     what.Append("exiting ");
450       if (gMC->IsTrackInside())      what.Append("inside ");
451       if (gMC->IsTrackDisappeared()) what.Append("disappeared ");
452       if (gMC->IsTrackStop())        what.Append("stopped ");
453       if (gMC->IsNewTrack())         what.Append("new ");
454       if (gMC->IsTrackAlive())       what.Append("alive ");
455       if (gMC->IsTrackOut())         what.Append("out ");
456       Int_t trackno = gAlice->GetMCApp()->GetCurrentTrackNumber();
457       Int_t mother  = gAlice->GetMCApp()->GetPrimary(trackno);
458       Warning("StepManager", "Track # %5d deposits a lot of energy\n" 
459               "  Volume:    %s\n" 
460               "  Momentum:  (%7.4f,%7.4f,%7.4f)\n"
461               "  PDG:       %d (%s)\n" 
462               "  Edep:      %-14.7f keV (mother %d)\n"
463               "  p/m:       %-7.4f/%-7.4f = %-14.7f\n"
464               "  Processes: %s\n"
465               "  What:      %s\n",
466               trackno, gMC->CurrentVolPath(), mom.X(), mom.Y(), mom.Z(),
467               iPart, (pdg ? pdg->GetName() : ""), edep, mother, mom.P(), mass, 
468               mom.P()/mass, processes.Data(), what.Data());
469     }
470     
471     if(gMC->IsTrackInside()) de=de+1000.*gMC->Edep();
472
473     if(gMC->IsTrackExiting() ||gMC->IsTrackDisappeared()|| gMC->IsTrackStop()){
474       hits[6]=de+1000.*gMC->Edep();
475       UShort_t detector = vol[0] / 2 + 1;
476       Char_t   ring     = (vol[0] % 2) == 0 ? 'I' : 'O';
477       UShort_t sector   = vol[1];
478       UShort_t strip    = vol[2];
479       gMC->TrackPosition(pos);
480       TVector3 cur(pos.Vect());
481       TVector3 old(hits[0], hits[1], hits[2]);
482       cur -= old;
483       Double_t len = cur.Mag();
484       AliFMDHit* h = new(lhits[fNhits++]) 
485         AliFMDHit(fIshunt,
486                   gAlice->GetMCApp()->GetCurrentTrackNumber(),
487                   detector, ring, sector, strip, 
488                   hits[0], hits[1], hits[2], hits[3], hits[4], hits[5], 
489                   hits[6], iPart, hits[8], len, 
490                   gMC->IsTrackDisappeared()||gMC->IsTrackStop());
491       if (hits[6] > 1 && p/mass > 1) fBad->Add(h);
492     } // IsTrackExiting()
493   }
494 }
495 //--------------------------------------------------------------------------
496 void 
497 AliFMDAlla::Response(Float_t Edep)
498 {
499   Float_t I=1.664*0.04*2.33/22400; // = 0.69e-6;
500   Float_t chargeOnly=Edep/I;
501   //Add noise ~500electrons
502   Int_t charge=500;
503   if (Edep>0)
504     charge=Int_t(gRandom->Gaus(chargeOnly,500));        
505 }   
506
507 //____________________________________________________________________
508 AliDigitizer* 
509 AliFMDAlla::CreateDigitizer(AliRunDigitizer* manager) const
510 {
511   // Create a digitizer object 
512   AliFMDDigitizerAlla* digitizer = new AliFMDDigitizerAlla(manager);
513   return digitizer;
514 }
515 //--------------------------------------------------------------------------
516 //
517 // EOF
518 //
519
520
521
522
523