]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMD.cxx
Major upgrade of AliRoot code
[u/mrichter/AliRoot.git] / PMD / AliPMD.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 Revision 1.12  2000/12/04 08:48:18  alibrary
19 Fixing problems in the HEAD
20
21 Revision 1.11  2000/11/17 10:15:24  morsch
22 Call to AliDetector::ResetHits() added to method  AliPMD::ResetHits()
23
24 Revision 1.10  2000/11/06 09:07:13  morsch
25 Set  fRecPoints to zero in default constructor.
26
27 Revision 1.9  2000/10/30 09:03:59  morsch
28 Prototype for PMD reconstructed hits (AliPMDRecPoint) added.
29
30 Revision 1.8  2000/10/20 06:24:40  fca
31 Put the PMD at the right position in the event display
32
33 Revision 1.7  2000/10/02 21:28:12  fca
34 Removal of useless dependecies via forward declarations
35
36 Revision 1.6  2000/01/19 17:17:06  fca
37 Introducing a list of lists of hits -- more hits allowed for detector now
38
39 Revision 1.5  1999/09/29 09:24:27  fca
40 Introduction of the Copyright and cvs Log
41
42 */
43
44 ///////////////////////////////////////////////////////////////////////////////
45 //
46 //                                                                           //
47 //  Photon Multiplicity Detector                                             //
48 //  This class contains the basic functions for the Photon Multiplicity      //
49 //  Detector. Functions specific to one particular geometry are              //
50 //  contained in the derived classes                                         //
51 //                                                                           //
52 //Begin_Html
53 /*
54 <img src="picts/AliPMDClass.gif">
55 </pre>
56 <br clear=left>
57 <font size=+2 color=red>
58 <p>The responsible person for this module is
59 <a href="mailto:sub@vecdec.veccal.ernet.in">Subhasis Chattopadhyay</a>.
60 </font>
61 <pre>
62 */
63 //End_Html
64 //                                                                           //
65 ///////////////////////////////////////////////////////////////////////////////
66
67 #include <TBRIK.h>
68 #include <TNode.h>
69 #include <TTree.h>
70 #include <TGeometry.h>
71 #include <TClonesArray.h>
72 #include <TFile.h>
73
74 #include "AliPMD.h"
75 #include "AliRun.h"
76 #include "AliMC.h" 
77 #include "AliConst.h" 
78 #include "AliPMDRecPoint.h"
79   
80 ClassImp(AliPMD)
81  
82 //_____________________________________________________________________________
83 AliPMD::AliPMD()
84 {
85   //
86   // Default constructor
87   //
88   fIshunt = 0;
89
90   // Always make the TClonesArray, otherwise the automatic streamer gets angry
91   fRecPoints  = new TClonesArray("AliPMDRecPoint",10000); 
92
93 }
94  
95 //_____________________________________________________________________________
96 AliPMD::AliPMD(const char *name, const char *title)
97   : AliDetector(name,title)
98 {
99   //
100   // Default constructor
101   //
102
103   // 
104   // Allocate the array of hits
105   fHits   = new TClonesArray("AliPMDhit",  405);
106   gAlice->AddHitList(fHits);
107
108   fRecPoints  = new TClonesArray("AliPMDRecPoint",10000); 
109   fNRecPoints = 0;
110   
111
112   fIshunt =  1;
113   
114   fPar[0] = 1;
115   fPar[1] = 1;
116   fPar[2] = 0.8;
117   fPar[3] = 0.02;
118   fIn[0]  = 6;
119   fIn[1]  = 20;
120   fIn[2]  = 600;
121   fIn[3]  = 27;
122   fIn[4]  = 27;
123   fGeo[0] = 0;
124   fGeo[1] = 0.2;
125   fGeo[2] = 4;
126   fPadSize[0] = 0.8;
127   fPadSize[1] = 1.0;
128   fPadSize[2] = 1.2;
129   fPadSize[3] = 1.5;
130 }
131
132 AliPMD::~AliPMD()
133 {
134   //
135   // Default constructor
136   //
137     delete fRecPoints;
138     fNRecPoints=0;
139 }
140
141 //_____________________________________________________________________________
142 void AliPMD::AddHit(Int_t track, Int_t *vol, Float_t *hits)
143 {
144   //
145   // Add a PMD hit
146   //
147   TClonesArray &lhits = *fHits;
148   AliPMDhit *newcell, *curcell;
149 //    printf("PMD++ Adding energy %f, prim %d, vol %d %d %d %d\n",
150 //       hits[3],gAlice->GetPrimary(track-1),vol[0],vol[1],vol[2],vol[3]);
151   newcell = new AliPMDhit(fIshunt, track, vol, hits);
152   Int_t i;
153   for (i=0; i<fNhits; i++) {
154     //
155     // See if this cell has already been hit
156     curcell=(AliPMDhit*) lhits[i];
157     if (*curcell==*newcell) {
158 //        printf("Cell with same numbers found\n") ; curcell->Print();
159       *curcell = *curcell+*newcell;
160 //        printf("Cell after addition\n") ; curcell->Print();
161       delete newcell;
162       return;
163     }
164   }
165   new(lhits[fNhits++]) AliPMDhit(newcell);
166   delete newcell;
167 }
168  
169 //_____________________________________________________________________________
170 void AliPMD::BuildGeometry()
171 {
172   //
173   // Build simple ROOT TNode geometry for event display
174   //
175
176   TNode *Node, *Top;
177   const int kColorPMD  = kRed;
178
179   //
180   Top=gAlice->GetGeometry()->GetNode("alice");
181
182   // PMD
183   new TBRIK("S_PMD","PMD box","void",300,300,5);
184   Top->cd();
185   Node = new TNode("PMD","PMD","S_PMD",0,0,-600,"");
186   Node->SetLineColor(kColorPMD);
187   fNodes->Add(Node);
188 }
189
190 //_____________________________________________________________________________
191 Int_t AliPMD::DistancetoPrimitive(Int_t , Int_t )
192 {
193   //
194   // Distance from mouse to detector on the screen
195   // dummy routine
196   //
197    return 9999;
198 }
199  
200 //_____________________________________________________________________________
201 void AliPMD::SetPAR(Float_t p1, Float_t p2, Float_t p3,Float_t p4)
202 {
203   //
204   // Set PMD parameters
205   //
206   fPar[0] = p1;
207   fPar[1] = p2;
208   fPar[2] = p3;
209   fPar[3] = p4;
210 }
211  
212 //_____________________________________________________________________________
213 void AliPMD::SetIN(Float_t p1, Float_t p2, Float_t p3,Float_t p4,Float_t p5)
214 {
215   //
216   // Set PMD parameters
217   //
218   fIn[0] = p1;
219   fIn[1] = p2;
220   fIn[2] = p3;
221   fIn[3] = p4;
222   fIn[4] = p5;
223 }
224  
225 //_____________________________________________________________________________
226 void AliPMD::SetGEO(Float_t p1, Float_t p2, Float_t p3)
227 {
228   //
229   // Set geometry parameters
230   //
231   fGeo[0] = p1;
232   fGeo[1] = p2;
233   fGeo[2] = p3;
234 }
235  
236 //_____________________________________________________________________________
237 void AliPMD::SetPadSize(Float_t p1, Float_t p2, Float_t p3,Float_t p4)
238 {
239   //
240   // Set pad size
241   //
242   fPadSize[0] = p1;
243   fPadSize[1] = p2;
244   fPadSize[2] = p3;
245   fPadSize[3] = p4;
246 }
247  
248 //_____________________________________________________________________________
249 void AliPMD::StepManager()
250 {
251   //
252   // Called at every step in PMD
253   //
254 }
255
256 void AliPMD::AddRecPoint(const AliPMDRecPoint &p)
257 {
258     //
259     // Add a PMD reconstructed hit to the list
260     //
261     TClonesArray &lrecpoints = *fRecPoints;
262     new(lrecpoints[fNRecPoints++]) AliPMDRecPoint(p);
263 }
264
265 void AliPMD::MakeBranch(Option_t* option, char *file)
266 {
267     // Create Tree branches for the PMD
268     
269     char *cR = strstr(option,"R");
270     
271     AliDetector::MakeBranch(option,file);
272
273     if (cR) {
274       printf("Make Branch - TreeR address %p\n",gAlice->TreeR());
275     
276       const Int_t kBufferSize = 4000;
277       char branchname[30];
278       
279       sprintf(branchname,"%sRecPoints",GetName());
280       if (fRecPoints   && gAlice->TreeR()) {
281         gAlice->MakeBranchInTree(gAlice->TreeR(), 
282                                  branchname, &fRecPoints, kBufferSize, file) ;
283       }
284    }    
285 }
286
287
288 void AliPMD::SetTreeAddress()
289 {
290   // Set branch address for the TreeR
291     char branchname[30];
292     AliDetector::SetTreeAddress();
293
294     TBranch *branch;
295     TTree *treeR = gAlice->TreeR();
296
297     sprintf(branchname,"%s",GetName());
298     if (treeR && fRecPoints) {
299         branch = treeR->GetBranch(branchname);
300         if (branch) branch->SetAddress(&fRecPoints);
301     }
302 }
303
304 void AliPMD::ResetHits()
305 {
306   //
307   // Reset number of hits and the hits array
308   //
309     AliDetector::ResetHits();
310     fNRecPoints   = 0;
311     if (fRecPoints)   fRecPoints->Clear();
312 }
313
314 ///////////////////////////////////////////////////////////////////////////////
315 //                                                                           //
316 //  Photon Multiplicity Detector Version 1                                   //
317 //                                                                           //
318 //Begin_Html
319 /*
320 <img src="picts/AliPMDv1Class.gif">
321 */
322 //End_Html
323 //                                                                           //
324 ///////////////////////////////////////////////////////////////////////////////
325
326
327
328 ClassImp(AliPMDhit)
329   
330 //_____________________________________________________________________________
331 AliPMDhit::AliPMDhit(Int_t shunt,Int_t track, Int_t *vol, Float_t *hits):
332   AliHit(shunt, track)
333 {
334   //
335   // Add a PMD hit
336   //
337   Int_t i;
338   for (i=0;i<5;i++) fVolume[i] = vol[i];
339   fX=hits[0];
340   fY=hits[1];
341   fZ=hits[2];
342   fEnergy=hits[3];
343 }
344