]> git.uio.no Git - u/mrichter/AliRoot.git/blob - CPV/AliCPVv0.cxx
Introduction of the Copyright and cvs Log
[u/mrichter/AliRoot.git] / CPV / AliCPVv0.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 /*
17 $Log$
18 */
19
20 /////////////////////////////////////////////////////////
21 //  Manager and hits classes for set:CPV version 0     //
22 //  Coarse geometry                                    //
23 //                                                     //
24 //  Author: Yuri Kharlov, IHEP, Protvino               //
25 //  e-mail: Yuri.Kharlov@cern.ch                       //
26 //  Last modified: 17 September 1999                   //
27 /////////////////////////////////////////////////////////
28  
29 // --- ROOT system ---
30 #include "TH1.h"
31 #include "TRandom.h"
32 #include "TFile.h"
33 #include "TTree.h"
34 #include "TBRIK.h"
35 #include "TNode.h"
36
37 // --- galice header files ---
38 #include "AliCPVv0.h"
39 #include "AliRun.h"
40 #include "AliMC.h" 
41
42 ClassImp(AliCPVv0)
43
44 //==============================================================================
45 //                            AliCPVv0
46 //==============================================================================
47
48 AliCPVv0::AliCPVv0()
49 {
50 }
51  
52 //______________________________________________________________________________
53
54 AliCPVv0::AliCPVv0(const char *name, const char *title)
55           : AliCPV(name, title)
56 {
57 }
58  
59 //______________________________________________________________________________
60
61 void AliCPVv0::CreateGeometry()
62 {
63
64   AliCPV *CPV_tmp = (AliCPV*)gAlice->GetModule("CPV");
65   if( NULL==CPV_tmp )
66   {
67     printf("There isn't CPV detector!\n");
68     return;
69   }
70
71   Int_t    rotation_matrix_number=0;
72   Float_t  par[3],
73            x,y,z;
74
75   // CPV creation
76
77   par[0] = GetPadZSize()   / 2 * GetNz();
78   par[1] = GetPadPhiSize() / 2 * GetNphi();
79   par[2] = GetThickness()  / 2;
80   gMC->Gsvolu("CPV","BOX ",GetCPVIdtmed(),par,3);
81
82   for( Int_t i=0; i<GetNofCradles(); i++ )
83   {
84     Float_t cradle_angle_pos = -90+(i-(GetNofCradles()-1)/2.) * GetAngle();
85
86     // Cradles are numerated in clock reversed order. (general way of angle increment)
87
88     Float_t r = GetRadius() + GetThickness()/2;
89     x = r*cos(cradle_angle_pos*kPI/180);
90     y = r*sin(cradle_angle_pos*kPI/180);
91     z = 0;
92     AliMatrix(rotation_matrix_number, 0,0 , 90,90+cradle_angle_pos , 90,180+cradle_angle_pos);
93     gMC->Gspos("CPV",i+1,"ALIC",x,y,z,rotation_matrix_number,"ONLY");
94   }
95   AddCPVCradles();
96
97 }
98
99 //______________________________________________________________________________
100
101 void AliCPVv0::StepManager()
102 {
103
104 //    if( gMC->IsTrackEntering() ) {
105 //      const char *VolumeName = gMC->CurrentVolName();
106 //      cout << "AliCPVv0::StepManager() entered to CPV to the volume " 
107 //           << VolumeName << "!\n";
108 //    }
109   
110   if( strcmp(gMC->CurrentVolName(),"CPV ")==0 && gMC->IsTrackEntering() )
111   {
112     // GEANT particle just have entered into CPV detector.
113
114     AliCPV &CPV = *(AliCPV*)gAlice->GetModule("CPV");
115
116     Int_t CradleNumber;
117     gMC->CurrentVolOffID(0,CradleNumber);
118     CradleNumber--;
119     AliCPVCradle  &cradle = CPV.GetCradle(CradleNumber);
120
121     // Current position of the hit in the cradle ref. system
122
123     TLorentzVector xyzt;
124     gMC -> TrackPosition(xyzt);
125     Float_t xyzm[3], xyzd[3], xyd[2];
126     for (Int_t i=0; i<3; i++) xyzm[i] = xyzt[i];
127     gMC -> Gmtod (xyzm, xyzd, 1);
128     for (Int_t i=0; i<2; i++) xyd[i]  = xyzd[i];
129
130     // Current momentum of the hit's track in the MARS ref. system
131     
132     TLorentzVector  pmom;
133     gMC -> TrackMomentum(pmom);
134     Int_t ipart = gMC->TrackPid();
135
136     // Store current particle in the list of Cradle particles.
137     cradle.AddHit(pmom,xyd,ipart);
138
139 //      if (gMC->TrackCharge()!=0.) CPV.Reconstruction(0,0);
140
141 //      CPV.Print("p");
142 //      CPV.Print("r");
143
144   }
145 }