]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALv0.cxx
This update is a step in the continuous development of EMCAL.
[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
16 /* $Id$ */
17
18 //_________________________________________________________________________
19 // Implementation version v0 of EMCAL Manager class 
20 // An object of this class does not produce hits nor digits
21 // It is the one to use if you do not want to produce outputs in TREEH or TREED
22 // This class places a Geometry of the EMCAL in the ALICE Detector as defined in AliEMCALGeometry.cxx                 
23 //*-- Author: Yves Schutz (SUBATECH)
24 //*-- and   : Sahal Yacoob (LBL / UCT)
25
26 // --- ROOT system ---
27 #include "TPGON.h"
28 #include "TTUBS.h"
29 #include "TNode.h"
30 #include "TRandom.h"
31 #include "TGeometry.h"
32
33
34 // --- Standard library ---
35
36 #include <stdio.h>
37 #include <string.h>
38 #include <stdlib.h>
39 #include <strstream.h>
40 #include <iostream.h>
41
42 // --- AliRoot header files ---
43
44 #include "AliEMCALv0.h"
45 #include "AliEMCALGeometry.h"
46 #include "AliConst.h"
47 #include "AliRun.h"
48 #include "AliMC.h"
49
50 ClassImp(AliEMCALv0)
51
52 //______________________________________________________________________
53 AliEMCALv0::AliEMCALv0(const char *name, const char *title):
54     AliEMCAL(name,title){
55     // Standard Constructor
56
57     if (strcmp(GetTitle(),"") != 0 )
58         fGeom =  AliEMCALGeometry::GetInstance(GetTitle(), "") ;
59
60 }
61 //______________________________________________________________________
62 void AliEMCALv0::BuildGeometry(){
63     // Display Geometry for display.C
64
65     const Int_t kColorArm1   = kBlue ;
66
67     // Difine the shape of the Calorimeter 
68   
69     new TTUBS("Envelop1", "Tubs that contains arm 1", "void", 
70               fGeom->GetEnvelop(0),     // rmin 
71               fGeom->GetEnvelop(1) +30 ,     // rmax
72               fGeom->GetEnvelop(2)/2.0, // half length in Z
73               fGeom->GetArm1PhiMin(),   // minimun phi angle
74               fGeom->GetArm1PhiMax()    // maximun phi angle
75         );
76     // Place the Node
77     TNode * envelop1node = new TNode("Envelop1", "Arm1 Envelop", "Envelop1") ;
78     envelop1node->SetLineColor(kColorArm1) ;
79     fNodes->Add(envelop1node) ;
80 }
81 //______________________________________________________________________
82 void AliEMCALv0::CreateGeometry(){
83     // Create the EMCAL geometry for Geant
84
85     AliEMCALv0 *emcaltmp = (AliEMCALv0*)gAlice->GetModule("EMCAL") ;
86
87     if ( emcaltmp == NULL ) {
88         Warning("CreateGeometry","detector not found!");
89         return;
90     } // end if
91     // Get pointer to the array containing media indices
92     Int_t *idtmed = fIdtmed->GetArray() - 1599 ;
93
94     // Create an Envelope within which to place the Detector 
95  
96     Float_t envelopA[5] ; 
97     envelopA[0] = fGeom->GetEnvelop(0) ;         // rmin
98     envelopA[1] = fGeom->GetEnvelop(1) + 30 ;    // rmax
99     envelopA[2] = fGeom->GetEnvelop(2) / 2.0 ;   // dz
100     envelopA[3] = fGeom->GetArm1PhiMin() ;       // minimun phi angle
101     envelopA[4] = fGeom->GetArm1PhiMax() ;       // maximun phi angle
102
103     gMC->Gsvolu("XEN1", "TUBS ", idtmed[1599], envelopA, 5) ; //filled with air
104
105     // Create the shapes of active material (LEAD/Aluminium/Scintillator) to be placed 
106
107     Float_t envelopB[10]; // First Layer of Aluminium
108     Float_t envelopC[10]; // Scintillator Layers
109     Float_t envelopD[10]; //  Lead Layers
110     
111     envelopC[0] = envelopD[0] =  envelopB[0] = fGeom->GetArm1PhiMin(); //starting position in Phi
112     envelopC[1] = envelopD[1] =  envelopB[1] = fGeom->GetArm1PhiMax() -  // Angular size of the Detector in Phi
113                                                fGeom->GetArm1PhiMin();
114     envelopC[2] = envelopD[2] =  envelopB[2] = fGeom->GetNPhi() ; // Number of Section in Phi
115     envelopD[3] = envelopC[3] = envelopB[3] = 2; // each section will be passed 2 z coordinates    
116
117     envelopB[4] = (fGeom->GetEnvelop(0) + fGeom->GetGap2Active()) /
118                   (tan(2*atan(exp(0.7)))) ;                          // z co-ordinate 1
119     envelopB[5] = fGeom->GetEnvelop(0) + fGeom->GetGap2Active(); //rmin at z1
120     envelopD[6] = envelopB[6] = envelopB[5] + 3.18;  //rmax at z1
121     envelopB[7] = (fGeom->GetEnvelop(0) + fGeom->GetGap2Active()) /
122                   (tan(2*atan(exp(-0.7)))) ;              // z co-ordinate 2
123     envelopB[8] = envelopB[5] ; //
124     envelopB[9] = envelopB[6] ; // radii are the same.
125
126     // filled shapes wit hactive material 
127     gMC->Gsvolu("XALU", "PGON", idtmed[1602], envelopB, 10); // Define Aluminium volume completely
128     gMC->Gsvolu("XPST", "PGON", idtmed[1601], 0, 0) ;  // The polystyrene layers will be defined when placed 
129     gMC->Gsvolu("XPBX", "PGON", idtmed[1600], 0, 0) ; //  as will the lead layers 
130     gMC->Gsdvn("XPHI", "XPST", fGeom->GetNPhi(), 2) ; //  Dividind eta polystyrene divisions into phi segments.
131     
132     Int_t idrotm = 1;
133     AliMatrix(idrotm, 90.0, 0., 90.0, 90.0, 0.0, 0.0) ;
134
135     // Position the Envelope in Alice  
136     gMC->Gspos("XEN1", 1, "ALIC", 0.0, 0.0, 0.0, idrotm, "MANY") ;
137     // Position Aluminium Layer in the Envelope 
138     gMC->Gspos("XALU", 1, "XEN1", 0.0, 0.0, 0.0 , idrotm, "ONLY") ;
139
140 // The loop below places the scintillator in Lead Layers alternately.
141  
142     for (int i = 0; i < (fGeom->GetNLayers()); i++ ){
143         envelopC[5] = envelopD[6] ; //rmin
144         envelopC[6] = envelopD[6] + ((i > +2)  ? 0.5 : 0.6)  ;  //rmax larger for first two layers (preshower)
145         envelopC[8] = envelopD[6] ; //rmin
146         envelopC[9] = envelopD[6] + ((i > 2 ) ? 0.5 : 0.6)  ;  //rmax larger for first two layers (preshower)
147         for (int j =0; j < (fGeom->GetNZ()) ; j++){
148             envelopC[4] = envelopD[6]/tan(2*atan(exp(0.7-(j*1.4/
149                                                       (fGeom->GetNZ()))))); //z begin  
150             envelopC[7] = envelopD[6]/tan(2*atan(exp(0.7-((j+1)*1.4/
151                                                       (fGeom->GetNZ()))))); // z end 
152             gMC->Gsposp("XPST",1+j+i*(fGeom->GetNZ()), "XEN1", 
153                         0.0, 0.0, 0.0 , idrotm, "ONLY", envelopC, 10); // Position and define layer
154         } // end for j
155         if (i < (fGeom->GetNLayers()-1)){
156             envelopD[5] = envelopC[6] ; //rmin
157             envelopD[6] = envelopC[6] + 0.5;  //rmax
158             envelopD[8] = envelopC[6] ; //rmin
159             envelopD[9] = envelopC[6] + 0.5;  //rmax
160             for (int j =0; j < (fGeom->GetNZ()) ; j++){
161                 envelopD[4] = envelopC[6]/tan(2*atan(exp(0.7-(j*1.4/
162                                                        (fGeom->GetNZ()))))); // z begin 
163                 envelopD[7] = envelopC[6]/tan(2*atan(exp(0.7-((j+1)*1.4/
164                                                         (fGeom->GetNZ()))))); // z end
165                 gMC->Gsposp("XPBX",1+ j+i*(fGeom->GetNZ()), "XEN1", 
166                             0.0, 0.0, 0.0 , idrotm, "MANY", envelopD, 10) ; // Position and Define Layer
167             } // end for j
168         } // end if i
169     }  // for i
170 }
171 //______________________________________________________________________
172 void AliEMCALv0::Init(void){
173     // Just prints an information message
174     Int_t i;
175
176     cout << endl;
177     for(i=0;i<35;i++) cout <<"*";
178     cout << " EMCAL_INIT ";
179     for(i=0;i<35;i++) cout << "*";
180     cout << endl;
181
182     // Here the EMCAL initialisation code (if any!)
183
184     if (fGeom!=0)  
185         cout << "AliEMCAL" << Version() << " : EMCAL geometry intialized for "
186              << fGeom->GetName() << endl ;
187     else
188         cout << "AliEMCAL" << Version() << 
189             " : EMCAL geometry initialization failed !" << endl ;
190     for(i=0;i<80;i++) printf("*");
191     cout << endl;
192 }