]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALv0.cxx
Simplification of Makefile and some small corrections
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALv0.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 /* $Id$ */
16
17 //_________________________________________________________________________
18 // Implementation version v0 of EMCAL Manager class 
19 // An object of this class does not produce hits nor digits
20 // It is the one to use if you do not want to produce outputs in TREEH or TREED
21 // This class places a Geometry of the EMCAL in the ALICE Detector as defined in AliEMCALGeometry.cxx                 
22 //*-- Author: Yves Schutz (SUBATECH)
23 //*-- and   : Sahal Yacoob (LBL / UCT)
24 //          : Aleksei Pavlinov (WSU)     SHASHLYK
25
26 // --- ROOT system ---
27
28 #include "TNode.h"
29 #include "TBRIK.h"
30 #include "TTRD1.h"
31 #include "TTRD2.h"
32 #include "TTRAP.h"
33 #include "TPGON.h"
34 #include "TTUBS.h"
35 #include "TGeometry.h"
36 #include "TVirtualMC.h"
37 #include "TArrayI.h"
38 #include "TROOT.h"
39 #include "TArrayF.h"
40 #include "TList.h"
41 #include "TVector2.h"
42
43 #include "AliEMCALShishKebabModule.h"
44 #include "AliEMCALShishKebabTrd1Module.h"
45
46 // --- Standard library ---
47
48 //#include <stdio.h>
49
50 // --- AliRoot header files ---
51
52 #include "AliEMCALv0.h"
53 #include "AliEMCALGeometry.h"
54 #include "AliRun.h"
55 #include "AliLog.h"
56
57 ClassImp(AliEMCALv0)
58
59 TArrayF ENVELOP1; // now global for BuildGeometry() - 30-aug-2004
60 int idAIR=1599, idPB = 1600, idSC = 1601, idSTEEL = 1603; // global
61 Int_t *idtmed=0, idrotm=0;
62 double sampleWidth=0.;
63 double parEMOD[5], smodPar0=0., smodPar1=0., smodPar2=0.;
64
65 //______________________________________________________________________
66 AliEMCALv0::AliEMCALv0(const char *name, const char *title):
67   AliEMCAL(name,title)
68 {
69   // ctor : title is used to identify the layout
70   AliEMCALGeometry *geom = GetGeometry() ; 
71   //geom->CreateListOfTrd1Modules(); 
72   fShishKebabModules = geom->GetShishKebabTrd1Modules(); 
73 }
74
75 //______________________________________________________________________
76 void AliEMCALv0::BuildGeometry()
77 {
78     // Display Geometry for display.C
79
80     const Int_t kColorArm1   = kBlue ;
81
82     AliEMCALGeometry * geom = GetGeometry();
83
84     TString gn(geom->GetName());
85     gn.ToUpper(); 
86
87     // Define the shape of the Calorimeter 
88     TNode * top = gAlice->GetGeometry()->GetNode("alice") ; // See AliceGeom/Nodes
89     TNode * envelopNode = 0;
90     char *envn = "Envelop1";
91     if(!gn.Contains("SHISH") || gn.Contains("TRD2")){
92       new TTUBS(envn, "Tubs that contains arm 1", "void", 
93               geom->GetEnvelop(0) -10, // rmin 
94               geom->GetEnvelop(1) +40 ,// rmax
95               geom->GetEnvelop(2)/2.0, // half length in Z
96               geom->GetArm1PhiMin(),   // minimum phi angle
97               geom->GetArm1PhiMax()    // maximum phi angle
98         );
99       top->cd();
100       envelopNode = new TNode(envn, "Arm1 Envelop", "Envelop1", 0., 0., 0., "") ;
101     } else {
102       if(gn.Contains("WSUC")) {
103         envelopNode = BuildGeometryOfWSUC();
104       } else { // Shish-kebab now for compact, twist and TRD1 cases (ALIC)
105         envn="Envelop2";
106         TPGON *pgon = new TPGON(envn, "PGON that contains arm 1", "void", 
107         geom->GetArm1PhiMin(),geom->GetArm1PhiMax()-geom->GetArm1PhiMin(),geom->GetNPhiSuperModule(), 2);
108       // define section
109         pgon->DefineSection(0, ENVELOP1[4],  ENVELOP1[5], ENVELOP1[6]);
110         pgon->DefineSection(1, ENVELOP1[7],  ENVELOP1[5], ENVELOP1[6]);
111         top->cd();
112         envelopNode = new TNode(envn, "Arm1 Envelop2", envn, 0., 0., 0., "") ;
113       }
114     }
115                                 
116     envelopNode->SetLineColor(kColorArm1) ;
117     fNodes->Add(envelopNode);
118 }
119
120 TNode *AliEMCALv0::BuildGeometryOfWSUC()
121 { // June 8, 2005; see directory geant3/TGeant3/G3toRoot.cxx
122   // enum EColor { kWhite, kBlack, kRed, kGreen, kBlue, kYellow, kMagenta, kCyan } - see $ROOTSYS/include/Gtypes.h
123    AliEMCALGeometry * g = GetGeometry(); 
124    TNode * top = gAlice->GetGeometry()->GetNode("alice") ; // See AliceGeom/Nodes
125    top->cd();
126
127    TNode *envelopNode = 0;
128    char *name = "";
129    /*
130     name = "WSUC";
131    new TBRIK(name, "WSUC(XEN1 in Geant)","void",ENVELOP1[0],ENVELOP1[1],ENVELOP1[2]);
132     envelopNode = new TNode(name, "envelope for WSUC", name, 0., 0., 0., "");
133    envelopNode->SetVisibility(0);
134    */
135
136    TNode *emod=0, *scmx=0;
137    name = "SMOD"; // super module
138    new TBRIK(name, "SMOD(SMOD in Geant)","void", smodPar0,smodPar1,smodPar2);
139    if(envelopNode) envelopNode->cd();
140    TNode *smod = new TNode(name, "SMOD", name, 0., 0., 0., "");   
141    smod->SetLineColor(kBlue) ;
142    if(envelopNode==0) envelopNode = smod;
143
144    name = "EMOD"; // see CreateEMOD 
145    TTRD1 *EMOD = new TTRD1(name, "EMOD(EMOD in Geant)","void", float(parEMOD[0]),
146    float(parEMOD[1]),float(parEMOD[2]),float(parEMOD[3]));
147
148    // SCMX = EMOD/4 for simplicity of drawing
149    name = "SCMX";
150    Float_t dz=0.,theta=0.,phi=0.,h1=0.,bl1=0.,tl1=0.,alpha1=0.,h2=0.,bl2=0.,tl2=0.,alpha2=0.;
151    h1     = EMOD->GetDy()/2.;
152    bl1    = EMOD->GetDx()/2.;
153    tl1    = bl1;
154    alpha1 = 0.;
155    h2     = EMOD->GetDy()/2.;
156    bl2    = EMOD->GetDx2()/2.;
157    tl2    = bl2;
158    alpha2 = 0.;
159
160    dz       = EMOD->GetDz();
161    double dr = TMath::Sqrt((h2-h1)*(h2-h1)+(bl2-bl1)*(bl2-bl1));
162    theta    = TMath::ATan2(dr,2.*dz) * TMath::RadToDeg(); 
163    phi      = 180.;
164
165    TTRAP *SCMX = new TTRAP(name, "SCMX(SCMX as in Geant)","void",
166              dz,theta,phi, h1,bl1,tl1,alpha1, h2,bl2,tl2,alpha2);
167    //   SCMX->Dump();
168    Float_t xShiftSCMX = (EMOD->GetDx() + EMOD->GetDx2())/4.;
169    Float_t yShiftSCMX = EMOD->GetDy()/2.;
170    printf(" xShiftSCMX %7.4f yShiftSCMX %7.4f \n",xShiftSCMX,yShiftSCMX);
171
172    name = "EMOD"; // see CreateEMOD 
173    smod->cd();
174    
175    AliEMCALShishKebabTrd1Module *mod=0;
176    Double_t  angle=90., xpos=0.,ypos=0.,zpos=0., xposSCMX=0.,yposSCMX=0.,zposSCMX=0.;
177    char rtmn[100], rtmt[100];
178    TRotMatrix *rtm=0, *rtmSCMX=0;
179    int numEmod=0;
180    for(int iz=0; iz<g->GetNZ(); iz++) {
181      mod   = (AliEMCALShishKebabTrd1Module*)fShishKebabModules->At(iz);
182      zpos = mod->GetPosZ()      - smodPar2;
183      ypos = mod->GetPosXfromR() - smodPar1;
184
185      angle = mod->GetThetaInDegree();
186      sprintf(rtmn,"rmEmod%5.1f",angle);
187      sprintf(rtmt,"rotation matrix for EMOD, iz=%i, angle = %6.3f",iz, angle);
188      if(iz==0) rtm = new TRotMatrix(rtmn, rtmt,0.,0., 90.,0., 90.,90.); // z'(x); y'(y); x'(z)
189      else      rtm = new TRotMatrix(rtmn, rtmt,90.-angle,270., 90.0,0.0, angle,90.);
190
191      TGeometry *tg = gAlice->GetGeometry();
192      for(int ix=0; ix<g->GetNPhi(); ix++) { // flat in phi
193        xpos = g->GetPhiModuleSize()*(2*ix+1 - g->GetNPhi())/2.;
194        sprintf(rtmt,"EMOD, iz %i, ix %i, angle %6.3f",iz,ix, angle);
195        TString namNode=name;
196        namNode += numEmod++;
197        smod->cd();
198        emod = new TNode(namNode.Data(), rtmt, (TShape*)EMOD, xpos,ypos,zpos,rtm);   
199        //       emod->SetLineColor(kGreen) ;
200        emod->SetVisibility(0); // SCMX will bi visible 
201        if(SCMX) { // 4(2x2) sensetive volume inside EMOD
202          emod->cd();
203          zposSCMX = 0.;
204          for(int jy=0; jy<2; jy++){ // division on y
205            yposSCMX = yShiftSCMX *(2*jy - 1); 
206            for(int jx=0; jx<2; jx++){ // division on x
207              Double_t theta1=90.,phi1=0., theta2=90.,phi2=90., theta3=0.,phi3=0 ; 
208              xposSCMX = xShiftSCMX *(2*jx - 1);              
209              namNode = "SCMX";
210              namNode += jy;
211              namNode += jx;
212              sprintf(rtmn,"rm%s",namNode.Data());
213              sprintf(rtmt,"rotation matrix for %s inside EMOD",namNode.Data());
214              rtmSCMX = tg->GetRotMatrix(rtmn);
215              if(jx == 1) {
216                phi1 = 180.; // x' = -x
217                phi2 = 270.; // y' = -y
218              }
219              if(rtmSCMX == 0) rtmSCMX = new TRotMatrix(rtmn,rtmt, theta1,phi1, theta2,phi2, theta3,phi3);
220              sprintf(rtmt,"%s inside %s", namNode.Data(), emod->GetName());
221              scmx = new TNode(namNode.Data(), rtmt, (TShape*)SCMX, xposSCMX,yposSCMX,zposSCMX,rtmSCMX);   
222              scmx->SetLineColor(kMagenta);
223            }
224          }
225        }
226      }
227    }
228    //   emod->Draw(); // for testing
229
230    return envelopNode;
231 }
232
233 //______________________________________________________________________
234 void AliEMCALv0::CreateGeometry()
235 {
236   // Create the EMCAL geometry for Geant
237   // Geometry of a tower
238   //|-----------------------------------------------------| XEN1
239   //| |                                                 | |
240   //| |    Al thickness = GetAlFrontThickness()         | |
241   //| |                                                 | |
242   //| |                                                 | |
243   //| |                                                 | |
244   //|  -------------------------------------------------  |
245   //| |    Air Gap = GetGap2Active()                    | |
246   //| |                                                 | |
247   //|  -------------------------------------------------  |
248   //| |    XU1 : XPST (ECAL e = GetECScintThick() )     | |
249   //|  -------------------------------------------------  |
250   //| |    XU1 : XPBX (ECAL e = GetECPbRadThick() )     | |
251   //|  -------------------------------------------------  |
252   //| |    XU1 : XPST (ECAL e = GetECScintThick()       | |
253   //|  -------------------------------------------------  |
254   //| |    XU1 : XPBX (ECAL e = GetECPbRadThick() )     | |
255   //|  -------------------------------------------------  |
256   //|    etc ..... GetNECLayers() - 1 times               |
257   //|  -------------------------------------------------  |
258   //| |    XUNLayer : XPST (ECAL e = GetECScintThick()  | |
259   //|  -------------------------------------------------  |
260
261     AliEMCALGeometry * geom = GetGeometry() ; 
262     TString gn(geom->GetName());
263     gn.ToUpper(); 
264
265     if(!(geom->IsInitialized())){
266         Error("CreateGeometry","EMCAL Geometry class has not been set up.");
267     } // end if
268
269     // Get pointer to the array containing media indices
270     idtmed = fIdtmed->GetArray() - 1599 ;
271
272     idrotm = 1;
273     //  gMC->Matrix(nmat, theta1, phi1, theta2, phi2, theta3, phi3) - see AliModule
274     AliMatrix(idrotm, 90.0, 0., 90.0, 90.0, 0.0, 0.0) ; 
275
276     // Create the EMCAL Mother Volume (a polygone) within which to place the Detector and named XEN1 
277
278     Float_t envelopA[10];
279     if(gn.Contains("TRD2")) { // TUBS
280        envelopA[0] = geom->GetEnvelop(0) - 10.; // rmin 
281        envelopA[1] = geom->GetEnvelop(1) + 12.; // rmax
282        //       envelopA[2] = geom->ZFromEtaR(geom->GetEnvelop(1), geom->GetArm1EtaMin());
283        envelopA[2] = 390.; // 6-feb-05
284        envelopA[3] = geom->GetArm1PhiMin();
285        envelopA[4] = geom->GetArm1PhiMax();
286        gMC->Gsvolu("XEN1", "TUBS", idtmed[1599], envelopA, 5) ;   // Tubs filled with air 
287        ENVELOP1.Set(5, envelopA);
288     // Position the EMCAL Mother Volume (XEN1) in Alice (ALIC)  
289        gMC->Gspos("XEN1", 1, "ALIC", 0.0, 0.0, 0.0, idrotm, "ONLY") ;
290     } else if(gn.Contains("TRD1") && gn.Contains("WSUC") ) { // TRD1 for WSUC facility
291       // 17-may-05 - just BOX
292       envelopA[0] = 26;
293       envelopA[1] = 15;
294       envelopA[2] = 30;
295       gMC->Gsvolu("XEN1", "BOX", idtmed[idSC], envelopA, 3) ;
296       ENVELOP1.Set(3);
297       for(int i=0; i<3; i++) ENVELOP1[i] = envelopA[i]; // 23-may-05  
298       // Position the EMCAL Mother Volume (XEN1) in WSUC  
299       gMC->Gspos("XEN1", 1, "WSUC", 0.0, 0.0, 0.0, idrotm, "ONLY") ;
300     } else { 
301       envelopA[0] = geom->GetArm1PhiMin();                         // minimum phi angle
302       envelopA[1] = geom->GetArm1PhiMax() - geom->GetArm1PhiMin(); // angular range in phi
303       envelopA[2] = geom->GetNPhi();                               // number of sections in phi
304       envelopA[3] = 2;                                             // 2 z coordinates
305       envelopA[4] = geom->ZFromEtaR(geom->GetEnvelop(1),
306                                    geom->GetArm1EtaMin());       // z coordinate 1
307     //add some padding for mother volume
308       envelopA[5] = geom->GetEnvelop(0) ;                          // rmin at z1
309       envelopA[6] = geom->GetEnvelop(1) ;                          // rmax at z1
310       envelopA[7] = geom->ZFromEtaR(geom->GetEnvelop(1),
311                                   geom->GetArm1EtaMax());        // z coordinate 2
312       envelopA[8] = envelopA[5] ;                                  // radii are the same.
313       envelopA[9] = envelopA[6] ;                                  // radii are the same.
314
315       if(gn.Contains("SHISH")) envelopA[2] = geom->GetNPhiSuperModule();
316
317       gMC->Gsvolu("XEN1", "PGON", idtmed[1599], envelopA, 10) ;   // Polygone filled with air 
318       ENVELOP1.Set(10, envelopA);
319       if (gDebug==2) {
320         printf("CreateGeometry: XEN1 = %f, %f\n", envelopA[5], envelopA[6]); 
321         printf("CreateGeometry: XU0 = %f, %f\n", envelopA[5], envelopA[6]); 
322       }
323     // Position the EMCAL Mother Volume (XEN1) in Alice (ALIC)  
324       gMC->Gspos(geom->GetNameOfEMCALEnvelope(), 1, "ALIC", 0.0, 0.0, 0.0, idrotm, "ONLY") ;
325     }
326
327     if(gn.Contains("SHISH")){
328       // COMPACT, TWIST, TRD2 or TRD1
329       CreateShishKebabGeometry();
330     }
331 }
332
333 //______________________________________________________________________
334 void AliEMCALv0::Init(void)
335 {
336     // Just prints an information message
337   
338   if(AliLog::GetGlobalDebugLevel()>0) { 
339     TString message("\n") ; 
340     message += "*****************************************\n" ;
341     
342     // Here the EMCAL initialisation code (if any!)
343     
344     AliEMCALGeometry * geom = GetGeometry() ; 
345     
346     if (geom!=0) {   
347       message += "AliEMCAL " ; 
348       message += Version() ; 
349       message += "EMCAL geometry initialized for " ; 
350       message += geom->GetName()  ;
351     }
352     else {
353       message += "AliEMCAL " ; 
354       message += Version() ;  
355       message += "EMCAL geometry initialization failed !" ; 
356     }
357     message += "\n*****************************************" ;
358     printf(message.Data() ) ; 
359   }
360 }
361
362 // 24-aug-04 by PAI
363 void AliEMCALv0::CreateShishKebabGeometry()
364 {  // TWIST, TRD1 and TRD2 
365   AliEMCALGeometry * g = GetGeometry(); 
366   TString gn(g->GetName()); gn.ToUpper(); 
367   // see AliModule::fIdtmed
368   //  idtmed = fIdtmed->GetArray() - 1599 ; // see AliEMCAL::::CreateMaterials()
369   //  int idAIR=1599, idPB = 1600, idSC = 1601, idSTEEL = 1603;
370   //  idAL = 1602;
371   Double_t par[10], xpos=0., ypos=0., zpos=0.;
372
373   CreateSmod(g->GetNameOfEMCALEnvelope());
374
375   CreateEmod("SMOD","EMOD"); // 18-may-05
376
377   if(gn.Contains("110DEG")) CreateEmod("SM10","EMOD"); // 12-oct-05
378
379   // Sensitive SC  (2x2 tiles)
380   double parSCM0[5], *dummy = 0, parTRAP[11];
381   Double_t trd1Angle = g->GetTrd1Angle()*TMath::DegToRad(), tanTmp = TMath::Tan(trd1Angle/2.);
382    if(!gn.Contains("TRD")) { // standard module
383     par[0] = (g->GetECPbRadThick()+g->GetECScintThick())*g->GetNECLayers()/2.;
384     par[1] = par[2] = g->GetPhiTileSize();   // Symetric case
385     gMC->Gsvolu("SCM0", "BOX", idtmed[idSC], par, 3); // 2x2 tiles
386     gMC->Gspos("SCM0", 1, "EMOD", 0., 0., 0., 0, "ONLY") ;
387     // Division to tile size
388     gMC->Gsdvn("SCM1","SCM0", g->GetNPHIdiv(), 2); // y-axis
389     gMC->Gsdvn("SCM2","SCM1", g->GetNETAdiv(), 3); // z-axis
390   // put LED to the SCM2 
391     par[0] = g->GetECPbRadThick()/2.;
392     par[1] = par[2] = g->GetPhiTileSize()/2.; // Symetric case
393     gMC->Gsvolu("PBTI", "BOX", idtmed[idPB], par, 3);
394
395     printf(" Pb tiles \n");
396     int nr=0;
397     ypos = zpos = 0.0;
398     xpos = -sampleWidth*g->GetNECLayers()/2. + g->GetECPbRadThick()/2.;
399     for(int ix=0; ix<g->GetNECLayers(); ix++){
400       gMC->Gspos("PBTI", ++nr, "SCM2", xpos, ypos, zpos, 0, "ONLY") ;
401     //    printf(" %i xpos %f \n", ix+1, xpos);
402       xpos += sampleWidth;
403     } 
404     printf(" Number of Pb tiles in SCM2 %i \n", nr);
405   } else if(gn.Contains("TRD1")) { // TRD1 - 30-sep-04
406     if(gn.Contains("MAY05")){
407       Double_t dzTmp = g->GetFrontSteelStrip()+g->GetPassiveScintThick();
408       parSCM0[0] = parEMOD[0] + tanTmp*dzTmp; // dx1
409       parSCM0[1] = parEMOD[1];                // dx2
410       parSCM0[2] = parEMOD[2];                // dy
411       for(int i=0; i<3; i++) parSCM0[i] -= g->GetLateralSteelStrip();
412       parSCM0[3] = parEMOD[3] - dzTmp/2.; // dz
413
414       gMC->Gsvolu("SCM0", "TRD1", idtmed[idAIR], parSCM0, 4);
415       gMC->Gspos("SCM0", 1, "EMOD", 0., 0., dzTmp/2., 0, "ONLY") ;
416     } else { // before MAY 2005
417       double wallThickness = g->GetPhiModuleSize()/2. -  g->GetPhiTileSize(); // Need check
418       for(int i=0; i<3; i++) parSCM0[i] = parEMOD[i] - wallThickness;
419       parSCM0[3] = parEMOD[3];
420       gMC->Gsvolu("SCM0", "TRD1", idtmed[idAIR], parSCM0, 4);
421       gMC->Gspos("SCM0", 1, "EMOD", 0., 0., 0., 0, "ONLY") ;
422     }
423
424     if(g->GetNPHIdiv()==2 && g->GetNETAdiv()==2) {
425     // Division to tile size - 1-oct-04
426       printf("<I> Divide SCM0 on y-axis %i\n", g->GetNETAdiv());
427       gMC->Gsdvn("SCMY","SCM0", g->GetNETAdiv(), 2); // y-axis
428     // Trapesoid 2x2
429       parTRAP[0] = parSCM0[3];    // dz
430       parTRAP[1] = TMath::ATan2((parSCM0[1]-parSCM0[0])/2.,2.*parSCM0[3])*180./TMath::Pi(); // theta
431       parTRAP[2] = 0.;           // phi
432       // bottom
433       parTRAP[3] = parSCM0[2]/2.; // H1
434       parTRAP[4] = parSCM0[0]/2.; // BL1
435       parTRAP[5] = parTRAP[4];    // TL1
436       parTRAP[6] = 0.0;           // ALP1
437       // top
438       parTRAP[7] = parSCM0[2]/2.; // H2
439       parTRAP[8] = parSCM0[1]/2.; // BL2
440       parTRAP[9] = parTRAP[8];    // TL2
441       parTRAP[10]= 0.0;           // ALP2
442       printf(" ** TRAP ** \n");
443       for(int i=0; i<11; i++) AliDebug(3, Form(" par[%2.2i] %9.4f\n", i, parTRAP[i]));
444
445       gMC->Gsvolu("SCMX", "TRAP", idtmed[idSC], parTRAP, 11);
446       xpos = +(parSCM0[1]+parSCM0[0])/4.;
447       gMC->Gspos("SCMX", 1, "SCMY", xpos, 0.0, 0.0, 0, "ONLY") ;
448
449       // Using rotation because SCMX should be the same due to Pb tiles
450       xpos = -xpos; 
451       AliMatrix(idrotm, 90.0,180., 90.0, 270.0, 0.0,0.0) ;
452       gMC->Gspos("SCMX", 2, "SCMY", xpos, 0.0, 0.0, idrotm, "ONLY");
453     // put LED to the SCM0 
454       AliEMCALShishKebabTrd1Module *mod = (AliEMCALShishKebabTrd1Module*)fShishKebabModules->At(0);
455       gMC->Gsvolu("PBTI", "BOX", idtmed[idPB], dummy, 0);
456
457       par[1] = parSCM0[2]/2;            // y 
458       par[2] = g->GetECPbRadThick()/2.; // z
459
460       int nr=0;
461       ypos = 0.0; 
462       zpos = -sampleWidth*g->GetNECLayers()/2. + g->GetECPbRadThick()/2.;
463       double xCenterSCMX =  (parTRAP[4] +  parTRAP[8])/2.;
464       if(!gn.Contains("NOPB")) { // for testing - 11-jul-05
465         printf(" Pb tiles \n");
466         for(int iz=0; iz<g->GetNECLayers(); iz++){
467           par[0] = (parSCM0[0] + mod->GetTanBetta()*sampleWidth*iz)/2.;
468           xpos   = par[0] - xCenterSCMX;
469           gMC->Gsposp("PBTI", ++nr, "SCMX", xpos, ypos, zpos, 0, "ONLY", par, 3) ;
470           AliDebug(3,Form(" %i xpos %f zpos %f par[0] %f \n", iz+1, xpos, zpos, par[0]));
471           zpos += sampleWidth;
472         } 
473         printf(" Number of Pb tiles in SCMX %i \n", nr);
474       }
475     } else if(g->GetNPHIdiv()==3 && g->GetNETAdiv()==3) {
476       Trd1Tower3X3(parSCM0);
477     } else if(g->GetNPHIdiv()==4 && g->GetNETAdiv()==4) {
478       Trd1Tower4X4();
479     }
480   } else if(gn.Contains("TRD2")) {    // TRD2 - 14-jan-05
481     //    Scm0InTrd2(g, parEMOD, parSCM0); // First dessin 
482     PbmoInTrd2(g, parEMOD, parSCM0); // Second dessin 
483   }
484 }
485
486 void AliEMCALv0::CreateSmod(const char* mother)
487 { // 18-may-05; mother="XEN1"; 
488   // child="SMOD" from first to 10th, "SM10" (11th and 12th) (TRD1 case)
489   // child="SMON" and "SMOP"("TRD2" case)
490   AliEMCALGeometry * g = GetGeometry(); 
491   TString gn(g->GetName()); gn.ToUpper();
492
493   Double_t par[3], parTubs[5], xpos=0., ypos=0., zpos=0., rpos=0., dphi=0., phi=0.0, phiRad=0.;
494   Double_t par1C = 0.;
495   //  ===== define Super Module from air - 14x30 module ==== ;
496   sampleWidth = double(g->GetECPbRadThick()+g->GetECScintThick());
497   printf("\n ## Super Module | sampleWidth %5.3f ## %s \n", sampleWidth, gn.Data());
498   par[0] = g->GetShellThickness()/2.;
499   par[1] = g->GetPhiModuleSize()*g->GetNPhi()/2.; 
500   par[2] = g->GetEtaModuleSize()*15.; 
501   idrotm=0;
502   int nphism = g->GetNumberOfSuperModules()/2; // 20-may-05
503   if(nphism>0) {
504     dphi = (g->GetArm1PhiMax() - g->GetArm1PhiMin())/nphism;
505     //    if(gn.Contains("110DEG")) dphi = (g->GetArm1PhiMax() - g->GetArm1PhiMin())/(nphism-1);
506     rpos = (g->GetEnvelop(0) + g->GetEnvelop(1))/2.;
507     printf(" rpos %8.2f : dphi %6.1f degree \n", rpos, dphi);
508   }
509
510   if (gn.Contains("TRD2")) { // tubs - 27-jan-05
511     parTubs[0] = g->GetTubsR();                       // rmin
512     parTubs[1] = parTubs[0] + g->GetShellThickness(); // rmax ?? 
513     parTubs[2] = 380./2.;                             // DZ half length in z; 11-oct-04 - for 26 division
514     parTubs[3] = -dphi/2.;                            // PHI1 starting angle of the segment;
515     parTubs[4] = +dphi/2.;                            // PHI2 ending angle of the segment;
516
517     gMC->Gsvolu("SMOP", "TUBS", idtmed[idAIR], parTubs, 5); // pozitive Z
518     gMC->Gsvolu("SMON", "TUBS", idtmed[idAIR], parTubs, 5); // negative Z
519
520     printf(" SMOP,N ** TUBS **\n"); 
521     printf("tmed %i | Rmin %7.2f Rmax %7.2f dz %7.2f phi1,2 (%7.2f,%7.2f)\n", 
522     idtmed[idAIR], parTubs[0],parTubs[1],parTubs[2], parTubs[3],parTubs[4]);
523     // have to add 1 cm plastic before EMOD - time solution 
524   } else if(gn.Contains("WSUC")) {
525     par[0] = g->GetPhiModuleSize()*g->GetNPhi()/2.; 
526     par[1] = g->GetShellThickness()/2.;
527     par[2] = g->GetEtaModuleSize()*g->GetNZ()/2. + 5; 
528
529     gMC->Gsvolu("SMOD", "BOX", idtmed[idAIR], par, 3);
530
531     printf("SMOD in WSUC : tmed %i | dx %7.2f dy %7.2f dz %7.2f (SMOD, BOX)\n", 
532     idtmed[idAIR], par[0],par[1],par[2]);
533     smodPar0 = par[0]; 
534     smodPar1 = par[1];
535     smodPar2 = par[2];
536     nphism   =  g->GetNumberOfSuperModules();
537   } else {
538     if     (gn.Contains("TWIST")) {
539       par[2] += 0.4;      // for 27 division
540     } else if(gn.Contains("TRD")) {
541       par[2]  = 350./2.; // 11-oct-04 - for 26 division
542       if(gn.Contains("TRD1")) {
543         printf(" par[0] %7.2f (old) \n",  par[0]);  
544         Float_t *parSM = g->GetSuperModulesPars(); 
545         for(int i=0; i<3; i++) par[i] = parSM[i];
546       }
547     }
548     gMC->Gsvolu("SMOD", "BOX", idtmed[idAIR], par, 3);
549     printf("tmed %i | dx %7.2f dy %7.2f dz %7.2f (SMOD, BOX)\n", 
550     idtmed[idAIR], par[0],par[1],par[2]);
551     smodPar0 = par[0]; 
552     smodPar2 = par[2];
553     if(gn.Contains("110DEG")) { // 12-oct-05
554       par1C = par[1];
555       par[1] /= 2.;
556       gMC->Gsvolu("SM10", "BOX", idtmed[idAIR], par, 3);
557       printf(" Super module with name \"SM10\" was created too par[1] = %f\n", par[1]);
558       par[1] = par1C;
559     }
560   // Steel plate
561     if(g->GetSteelFrontThickness() > 0.0) { // 28-mar-05
562       par[0] = g->GetSteelFrontThickness()/2.;
563       gMC->Gsvolu("STPL", "BOX", idtmed[idSTEEL], par, 3);
564       printf("tmed %i | dx %7.2f dy %7.2f dz %7.2f (STPL) \n", idtmed[idSTEEL], par[0],par[1],par[2]);
565       xpos = -(g->GetShellThickness() - g->GetSteelFrontThickness())/2.;
566       gMC->Gspos("STPL", 1, "SMOD", xpos, 0.0, 0.0, 0, "ONLY") ;
567     }
568   }
569
570   int nr=0, nrsmod=0, i0=0;
571   if(gn.Contains("TEST")) {nphism = 1;} // just only 2 super modules;
572
573   // Turn whole super module
574   int TurnSupMod = 1; // should be ONE; for testing =0
575   for(int i=i0; i<nphism; i++) {
576     if (gn.Contains("TRD2")) {      // tubs - 27-jan-05
577       if(i==i0) {
578         printf("** TRD2 ** ");
579         if(TurnSupMod==1) printf(" No 3 degree rotation !!! ");
580         printf("\n");
581       }
582       Double_t phic=0., phicRad=0.; // phi angle of arc center
583       phic    = g->GetArm1PhiMin() + dphi*(2*i+1)/2.; //
584       phicRad = phic*TMath::DegToRad();
585       phi     = phic - g->GetTubsTurnAngle();
586       phiRad  = phi*TMath::DegToRad();
587       if(TurnSupMod==1) {
588         TVector2  vc;     // position of super module center
589         vc.SetMagPhi(parTubs[0], phicRad);
590         TVector2  vcTurn; // position of super module center with turn
591         vcTurn.SetMagPhi(parTubs[0], phiRad);
592         TVector2 vcShift = vc - vcTurn;
593         phic = phi;
594
595         xpos = vcShift.X();
596         ypos = vcShift.Y();
597       } else { // 1-mar-05 ; just for testing - no turn od SMOD; looks good
598         xpos = ypos = 0.0;
599       }
600       zpos = parTubs[2];
601       AliMatrix(idrotm, 90.0, phic, 90.0, 90.0+phic, 0.0, 0.0);
602
603       gMC->Gspos("SMOP", ++nr, mother, xpos, ypos, zpos, idrotm, "ONLY") ;
604       printf("SMOP %2i | %2i idrotm %3i phi %6.1f(%5.3f) xpos %7.2f ypos %7.2f zpos %7.2f \n", 
605       i, nr, idrotm, phic, phicRad, xpos, ypos, zpos);
606       printf(" phiy(90+phic)  %6.1f \n", 90. + phic);
607
608       if(!gn.Contains("TEST1") && g->GetNumberOfSuperModules() > 1){
609         //        double  phiy = 90. + phic + 180.;
610         //        if(phiy>=360.) phiy -= 360.;
611         //        printf(" phiy  %6.1f \n", phiy);
612         //        AliMatrix(idrotm, 90.0, phic, 90.0, phiy, 180.0, 0.0);
613         gMC->Gspos("SMON", nr, mother, xpos, ypos, -zpos, idrotm, "ONLY") ;
614         printf("SMON %2i | %2i idrotm %3i phi %6.1f(%5.3f) xpos %7.2f ypos %7.2f zpos %7.2f \n", 
615         i, nr, idrotm, phic, phicRad, xpos, ypos, -zpos);
616       }
617     } else if(gn.Contains("WSUC")) {
618       xpos = ypos = zpos = 0.0;
619       idrotm = 0;
620       gMC->Gspos("SMOD", 1, mother, xpos, ypos, zpos, idrotm, "ONLY") ;
621       printf(" idrotm %3i phi %6.1f(%5.3f) xpos %7.2f ypos %7.2f zpos %7.2f \n", 
622       idrotm, phi, phiRad, xpos, ypos, zpos);
623       nr++;
624     } else { // TRD1 
625       TString smName("SMOD"); // 12-oct-05
626       if(i==5 && gn.Contains("110DEG")) {
627         smName = "SM10";
628         nrsmod = nr;
629         nr     = 0;
630       }
631       phi    = g->GetArm1PhiMin() + dphi*(2*i+1)/2.; // phi= 70, 90, 110, 130, 150, 170
632       phiRad = phi*TMath::Pi()/180.;
633
634       AliMatrix(idrotm, 90.0, phi, 90.0, 90.0+phi, 0.0, 0.0);
635
636       xpos = rpos * TMath::Cos(phiRad);
637       ypos = rpos * TMath::Sin(phiRad);
638       zpos = smodPar2; // 21-sep-04
639       if(i==5 && gn.Contains("110DEG")) {
640         xpos += (par1C/2. * TMath::Sin(phiRad)); 
641         ypos -= (par1C/2. * TMath::Cos(phiRad)); 
642       }
643       
644       // 1th module in z-direction;
645       gMC->Gspos(smName.Data(), ++nr, mother, xpos, ypos, zpos, idrotm, "ONLY") ;
646       AliDebug(3, Form(" %s : %2i idrotm %3i phi %6.1f(%5.3f) xpos %7.2f ypos %7.2f zpos %7.2f : i %i \n", 
647                        smName.Data(), nr, idrotm, phi, phiRad, xpos, ypos, zpos, i));
648       // 2th module in z-direction;
649       if(gn.Contains("TWIST") || gn.Contains("TRD")) {
650       // turn arround X axis; 0<phi<360
651         double phiy = 90. + phi + 180.;
652         if(phiy>=360.) phiy -= 360.;
653  
654         AliMatrix(idrotm, 90.0, phi, 90.0, phiy, 180.0, 0.0);
655         gMC->Gspos(smName.Data(), ++nr, mother, xpos, ypos, -zpos, idrotm, "ONLY");
656         AliDebug(3, Form(" %s : %2i idrotm %3i phiy %6.1f  xpos %7.2f ypos %7.2f zpos %7.2f \n", 
657                          smName.Data(), nr, idrotm, phiy, xpos, ypos, -zpos));
658       } else {
659         gMC->Gspos("SMOD", ++nr, mother, xpos, ypos, -zpos, idrotm, "ONLY");
660       }
661     }
662   }
663   printf(" Number of Super Modules %i \n", nr+nrsmod);
664 }
665
666 void AliEMCALv0::CreateEmod(const char* mother, const char* child)
667 { // 17-may-05; mother="SMOD"; child="EMOD"
668   AliEMCALGeometry * g = GetGeometry(); 
669   TString gn(g->GetName()); gn.ToUpper(); 
670   // Module definition
671   Double_t par[10], parTubs[5], xpos=0., ypos=0., zpos=0., rpos=0.;
672   Double_t parSCPA[5], zposSCPA=0.; // passive SC - 13-MAY-05, TRD1 case
673   Double_t trd1Angle = g->GetTrd1Angle()*TMath::DegToRad(), tanTrd1 = TMath::Tan(trd1Angle/2.);
674   Double_t tanTrd2y  = TMath::Tan(g->GetTrd2AngleY()*TMath::DegToRad()/2.);
675   int nr=0;
676   idrotm=0;
677   if(!gn.Contains("TRD")) { // standard module
678     par[0] = (sampleWidth*g->GetNECLayers())/2.; 
679     par[1] = par[2] = g->GetPhiModuleSize()/2.;
680     gMC->Gsvolu(child, "BOX", idtmed[idSTEEL], par, 3);
681
682   } else if (gn.Contains("TRD1")){ // TRD1 system coordinate iz differnet
683     if(strcmp(mother,"SMOD")==0) {
684       parEMOD[0] = g->GetEtaModuleSize()/2.;   // dx1
685       parEMOD[1] = g->Get2Trd1Dx2()/2.;        // dx2
686       parEMOD[2] = g->GetPhiModuleSize()/2.;;  // dy
687       parEMOD[3] = g->GetLongModuleSize()/2.;  // dz
688       gMC->Gsvolu(child, "TRD1", idtmed[idSTEEL], parEMOD, 4);
689       if(gn.Contains("WSUC") || gn.Contains("MAY05")){
690         parSCPA[0] = g->GetEtaModuleSize()/2. + tanTrd1*g->GetFrontSteelStrip();   // dx1
691         parSCPA[1] = parSCPA[0]               + tanTrd1*g->GetPassiveScintThick(); // dx2
692         parSCPA[2] = g->GetPhiModuleSize()/2.;     // dy
693         parSCPA[3] = g->GetPassiveScintThick()/2.; // dz
694         gMC->Gsvolu("SCPA", "TRD1", idtmed[idSC], parSCPA, 4);
695         zposSCPA   = -parEMOD[3] + g->GetFrontSteelStrip() + g->GetPassiveScintThick()/2.;
696         gMC->Gspos ("SCPA", ++nr, child, 0.0, 0.0, zposSCPA, 0, "ONLY");
697       }
698     }
699   } else if (gn.Contains("TRD2")){ // TRD2 as for TRD1 - 27-jan-05
700     parEMOD[0] = g->GetEtaModuleSize()/2.;   // dx1
701     parEMOD[1] = g->Get2Trd1Dx2()/2.;        // dx2
702     parEMOD[2] = g->GetPhiModuleSize()/2.;   // dy1
703     parEMOD[3] = parEMOD[2] + tanTrd2y*g->GetLongModuleSize();// dy2
704     parEMOD[4] = g->GetLongModuleSize()/2.;  // dz
705     gMC->Gsvolu(child, "TRD2", idtmed[idSTEEL], parEMOD, 5);
706   }
707
708   nr   = 0;
709   if(gn.Contains("TWIST")) { // 13-sep-04
710     fShishKebabModules = new TList;
711     AliEMCALShishKebabModule *mod=0, *mTmp; // current module
712     for(int iz=0; iz<g->GetNZ(); iz++) {
713       //for(int iz=0; iz<4; iz++) {
714       if(iz==0) {
715         mod    = new AliEMCALShishKebabModule();
716         idrotm = 0;
717       } else {
718         mTmp = new AliEMCALShishKebabModule(*mod);
719         mod  = mTmp;
720         Double_t  angle = mod->GetThetaInDegree();
721         AliMatrix(idrotm, angle,0., 90.0,90.0, 90.-angle, 180.);
722       }
723       fShishKebabModules->Add(mod);
724
725       xpos = mod->GetPosXfromR() + g->GetSteelFrontThickness() - smodPar0;
726       zpos = mod->GetPosZ() - smodPar2;
727       for(int iy=0; iy<g->GetNPhi(); iy++) {
728         ypos = g->GetPhiModuleSize()*(2*iy+1 - g->GetNPhi())/2.;
729         gMC->Gspos(child, ++nr, mother, xpos, ypos, zpos, idrotm, "ONLY") ;
730         //        printf(" %3i(%2i,2i) xpos %7.2f ypos %7.2f zpos %7.2f \n", nr,iy,iz, xpos, ypos, zpos);
731       }
732     }    
733   } else if(gn.Contains("TRD")) { // 30-sep-04; 27-jan-05 - as for TRD1 as for TRD2
734     // X->Z(0, 0); Y->Y(90, 90); Z->X(90, 0)
735     AliEMCALShishKebabTrd1Module *mod=0, *mTmp; // current module
736
737     for(int iz=0; iz<g->GetNZ(); iz++) {
738       Double_t  angle=90., phiOK=0;
739       if(gn.Contains("TRD1")) {
740         mod = (AliEMCALShishKebabTrd1Module*)fShishKebabModules->At(iz);
741         angle = mod->GetThetaInDegree();
742         if(!gn.Contains("WSUC")) { // ALICE 
743           if(iz==0) AliMatrix(idrotm, 0.,0., 90.,90., 90.,0.); // z'(x); y'(y); x'(z)
744           else      AliMatrix(idrotm, 90.-angle,180., 90.0,90.0, angle, 0.);
745           phiOK = mod->GetCenterOfModule().Phi()*180./TMath::Pi(); 
746           //          printf(" %2i | angle | %6.3f - %6.3f = %6.3f(eta %5.3f)\n", 
747           //iz+1, angle, phiOK, angle-phiOK, mod->GetEtaOfCenterOfModule());
748           xpos = mod->GetPosXfromR() + g->GetSteelFrontThickness() - smodPar0;
749           zpos = mod->GetPosZ() - smodPar2;
750
751           int iyMax = g->GetNPhi();
752           if(strcmp(mother,"SMOD") && gn.Contains("110DEG")) {
753             iyMax /= 2;
754           }
755           for(int iy=0; iy<iyMax; iy++) { // flat in phi
756             ypos = g->GetPhiModuleSize()*(2*iy+1 - iyMax)/2.;
757             gMC->Gspos(child, ++nr, mother, xpos, ypos, zpos, idrotm, "ONLY") ;
758         //printf(" %2i xpos %7.2f ypos %7.2f zpos %7.2f idrotm %i\n", nr, xpos, ypos, zpos, idrotm);
759             AliDebug(3,Form("%3.3i(%2.2i,%2.2i) ", nr,iy+1,iz+1));
760           }
761           printf("\n");
762         } else {
763           if(iz==0) AliMatrix(idrotm, 0.,0., 90.,0., 90.,90.); // (x')z; y'(x); z'(y)
764           else      AliMatrix(idrotm, 90-angle,270., 90.0,0.0, angle,90.);
765           phiOK = mod->GetCenterOfModule().Phi()*180./TMath::Pi(); 
766           printf(" %2i | angle -phiOK | %6.3f - %6.3f = %6.3f(eta %5.3f)\n", 
767           iz+1, angle, phiOK, angle-phiOK, mod->GetEtaOfCenterOfModule());
768           zpos = mod->GetPosZ()      - smodPar2;
769           ypos = mod->GetPosXfromR() - smodPar1;
770           printf(" zpos %7.2f ypos %7.2f idrotm %i\n xpos ", zpos, xpos, idrotm);
771           for(int ix=0; ix<g->GetNPhi(); ix++) { // flat in phi
772             xpos = g->GetPhiModuleSize()*(2*ix+1 - g->GetNPhi())/2.;
773             gMC->Gspos(child, ++nr, mother, xpos, ypos, zpos, idrotm, "ONLY") ;
774             printf(" %7.2f ", xpos);
775           }
776           printf("\n");
777         }
778       } else if(gn.Contains("TRD2")){ // 1-feb-05 - TRD2;  curve in phi
779         double angEtaRow = 0.;
780         double theta1=0.,phi1=0., theta2=0.,phi2=0., theta3=0.,phi3=0.;
781         angle=90.;
782         if(iz==0) {
783           mod   = new AliEMCALShishKebabTrd1Module();
784         } else {
785           mTmp  = new AliEMCALShishKebabTrd1Module(*mod);
786           mod   = mTmp;
787           angle = mod->GetThetaInDegree();
788         }
789
790         fShishKebabModules->Add(mod);
791         phiOK = mod->GetCenterOfModule().Phi()*180./TMath::Pi(); 
792         printf(" %i | theta | %6.3f - %6.3f = %6.3f\n", iz+1, angle, phiOK, angle-phiOK);
793
794         zpos = mod->GetPosZ() - parTubs[2];
795         rpos = parTubs[0] + mod->GetPosXfromR();
796
797         angle     = mod->GetThetaInDegree();
798         Double_t stepAngle = (parTubs[4] -  parTubs[3])/g->GetNPhi(); // 11-mar-04
799         for(int iy=0; iy<g->GetNPhi(); iy++) {
800           angEtaRow = parTubs[3] + stepAngle*(0.5+double(iy));
801           //          angEtaRow = 0;
802           theta1  = 90. +  angle; phi1 = angEtaRow;      // x' ~-z;
803           theta2  = 90.;          phi2 = 90. + angEtaRow;// y' ~ y;
804           theta3  = angle;        phi3 = angEtaRow;      // z' ~ x;
805           if(phi3 < 0.0) phi3 += 360.; 
806           AliMatrix(idrotm, theta1,phi1, theta2,phi2, theta3,phi3);
807
808           xpos = rpos * TMath::Cos(angEtaRow*TMath::DegToRad());
809           ypos = rpos * TMath::Sin(angEtaRow*TMath::DegToRad());
810           gMC->Gspos(child, ++nr, "SMOP", xpos, ypos, zpos, idrotm, "ONLY") ;
811           // SMON; 
812           phi1    = 180 + angEtaRow;
813           theta3  = 180.-theta3;  phi3 = angEtaRow;
814           AliMatrix(idrotm, theta1,phi1, theta2,phi2, theta3,phi3);
815           gMC->Gspos(child,  nr, "SMON", xpos, ypos, -zpos, idrotm, "ONLY") ;
816           if(0) {
817             printf(" angEtaRow(phi) %7.2f |  angle(eta) %7.2f \n",  angEtaRow, angle);
818             printf("iy=%2i xpos %7.2f ypos %7.2f zpos %7.2f idrotm %i\n", iy, xpos, ypos, zpos, idrotm);
819           }
820         } // for(int iy=0; iy<g->GetNPhi(); iy++)
821       }
822     } 
823   } else {
824     xpos = g->GetSteelFrontThickness()/2.;
825     for(int iz=0; iz<g->GetNZ(); iz++) {
826       zpos = -smodPar2 + g->GetEtaModuleSize()*(2*iz+1)/2.;
827       for(int iy=0; iy<g->GetNPhi(); iy++) {
828         ypos = g->GetPhiModuleSize()*(2*iy+1 - g->GetNPhi())/2.;
829         gMC->Gspos(child, ++nr, mother, xpos, ypos, zpos, 0, "ONLY") ;
830       //printf(" %2i xpos %7.2f ypos %7.2f zpos %7.2f \n", nr, xpos, ypos, zpos);
831       }
832     }
833   }
834   printf(" Number of modules in Super Module %i \n", nr);
835 }
836
837 // 8-dec-04 by PAI
838 void AliEMCALv0::Trd1Tower3X3(const double parSCM0[4])
839 {// PB should be for whole SCM0 - ?
840   double parTRAP[11], *dummy=0;
841   AliEMCALGeometry * g = GetGeometry(); 
842   TString gn(g->GetName()), scmx; 
843   gn.ToUpper(); 
844  // Division to tile size 
845   Info("Trd1Tower3X3()","<I> Divide SCM0 on y-axis %i", g->GetNETAdiv());
846   gMC->Gsdvn("SCMY","SCM0", g->GetNETAdiv(), 2); // y-axis
847   double dx1=parSCM0[0], dx2=parSCM0[1], dy=parSCM0[2], dz=parSCM0[3];
848   double ndiv=3., xpos=0.0;
849   // should be defined once
850   gMC->Gsvolu("PBTI", "BOX", idtmed[idPB], dummy, 0);
851   if(gn.Contains("TEST")==0) { // one name for all trapesoid
852     scmx = "SCMX"; 
853     gMC->Gsvolu(scmx.Data(), "TRAP", idtmed[idSC], dummy, 0);
854   }
855
856   
857   for(int ix=1; ix<=3; ix++) { // 3X3
858     // ix=1
859     parTRAP[0] = dz;
860     parTRAP[1] = TMath::ATan2((dx2-dx1)/2.,2.*dz)*TMath::RadToDeg(); // theta
861     parTRAP[2] = 0.;           // phi
862     // bottom
863     parTRAP[3] = dy/ndiv;      // H1
864     parTRAP[4] = dx1/ndiv;     // BL1
865     parTRAP[5] = parTRAP[4];   // TL1
866     parTRAP[6] = 0.0;          // ALP1
867     // top
868     parTRAP[7] = dy/ndiv;      // H2
869     parTRAP[8] = dx2/ndiv;     // BL2
870     parTRAP[9] = parTRAP[8];   // TL2
871     parTRAP[10]= 0.0;          // ALP2
872     xpos = +(dx1+dx2)/3.;      // 6 or 3
873
874     if      (ix==3) {
875       parTRAP[1] = -parTRAP[1];
876       xpos = -xpos;
877     } else if(ix==2) { // central part is box but we treat as trapesoid due to numbering
878       parTRAP[1] = 0.;
879       parTRAP[8] = dx1/ndiv;     // BL2
880       parTRAP[9] = parTRAP[8];   // TL2
881       xpos = 0.0;
882     }
883     printf(" ** TRAP ** xpos %9.3f\n", xpos);
884     for(int i=0; i<11; i++) printf(" par[%2.2i] %9.4f\n", i, parTRAP[i]);
885     if(gn.Contains("TEST")){
886       scmx = "SCX"; scmx += ix;
887       gMC->Gsvolu(scmx.Data(), "TRAP", idtmed[idSC], parTRAP, 11);
888       gMC->Gspos(scmx.Data(), 1, "SCMY", xpos, 0.0, 0.0, 0, "ONLY") ;
889     } else {
890       gMC->Gsposp(scmx.Data(), ix, "SCMY", xpos, 0.0, 0.0, 0, "ONLY", parTRAP, 11) ;
891     }
892     PbInTrap(parTRAP, scmx);
893   }
894
895   Info("Trd1Tower3X3()", "Ver. 1.0 : was tested.");
896 }
897
898 // 8-dec-04 by PAI
899 void AliEMCALv0::PbInTrap(const double parTRAP[11], TString n)
900 {// see for example CreateShishKebabGeometry(); just for case TRD1
901   static int nr=0;
902   printf(" Pb tiles : nrstart %i\n", nr);
903   AliEMCALGeometry * g = GetGeometry(); 
904
905   double par[3];
906   double sampleWidth = double(g->GetECPbRadThick()+g->GetECScintThick());
907   double xpos = 0.0, ypos = 0.0;
908   double zpos = -sampleWidth*g->GetNECLayers()/2. + g->GetECPbRadThick()/2.;
909
910   double coef = (parTRAP[8] -  parTRAP[4]) / (2.*parTRAP[0]);
911   double xCenterSCMX =  (parTRAP[4] +  parTRAP[8])/2.; // ??
912   //  double tan = TMath::Tan(parTRAP[1]*TMath::DegToRad());
913
914   par[1] = parTRAP[3];              // y 
915   par[2] = g->GetECPbRadThick()/2.; // z
916   for(int iz=0; iz<g->GetNECLayers(); iz++){
917     par[0] = parTRAP[4] + coef*sampleWidth*iz;
918     xpos   = par[0] - xCenterSCMX;
919     if(parTRAP[1] < 0.) xpos = -xpos;
920     gMC->Gsposp("PBTI", ++nr, n.Data(), xpos, ypos, zpos, 0, "ONLY", par, 3) ;
921     printf(" %i xpos %9.3f zpos %9.3f par[0] %9.3f |", iz+1, xpos, zpos, par[0]);
922     zpos += sampleWidth;
923     if(iz%2>0) printf("\n");
924   } 
925   printf(" Number of Pb tiles in SCMX %i coef %9.7f \n", nr, coef);
926   printf(" par[1] %9.3f  par[2] %9.3f ypos %9.3f \n", par[1], par[2], ypos); 
927   Info("PbInTrap", "Ver. 1.0 : was tested.");
928 }
929
930 // 8-dec-04 by PAI
931 void AliEMCALv0::Trd1Tower4X4()
932 {// see for example CreateShishKebabGeometry()
933 }
934 // 3-feb-05
935 void AliEMCALv0::Scm0InTrd2(const AliEMCALGeometry * g, const Double_t parEMOD[5], Double_t parSCM0[5])
936 {
937   // Passive material inside the detector
938   double wallThickness = g->GetPhiModuleSize()/2. -  g->GetPhiTileSize(); //Need check
939   printf(" wall thickness %7.5f \n", wallThickness);
940   for(int i=0; i<4; i++) { // on pictures sometimes I can not see 0 -> be carefull!!
941     parSCM0[i] = parEMOD[i] - wallThickness;
942     printf(" %i parSCMO %7.3f parEMOD %7.3f : dif %7.3f \n", i, parSCM0[i],parEMOD[i], parSCM0[i]-parEMOD[i]);
943   }
944   parSCM0[4] = parEMOD[4];
945   gMC->Gsvolu("SCM0", "TRD2", idtmed[idSC], parSCM0, 5); // idAIR -> idSC
946   gMC->Gspos("SCM0", 1, "EMOD", 0., 0., 0., 0, "ONLY") ;
947   // Division 
948   if(g->GetNPHIdiv()==2 && g->GetNETAdiv()==2) {
949     Division2X2InScm0(g, parSCM0);
950   } else {
951     Info("Scm0InTrd2"," no division SCM0 in this geometry |%s|\n", g->GetName());
952     assert(0);
953   }
954 }
955
956 void AliEMCALv0::Division2X2InScm0(const AliEMCALGeometry * g, const Double_t parSCM0[5])
957 {
958   Double_t parTRAP[11], xpos=0.,ypos=0., dx1=0.0,dx2=0.,dy1=0.0,dy2=0.,dz=0.
959   ,dr1=0.0,dr2=0.;
960   idrotm=0;
961
962   Info("Division2X2InScm0","Divide SCM0 on y-axis %i\n", g->GetNETAdiv());
963   TString n("SCMX"), overLapFlagSCMY("ONLY"), overLapFlagSCMX("ONLY");
964   n = "SCM0"; // for testing - 14-mar-05
965   if(n=="SCM0"){
966     PbInTrapForTrd2(parSCM0, n);
967     // overLapFlagSCMY=overLapFlagSCMX="MANY"; // do not work
968     return;
969   }
970
971   dy1 = parSCM0[2] , dy2 = parSCM0[3], dz = parSCM0[4];
972
973   parTRAP[0] = parSCM0[4];    // dz
974   parTRAP[1] = TMath::ATan2((dy2-dy1)/2.,2.*dz)*TMath::RadToDeg();
975   parTRAP[2] = 90.;           // phi
976   // bottom
977   parTRAP[3] = parSCM0[2]/2.; // H1
978   parTRAP[4] = parSCM0[0];    // BL1
979   parTRAP[5] = parTRAP[4];    // TL1
980   parTRAP[6] = 0.0;           // ALP1
981   // top
982   parTRAP[7] = parSCM0[3]/2.; // H2
983   parTRAP[8] = parSCM0[1];    // BL2
984   parTRAP[9] = parTRAP[8];    // TL2
985   parTRAP[10]= 0.0;           // ALP2
986   printf(" ** SCMY ** \n");
987   for(int i=0; i<11; i++) printf(" par[%2.2i] %9.4f\n", i, parTRAP[i]);
988
989   idrotm=0;
990   gMC->Gsvolu("SCMY", "TRAP", idtmed[idSC], parTRAP, 11); // idAIR -> idSC
991   ypos = +(parTRAP[3]+parTRAP[7])/2.; //
992   printf(" Y shift SCMY inside SCM0 : %7.3f : opt %s\n", ypos, overLapFlagSCMY.Data()); 
993   gMC->Gspos("SCMY", 1, "SCM0", 0.0, ypos, 0.0, idrotm,  overLapFlagSCMY.Data()) ;
994   // Rotation SCMY around z-axis on 180 degree; x'=-x; y'=-y and z=z
995   AliMatrix(idrotm, 90.0,180., 90.0, 270.0, 0.0,0.0) ;
996   // We may have problem with numeration due to rotation - 4-feb-05
997   gMC->Gspos("SCMY", 2, "SCM0", 0.0, -ypos, 0.0, idrotm,  overLapFlagSCMY.Data()); 
998
999   Info("Division2X2InScm0","Divide SCMY on x-axis %i\n", g->GetNPHIdiv());
1000   dx1 = parSCM0[0]; 
1001   dx2 = parSCM0[1]; 
1002   dr1=TMath::Sqrt(dx1*dx1+dy1*dy1);
1003   dr2=TMath::Sqrt(dx2*dx2+dy2*dy2);
1004
1005   parTRAP[0] = parSCM0[4];    // dz
1006   parTRAP[1] = TMath::ATan2((dr2-dr1)/2.,2.*dz)*TMath::RadToDeg(); // 
1007   parTRAP[2] = 45.;           // phi
1008   // bottom
1009   parTRAP[3] = parSCM0[2]/2.; // H1
1010   parTRAP[4] = parSCM0[0]/2.; // BL1
1011   parTRAP[5] = parTRAP[4];    // TL1
1012   parTRAP[6] = 0.0;           // ALP1
1013   // top
1014   parTRAP[7] = parSCM0[3]/2.; // H2
1015   parTRAP[8] = parSCM0[1]/2;  // BL2
1016   parTRAP[9] = parTRAP[8];    // TL2
1017   parTRAP[10]= 0.0;           // ALP2
1018   printf(" ** SCMX ** \n");
1019   for(int i=0; i<11; i++) printf(" par[%2.2i] %9.4f\n", i, parTRAP[i]);
1020
1021   idrotm=0;
1022   gMC->Gsvolu("SCMX", "TRAP", idtmed[idSC], parTRAP, 11);
1023   xpos = (parTRAP[4]+parTRAP[8])/2.;
1024   printf(" X shift SCMX inside SCMX : %7.3f : opt %s\n", xpos, overLapFlagSCMX.Data()); 
1025   gMC->Gspos("SCMX", 1, "SCMY", xpos, 0.0, 0.0, idrotm,  overLapFlagSCMX.Data()) ;
1026   //  AliMatrix(idrotm, 90.0,270., 90.0, 0.0, 0.0,0.0); // x'=-y; y'=x; z'=z
1027   AliMatrix(idrotm, 90.0,90., 90.0, -180.0, 0.0,0.0);     // x'=y;  y'=-x; z'=z
1028   gMC->Gspos("SCMX", 2, "SCMY", -xpos, 0.0, 0.0, idrotm,  overLapFlagSCMX.Data()) ;
1029   // PB:
1030   if(n=="SCMX" && overLapFlagSCMY == "ONLY") {
1031     PbInTrapForTrd2(parTRAP, n);
1032   }
1033
1034
1035 // 4-feb-05 by PAI
1036 void AliEMCALv0::PbInTrapForTrd2(const double *parTRAP, TString name)
1037 {// TRD2 cases
1038   Double_t *dummy=0;
1039   TString pbShape("BOX"), pbtiChonly("ONLY");
1040   if(name=="SCM0") {
1041     pbShape    = "TRD2";
1042     //    pbtiChonly = "MANY";
1043   }
1044   gMC->Gsvolu("PBTI", pbShape.Data(), idtmed[idPB], dummy, 0);
1045
1046   int nr=0;
1047   Info("PbInTrapForTrd2"," Pb tiles inside %s: shape %s :pbtiChonly %s\n nrstart %i\n", 
1048   name.Data(), pbShape.Data(), pbtiChonly.Data(), nr);
1049   AliEMCALGeometry * g = GetGeometry(); 
1050
1051   double par[5], parPB[5], parSC[5];
1052   double sampleWidth = double(g->GetECPbRadThick()+g->GetECScintThick());
1053   double xpos = 0.0, ypos = 0.0;
1054   double zpos = -sampleWidth*g->GetNECLayers()/2. + g->GetECPbRadThick()/2.;
1055   if(name == "SCMX") { // common trapezoid - 11 parameters
1056     double coef = (parTRAP[8] -  parTRAP[4]) / (2.*parTRAP[0]);
1057     double xCenterSCMX =  (parTRAP[4] +  parTRAP[8])/2.; // the same for y
1058     printf(" xCenterSCMX %8.5f : coef %8.7f \n", xCenterSCMX, coef);
1059
1060     par[2] = g->GetECPbRadThick()/2.; // z
1061     for(int iz=0; iz<g->GetNECLayers(); iz++){
1062       par[0] = parTRAP[4] + coef*sampleWidth*iz;
1063       par[1] = par[0];
1064       xpos   = ypos = par[0] - xCenterSCMX;
1065     //if(parTRAP[1] < 0.) xpos = -xpos;
1066       gMC->Gsposp("PBTI", ++nr, name.Data(), xpos, ypos, zpos, 0, "ONLY", par, 3) ;
1067       printf(" %2.2i xpos %8.5f zpos %6.3f par[0,1] %6.3f |", iz+1, xpos, zpos, par[0]);
1068       if(iz%2>0) printf("\n");
1069       zpos += sampleWidth;
1070     } 
1071     printf(" Number of Pb tiles in SCMX %i coef %9.7f \n", nr, coef);
1072     printf(" par[1] %9.5f  par[2] %9.5f ypos %9.5f \n", par[1], par[2], ypos); 
1073   } else if(name == "SCM0") { // 1-mar-05 ; TRD2 - 5 parameters
1074     printf(" SCM0 par = ");
1075     for(int i=0; i<5; i++) printf(" %9.5f ", parTRAP[i]);
1076     printf("\n zpos %f \n",zpos);
1077
1078     double tanx = (parTRAP[1] -  parTRAP[0]) / (2.*parTRAP[4]); //  tanx =  tany now
1079     double tany = (parTRAP[3] -  parTRAP[2]) / (2.*parTRAP[4]), ztmp=0.;
1080     parPB[4] = g->GetECPbRadThick()/2.;
1081     parSC[2] = g->GetECScintThick()/2.;
1082     for(int iz=0; iz<g->GetNECLayers(); iz++){
1083       ztmp     = sampleWidth*double(iz);
1084       parPB[0] = parTRAP[0] + tanx*ztmp;
1085       parPB[1] = parPB[0]   + tanx*g->GetECPbRadThick();
1086       parPB[2] = parTRAP[2] + tany*ztmp;
1087       parPB[3] = parPB[2]   + tany*g->GetECPbRadThick();
1088       gMC->Gsposp("PBTI", ++nr, name.Data(), xpos, ypos, zpos, 0, pbtiChonly.Data(), parPB, 5) ;
1089       printf("\n PBTI %2i | zpos %6.3f | par = ", nr, zpos);
1090       /*
1091       for(int i=0; i<5; i++) printf(" %9.5f ", parPB[i]);
1092       // individual SC tile
1093       parSC[0] = parPB[0];
1094       parSC[1] = parPB[1];
1095       gMC->Gsposp("SCTI", nr, name.Data(), xpos, ypos, zpos+g->GetECScintThick(), 
1096       0, pbtiChonly.Data(), parSC, 3) ;
1097       printf("\n SCTI     zpos %6.3f | par = ", zpos+g->GetECScintThick());
1098       for(int i=0; i<3; i++) printf(" %9.5f ", parPB[i]);
1099       */
1100       zpos  += sampleWidth;
1101     }
1102     printf("\n");
1103   }
1104   Info("PbInTrapForTrd2", "Ver. 0.03 : was tested.");
1105 }
1106
1107 // 15-mar-05
1108 void AliEMCALv0::PbmoInTrd2(const AliEMCALGeometry * g, const Double_t parEMOD[5], Double_t parPBMO[5])
1109 {
1110   Info("PbmoInTrd2"," started : geometry %s ", g->GetName());
1111   double wallThickness = g->GetPhiModuleSize()/2. -  g->GetPhiTileSize();
1112   printf(" wall thickness %7.5f \n", wallThickness);
1113   for(int i=0; i<4; i++) {
1114     parPBMO[i] = parEMOD[i] - wallThickness;
1115     printf(" %i parPBMO %7.3f parEMOD %7.3f : dif %7.3f \n", 
1116     i, parPBMO[i],parEMOD[i], parPBMO[i]-parEMOD[i]);
1117   }
1118   parPBMO[4] = parEMOD[4];
1119   gMC->Gsvolu("PBMO", "TRD2", idtmed[idPB], parPBMO, 5);
1120   gMC->Gspos("PBMO", 1, "EMOD", 0., 0., 0., 0, "ONLY") ;
1121   // Division 
1122   if(g->GetNPHIdiv()==2 && g->GetNETAdiv()==2) {
1123     Division2X2InPbmo(g, parPBMO);
1124     printf(" PBMO, division 2X2 | geometry |%s|\n", g->GetName());
1125   } else {
1126     printf(" no division PBMO in this geometry |%s|\n", g->GetName());
1127     assert(0);
1128   }
1129 }
1130
1131 void AliEMCALv0::Division2X2InPbmo(const AliEMCALGeometry * g, const Double_t parPBMO[5])
1132 {
1133   Info("Division2X2InPbmo"," started : geometry %s ", g->GetName());
1134   //Double_t *dummy=0;
1135   //  gMC->Gsvolu("SCTI", "BOX", idtmed[idSC], dummy, 0);
1136
1137   double parSC[3];
1138   double sampleWidth = double(g->GetECPbRadThick()+g->GetECScintThick());
1139   double xpos = 0.0, ypos = 0.0, zpos = 0.0, ztmp=0;;
1140   double tanx = (parPBMO[1] -  parPBMO[0]) / (2.*parPBMO[4]); //  tanx =  tany now
1141   double tany = (parPBMO[3] -  parPBMO[2]) / (2.*parPBMO[4]);
1142   char name[10], named[10], named2[10];
1143
1144   printf(" PBMO par = ");
1145   for(int i=0; i<5; i++) printf(" %9.5f ", parPBMO[i]);
1146   printf("\n");
1147
1148   parSC[2] = g->GetECScintThick()/2.;
1149   zpos = -sampleWidth*g->GetNECLayers()/2. + g->GetECPbRadThick() + g->GetECScintThick()/2.;
1150   printf(" parSC[2] %9.5f \n", parSC[2]);
1151   for(int iz=0; iz<g->GetNECLayers(); iz++){
1152     ztmp     = g->GetECPbRadThick() + sampleWidth*double(iz); // Z for previous PB
1153     parSC[0] =  parPBMO[0] + tanx*ztmp;
1154     parSC[1] =  parPBMO[2] + tany*ztmp;
1155
1156     sprintf(name,"SC%2.2i", iz+1);
1157     gMC->Gsvolu(name, "BOX", idtmed[idSC], parSC, 3);
1158     gMC->Gspos(name, 1, "PBMO", xpos, ypos, zpos, 0, "ONLY") ;
1159     printf("%s | zpos %6.3f | parSC[0,1]=(%7.5f,%7.5f) -> ", 
1160     name, zpos, parSC[0], parSC[1]);
1161     
1162     sprintf(named,"SY%2.2i", iz+1);
1163     printf(" %s -> ", named);
1164     gMC->Gsdvn(named,name, 2, 2);
1165
1166     sprintf(named2,"SX%2.2i", iz+1);
1167     printf(" %s \n", named2);
1168     gMC->Gsdvn(named2,named, 2, 1);
1169
1170     zpos    += sampleWidth;
1171   }
1172 }
1173
1174 AliEMCALShishKebabTrd1Module* AliEMCALv0::GetShishKebabModule(Int_t neta)
1175 { // 28-oct-05
1176   AliEMCALShishKebabTrd1Module* trd1=0;
1177   if(fShishKebabModules && neta>=0 && neta<fShishKebabModules->GetSize()) {
1178     trd1 = (AliEMCALShishKebabTrd1Module*)fShishKebabModules->At(neta);
1179   }
1180   return trd1;
1181 }