1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 Revision 1.22 2000/09/12 18:11:13 fca
19 zero hits area before using
21 Revision 1.21 2000/07/21 10:21:07 morsch
22 fNrawch = 0; and fNrechits = 0; in the default constructor.
24 Revision 1.20 2000/07/10 15:28:39 fca
25 Correction of the inheritance scheme
27 Revision 1.19 2000/06/30 16:29:51 dibari
28 Added kDebugLevel variable to control output size on demand
30 Revision 1.18 2000/06/12 15:15:46 jbarbosa
33 Revision 1.17 2000/06/09 14:58:37 jbarbosa
34 New digitisation per particle type
36 Revision 1.16 2000/04/19 12:55:43 morsch
37 Newly structured and updated version (JB, AM)
42 ////////////////////////////////////////////////
43 // Manager and hits classes for set:RICH //
44 ////////////////////////////////////////////////
52 #include <TObjArray.h>
55 #include <TParticle.h>
60 #include "AliRICHSegmentation.h"
61 #include "AliRICHHit.h"
62 #include "AliRICHCerenkov.h"
63 #include "AliRICHPadHit.h"
64 #include "AliRICHDigit.h"
65 #include "AliRICHTransientDigit.h"
66 #include "AliRICHRawCluster.h"
67 #include "AliRICHRecHit.h"
68 #include "AliRICHHitMapA1.h"
69 #include "AliRICHClusterFinder.h"
74 #include "AliPoints.h"
75 #include "AliCallf77.h"
79 // Static variables for the pad-hit iterator routines
80 static Int_t sMaxIterPad=0;
81 static Int_t sCurIterPad=0;
82 static TClonesArray *fClusters2;
83 static TClonesArray *fHits2;
88 //___________________________________________
91 // Default constructor for RICH manager class
105 //___________________________________________
106 AliRICH::AliRICH(const char *name, const char *title)
107 : AliDetector(name,title)
111 <img src="gif/alirich.gif">
115 fHits = new TClonesArray("AliRICHHit",1000 );
116 gAlice->AddHitList(fHits);
117 fPadHits = new TClonesArray("AliRICHPadHit",100000);
118 fCerenkovs = new TClonesArray("AliRICHCerenkov",1000);
119 gAlice->AddHitList(fCerenkovs);
120 //gAlice->AddHitList(fHits);
125 fNdch = new Int_t[kNCH];
127 fDchambers = new TObjArray(kNCH);
129 fRecHits = new TObjArray(kNCH);
133 for (i=0; i<kNCH ;i++) {
134 (*fDchambers)[i] = new TClonesArray("AliRICHDigit",10000);
138 fNrawch = new Int_t[kNCH];
140 fRawClusters = new TObjArray(kNCH);
141 //printf("Created fRwClusters with adress:%p",fRawClusters);
143 for (i=0; i<kNCH ;i++) {
144 (*fRawClusters)[i] = new TClonesArray("AliRICHRawCluster",10000);
148 fNrechits = new Int_t[kNCH];
150 for (i=0; i<kNCH ;i++) {
151 (*fRecHits)[i] = new TClonesArray("AliRICHRecHit",1000);
153 //printf("Created fRecHits with adress:%p",fRecHits);
156 SetMarkerColor(kRed);
159 AliRICH::AliRICH(const AliRICH& RICH)
165 //___________________________________________
169 // Destructor of RICH manager class
177 //___________________________________________
178 void AliRICH::AddHit(Int_t track, Int_t *vol, Float_t *hits)
182 // Adds a hit to the Hits list
185 TClonesArray &lhits = *fHits;
186 new(lhits[fNhits++]) AliRICHHit(fIshunt,track,vol,hits);
188 //_____________________________________________________________________________
189 void AliRICH::AddCerenkov(Int_t track, Int_t *vol, Float_t *cerenkovs)
193 // Adds a RICH cerenkov hit to the Cerenkov Hits list
196 TClonesArray &lcerenkovs = *fCerenkovs;
197 new(lcerenkovs[fNcerenkovs++]) AliRICHCerenkov(fIshunt,track,vol,cerenkovs);
198 //printf ("Done for Cerenkov %d\n\n\n\n",fNcerenkovs);
200 //___________________________________________
201 void AliRICH::AddPadHit(Int_t *clhits)
205 // Add a RICH pad hit to the list
208 TClonesArray &lPadHits = *fPadHits;
209 new(lPadHits[fNPadHits++]) AliRICHPadHit(clhits);
211 //_____________________________________________________________________________
212 void AliRICH::AddDigits(Int_t id, Int_t *tracks, Int_t *charges, Int_t *digits)
216 // Add a RICH digit to the list
219 TClonesArray &ldigits = *((TClonesArray*)(*fDchambers)[id]);
220 new(ldigits[fNdch[id]++]) AliRICHDigit(tracks,charges,digits);
223 //_____________________________________________________________________________
224 void AliRICH::AddRawCluster(Int_t id, const AliRICHRawCluster& c)
227 // Add a RICH digit to the list
230 TClonesArray &lrawcl = *((TClonesArray*)(*fRawClusters)[id]);
231 new(lrawcl[fNrawch[id]++]) AliRICHRawCluster(c);
234 //_____________________________________________________________________________
235 void AliRICH::AddRecHit(Int_t id, Float_t *rechit, Float_t *photons, Int_t *padsx, Int_t* padsy)
239 // Add a RICH reconstructed hit to the list
242 TClonesArray &lrec = *((TClonesArray*)(*fRecHits)[id]);
243 new(lrec[fNrechits[id]++]) AliRICHRecHit(id,rechit,photons,padsx,padsy);
246 //___________________________________________
247 void AliRICH::BuildGeometry()
252 // Builds a TNode geometry for event display
256 const int kColorRICH = kGreen;
258 top=gAlice->GetGeometry()->GetNode("alice");
261 new TBRIK("S_RICH","S_RICH","void",71.09999,11.5,73.15);
264 Float_t pos1[3]={0,471.8999,165.2599};
265 //Chamber(0).SetChamberTransform(pos1[0],pos1[1],pos1[2],
266 new TRotMatrix("rot993","rot993",90,0,70.69,90,19.30999,-90);
267 node = new TNode("RICH1","RICH1","S_RICH",pos1[0],pos1[1],pos1[2],"rot993");
270 node->SetLineColor(kColorRICH);
274 Float_t pos2[3]={171,470,0};
275 //Chamber(1).SetChamberTransform(pos2[0],pos2[1],pos2[2],
276 new TRotMatrix("rot994","rot994",90,-20,90,70,0,0);
277 node = new TNode("RICH2","RICH2","S_RICH",pos2[0],pos2[1],pos2[2],"rot994");
280 node->SetLineColor(kColorRICH);
283 Float_t pos3[3]={0,500,0};
284 //Chamber(2).SetChamberTransform(pos3[0],pos3[1],pos3[2],
285 new TRotMatrix("rot995","rot995",90,0,90,90,0,0);
286 node = new TNode("RICH3","RICH3","S_RICH",pos3[0],pos3[1],pos3[2],"rot995");
289 node->SetLineColor(kColorRICH);
292 Float_t pos4[3]={-171,470,0};
293 //Chamber(3).SetChamberTransform(pos4[0],pos4[1],pos4[2],
294 new TRotMatrix("rot996","rot996",90,20,90,110,0,0);
295 node = new TNode("RICH4","RICH4","S_RICH",pos4[0],pos4[1],pos4[2],"rot996");
298 node->SetLineColor(kColorRICH);
301 Float_t pos5[3]={161.3999,443.3999,-165.3};
302 //Chamber(4).SetChamberTransform(pos5[0],pos5[1],pos5[2],
303 new TRotMatrix("rot997","rot997",90,340,108.1999,70,18.2,70);
304 node = new TNode("RICH5","RICH5","S_RICH",pos5[0],pos5[1],pos5[2],"rot997");
306 node->SetLineColor(kColorRICH);
309 Float_t pos6[3]={0., 471.9, -165.3,};
310 //Chamber(5).SetChamberTransform(pos6[0],pos6[1],pos6[2],
311 new TRotMatrix("rot998","rot998",90,0,109.3099,90,19.30999,90);
312 node = new TNode("RICH6","RICH6","S_RICH",pos6[0],pos6[1],pos6[2],"rot998");
315 node->SetLineColor(kColorRICH);
318 Float_t pos7[3]={-161.399,443.3999,-165.3};
319 //Chamber(6).SetChamberTransform(pos7[0],pos7[1],pos7[2],
320 new TRotMatrix("rot999","rot999",90,20,108.1999,110,18.2,110);
321 node = new TNode("RICH7","RICH7","S_RICH",pos7[0],pos7[1],pos7[2],"rot999");
322 node->SetLineColor(kColorRICH);
327 //___________________________________________
328 void AliRICH::CreateGeometry()
331 // Create the geometry for RICH version 1
333 // Modified by: N. Colonna (INFN - BARI, Nicola.Colonna@ba.infn.it)
334 // R.A. Fini (INFN - BARI, Rosanna.Fini@ba.infn.it)
335 // R.A. Loconsole (Bari University, loco@riscom.ba.infn.it)
339 <img src="picts/AliRICHv1.gif">
344 <img src="picts/AliRICHv1Tree.gif">
348 AliRICH *pRICH = (AliRICH *) gAlice->GetDetector("RICH");
349 AliRICHSegmentation* segmentation;
350 AliRICHGeometry* geometry;
351 AliRICHChamber* iChamber;
353 iChamber = &(pRICH->Chamber(0));
354 segmentation=iChamber->GetSegmentationModel(0);
355 geometry=iChamber->GetGeometryModel();
358 distance = geometry->GetFreonThickness()/2 + geometry->GetQuartzThickness() + geometry->GetGapThickness();
359 geometry->SetRadiatorToPads(distance);
362 Int_t *idtmed = fIdtmed->GetArray()-999;
369 // --- Define the RICH detector
370 // External aluminium box
372 par[1] = 11.5; //Original Settings
377 gMC->Gsvolu("RICH", "BOX ", idtmed[1009], par, 3);
379 // Sensitive part of the whole RICH
381 par[1] = 11.5; //Original Settings
386 gMC->Gsvolu("SRIC", "BOX ", idtmed[1000], par, 3);
390 par[1] = .188; //Original Settings
395 gMC->Gsvolu("HONE", "BOX ", idtmed[1001], par, 3);
399 par[1] = .025; //Original Settings
404 gMC->Gsvolu("ALUM", "BOX ", idtmed[1009], par, 3);
407 par[0] = geometry->GetQuartzWidth()/2;
408 par[1] = geometry->GetQuartzThickness()/2;
409 par[2] = geometry->GetQuartzLength()/2;
411 par[1] = .25; //Original Settings
413 /*par[0] = geometry->GetQuartzWidth()/2;
414 par[1] = geometry->GetQuartzThickness()/2;
415 par[2] = geometry->GetQuartzLength()/2;*/
416 //printf("\n\n\n\n\n\n\n\\n\n\n\n Gap Thickness: %f %f %f\n\n\n\n\n\n\n\n\n\n\n\n\n\n",par[0],par[1],par[2]);
417 gMC->Gsvolu("QUAR", "BOX ", idtmed[1002], par, 3);
419 // Spacers (cylinders)
422 par[2] = geometry->GetFreonThickness()/2;
423 gMC->Gsvolu("SPAC", "TUBE", idtmed[1002], par, 3);
427 par[1] = .2; //Original Settings
432 gMC->Gsvolu("OQUA", "BOX ", idtmed[1007], par, 3);
434 // Frame of opaque quartz
435 par[0] = geometry->GetOuterFreonWidth()/2;
436 par[1] = geometry->GetFreonThickness()/2;
437 par[2] = geometry->GetOuterFreonLength()/2 + 1;
439 par[1] = .5; //Original Settings
444 gMC->Gsvolu("OQF1", "BOX ", idtmed[1007], par, 3);
446 par[0] = geometry->GetInnerFreonWidth()/2;
447 par[1] = geometry->GetFreonThickness()/2;
448 par[2] = geometry->GetInnerFreonLength()/2 + 1;
449 gMC->Gsvolu("OQF2", "BOX ", idtmed[1007], par, 3);
451 // Little bar of opaque quartz
453 par[1] = geometry->GetQuartzThickness()/2;
454 par[2] = geometry->GetInnerFreonLength()/2 - 2.4;
456 par[1] = .25; //Original Settings
461 gMC->Gsvolu("BARR", "BOX ", idtmed[1007], par, 3);
464 par[0] = geometry->GetOuterFreonWidth()/2;
465 par[1] = geometry->GetFreonThickness()/2;
466 par[2] = geometry->GetOuterFreonLength()/2;
468 par[1] = .5; //Original Settings
473 gMC->Gsvolu("FRE1", "BOX ", idtmed[1003], par, 3);
475 par[0] = geometry->GetInnerFreonWidth()/2;
476 par[1] = geometry->GetFreonThickness()/2;
477 par[2] = geometry->GetInnerFreonLength()/2;
478 gMC->Gsvolu("FRE2", "BOX ", idtmed[1003], par, 3);
482 par[1] = geometry->GetGapThickness()/2;
483 //printf("\n\n\n\n\n\n\n\\n\n\n\n Gap Thickness: %f\n\n\n\n\n\n\n\n\n\n\n\n\n\n",par[1]);
485 gMC->Gsvolu("META", "BOX ", idtmed[1004], par, 3);
489 par[1] = geometry->GetProximityGapThickness()/2;
490 //printf("\n\n\n\n\n\n\n\\n\n\n\n Gap Thickness: %f\n\n\n\n\n\n\n\n\n\n\n\n\n\n",par[1]);
492 gMC->Gsvolu("GAP ", "BOX ", idtmed[1008], par, 3);
498 gMC->Gsvolu("CSI ", "BOX ", idtmed[1005], par, 3);
504 gMC->Gsvolu("GRID", "TUBE", idtmed[1006], par, 3);
506 // --- Places the detectors defined with GSVOLU
507 // Place material inside RICH
508 gMC->Gspos("SRIC", 1, "RICH", 0., 0., 0., 0, "ONLY");
510 gMC->Gspos("ALUM", 1, "SRIC", 0., 1.276 - geometry->GetGapThickness()/2 - geometry->GetQuartzThickness() - geometry->GetFreonThickness()- .4 -.05 - .376 -.025, 0., 0, "ONLY");
511 gMC->Gspos("HONE", 1, "SRIC", 0., 1.276- geometry->GetGapThickness()/2 - geometry->GetQuartzThickness() - geometry->GetFreonThickness()- .4 -.05 - .188, 0., 0, "ONLY");
512 gMC->Gspos("ALUM", 2, "SRIC", 0., 1.276 - geometry->GetGapThickness()/2 - geometry->GetQuartzThickness() - geometry->GetFreonThickness()- .4 - .025, 0., 0, "ONLY");
513 gMC->Gspos("OQUA", 1, "SRIC", 0., 1.276 - geometry->GetGapThickness()/2 - geometry->GetQuartzThickness() - geometry->GetFreonThickness()- .2, 0., 0, "ONLY");
515 AliMatrix(idrotm[1019], 0., 0., 90., 0., 90., 90.);
517 Int_t nspacers = (Int_t)(TMath::Abs(geometry->GetInnerFreonLength()/14.4));
518 //printf("\n\n\n\n\n\n\n\n\n\n\n\n\n\n\n\n Spacers:%d\n\n\n\n\n\n\n\n\n\n\n\n\n\n\n\n",nspacers);
520 //printf("Nspacers: %d", nspacers);
522 //for (i = 1; i <= 9; ++i) {
523 //zs = (5 - i) * 14.4; //Original settings
524 for (i = 0; i < nspacers; i++) {
525 zs = (TMath::Abs(nspacers/2) - i) * 14.4;
526 gMC->Gspos("SPAC", i, "FRE1", 6.7, 0., zs, idrotm[1019], "ONLY"); //Original settings
527 //gMC->Gspos("SPAC", i, "FRE1", zs, 0., 6.7, idrotm[1019], "ONLY");
529 //for (i = 10; i <= 18; ++i) {
530 //zs = (14 - i) * 14.4; //Original settings
531 for (i = nspacers; i < nspacers*2; ++i) {
532 zs = (nspacers + TMath::Abs(nspacers/2) - i) * 14.4;
533 gMC->Gspos("SPAC", i, "FRE1", -6.7, 0., zs, idrotm[1019], "ONLY"); //Original settings
534 //gMC->Gspos("SPAC", i, "FRE1", zs, 0., -6.7, idrotm[1019], "ONLY");
537 //for (i = 1; i <= 9; ++i) {
538 //zs = (5 - i) * 14.4; //Original settings
539 for (i = 0; i < nspacers; i++) {
540 zs = (TMath::Abs(nspacers/2) - i) * 14.4;
541 gMC->Gspos("SPAC", i, "FRE2", 6.7, 0., zs, idrotm[1019], "ONLY"); //Original settings
542 //gMC->Gspos("SPAC", i, "FRE2", zs, 0., 6.7, idrotm[1019], "ONLY");
544 //for (i = 10; i <= 18; ++i) {
545 //zs = (5 - i) * 14.4; //Original settings
546 for (i = nspacers; i < nspacers*2; ++i) {
547 zs = (nspacers + TMath::Abs(nspacers/2) - i) * 14.4;
548 gMC->Gspos("SPAC", i, "FRE2", -6.7, 0., zs, idrotm[1019], "ONLY"); //Original settings
549 //gMC->Gspos("SPAC", i, "FRE2", zs, 0., -6.7, idrotm[1019], "ONLY");
552 /*gMC->Gspos("FRE1", 1, "OQF1", 0., 0., 0., 0, "ONLY");
553 gMC->Gspos("FRE2", 1, "OQF2", 0., 0., 0., 0, "ONLY");
554 gMC->Gspos("OQF1", 1, "SRIC", 31.3, -4.724, 41.3, 0, "ONLY");
555 gMC->Gspos("OQF2", 2, "SRIC", 0., -4.724, 0., 0, "ONLY");
556 gMC->Gspos("OQF1", 3, "SRIC", -31.3, -4.724, -41.3, 0, "ONLY");
557 gMC->Gspos("BARR", 1, "QUAR", -21.65, 0., 0., 0, "ONLY"); //Original settings
558 gMC->Gspos("BARR", 2, "QUAR", 21.65, 0., 0., 0, "ONLY"); //Original settings
559 gMC->Gspos("QUAR", 1, "SRIC", 0., -3.974, 0., 0, "ONLY");
560 gMC->Gspos("GAP ", 1, "META", 0., 4.8, 0., 0, "ONLY");
561 gMC->Gspos("META", 1, "SRIC", 0., 1.276, 0., 0, "ONLY");
562 gMC->Gspos("CSI ", 1, "SRIC", 0., 6.526, 0., 0, "ONLY");*/
565 gMC->Gspos("FRE1", 1, "OQF1", 0., 0., 0., 0, "ONLY");
566 gMC->Gspos("FRE2", 1, "OQF2", 0., 0., 0., 0, "ONLY");
567 gMC->Gspos("OQF1", 1, "SRIC", geometry->GetOuterFreonWidth()/2 + geometry->GetInnerFreonWidth()/2, 1.276 - geometry->GetGapThickness()/2- geometry->GetQuartzThickness() -geometry->GetFreonThickness()/2, 0., 0, "ONLY"); //Original settings (31.3)
568 gMC->Gspos("OQF2", 2, "SRIC", 0., 1.276 - geometry->GetGapThickness()/2 - geometry->GetQuartzThickness() - geometry->GetFreonThickness()/2, 0., 0, "ONLY"); //Original settings
569 gMC->Gspos("OQF1", 3, "SRIC", - (geometry->GetOuterFreonWidth()/2 + geometry->GetInnerFreonWidth()/2), 1.276 - geometry->GetGapThickness()/2 - geometry->GetQuartzThickness() - geometry->GetFreonThickness()/2, 0., 0, "ONLY"); //Original settings (-31.3)
570 gMC->Gspos("BARR", 1, "QUAR", -21.65, 0., 0., 0, "ONLY"); //Original settings
571 gMC->Gspos("BARR", 2, "QUAR", 21.65, 0., 0., 0, "ONLY"); //Original settings
572 gMC->Gspos("QUAR", 1, "SRIC", 0., 1.276 - geometry->GetGapThickness()/2 - geometry->GetQuartzThickness()/2, 0., 0, "ONLY");
573 gMC->Gspos("GAP ", 1, "META", 0., geometry->GetGapThickness()/2 - geometry->GetProximityGapThickness()/2 - 0.0001, 0., 0, "ONLY");
574 gMC->Gspos("META", 1, "SRIC", 0., 1.276, 0., 0, "ONLY");
575 gMC->Gspos("CSI ", 1, "SRIC", 0., 1.276 + geometry->GetGapThickness()/2 + .25, 0., 0, "ONLY");
577 //printf("Position of the gap: %f to %f\n", 1.276 + geometry->GetGapThickness()/2 - geometry->GetProximityGapThickness()/2 - .2, 1.276 + geometry->GetGapThickness()/2 - geometry->GetProximityGapThickness()/2 + .2);
579 // Place RICH inside ALICE apparatus
581 AliMatrix(idrotm[1000], 90., 0., 70.69, 90., 19.31, -90.);
582 AliMatrix(idrotm[1001], 90., -20., 90., 70., 0., 0.);
583 AliMatrix(idrotm[1002], 90., 0., 90., 90., 0., 0.);
584 AliMatrix(idrotm[1003], 90., 20., 90., 110., 0., 0.);
585 AliMatrix(idrotm[1004], 90., 340., 108.2, 70., 18.2, 70.);
586 AliMatrix(idrotm[1005], 90., 0., 109.31, 90., 19.31, 90.);
587 AliMatrix(idrotm[1006], 90., 20., 108.2, 110., 18.2, 110.);
589 gMC->Gspos("RICH", 1, "ALIC", 0., 471.9, 165.26, idrotm[1000], "ONLY");
590 gMC->Gspos("RICH", 2, "ALIC", 171., 470., 0., idrotm[1001], "ONLY");
591 gMC->Gspos("RICH", 3, "ALIC", 0., 500., 0., idrotm[1002], "ONLY");
592 gMC->Gspos("RICH", 4, "ALIC", -171., 470., 0., idrotm[1003], "ONLY");
593 gMC->Gspos("RICH", 5, "ALIC", 161.4, 443.4, -165.3, idrotm[1004], "ONLY");
594 gMC->Gspos("RICH", 6, "ALIC", 0., 471.9, -165.3, idrotm[1005], "ONLY");
595 gMC->Gspos("RICH", 7, "ALIC", -161.4, 443.4, -165.3, idrotm[1006], "ONLY");
600 //___________________________________________
601 void AliRICH::CreateMaterials()
604 // *** DEFINITION OF AVAILABLE RICH MATERIALS ***
605 // ORIGIN : NICK VAN EIJNDHOVEN
606 // Modified by: N. Colonna (INFN - BARI, Nicola.Colonna@ba.infn.it)
607 // R.A. Fini (INFN - BARI, Rosanna.Fini@ba.infn.it)
608 // R.A. Loconsole (Bari University, loco@riscom.ba.infn.it)
610 Int_t isxfld = gAlice->Field()->Integ();
611 Float_t sxmgmx = gAlice->Field()->Max();
614 /************************************Antonnelo's Values (14-vectors)*****************************************/
616 Float_t ppckov[14] = { 5.63e-9,5.77e-9,5.9e-9,6.05e-9,6.2e-9,6.36e-9,6.52e-9,
617 6.7e-9,6.88e-9,7.08e-9,7.3e-9,7.51e-9,7.74e-9,8e-9 };
618 Float_t rIndexQuarz[14] = { 1.528309,1.533333,
619 1.538243,1.544223,1.550568,1.55777,
620 1.565463,1.574765,1.584831,1.597027,
621 1.611858,1.6277,1.6472,1.6724 };
622 Float_t rIndexOpaqueQuarz[14] = { 1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1. };
623 Float_t rIndexMethane[14] = { 1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1. };
624 Float_t rIndexGrid[14] = { 1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1. };
625 Float_t abscoFreon[14] = { 179.0987,179.0987,
626 179.0987,179.0987,179.0987,142.92,56.65,13.95,10.43,7.07,2.03,.5773,.33496,0. };
627 //Float_t abscoFreon[14] = { 1e-5,1e-5,1e-5,1e-5,1e-5,1e-5,1e-5,1e-5,1e-5,
628 // 1e-5,1e-5,1e-5,1e-5,1e-5 };
629 Float_t abscoQuarz[14] = { 64.035,39.98,35.665,31.262,27.527,22.815,21.04,17.52,
630 14.177,9.282,4.0925,1.149,.3627,.10857 };
631 Float_t abscoOpaqueQuarz[14] = { 1e-5,1e-5,1e-5,1e-5,1e-5,1e-5,1e-5,1e-5,1e-5,
632 1e-5,1e-5,1e-5,1e-5,1e-5 };
633 Float_t abscoCsI[14] = { 1e-4,1e-4,1e-4,1e-4,1e-4,1e-4,1e-4,1e-4,1e-4,1e-4,
634 1e-4,1e-4,1e-4,1e-4 };
635 Float_t abscoMethane[14] = { 1e6,1e6,1e6,1e6,1e6,1e6,1e6,1e6,1e6,1e6,1e6,
637 Float_t abscoGrid[14] = { 1e-4,1e-4,1e-4,1e-4,1e-4,1e-4,1e-4,1e-4,1e-4,1e-4,
638 1e-4,1e-4,1e-4,1e-4 };
639 Float_t efficAll[14] = { 1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1. };
640 Float_t efficCsI[14] = { 6e-4,.005,.0075,.01125,.045,.117,.135,.16575,
641 .17425,.1785,.1836,.1904,.1938,.221 };
642 Float_t efficGrid[14] = { 1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1. };
646 /**********************************End of Antonnelo's Values**********************************/
648 /**********************************Values from rich_media.f (31-vectors)**********************************/
651 //Photons energy intervals
655 ppckov[i] = (Float_t(i)*0.1+5.5)*1e-9;
656 //printf ("Energy intervals: %e\n",ppckov[i]);
660 //Refraction index for quarz
661 Float_t rIndexQuarz[26];
668 Float_t ene=ppckov[i]*1e9;
669 Float_t a=f1/(e1*e1 - ene*ene);
670 Float_t b=f2/(e2*e2 - ene*ene);
671 rIndexQuarz[i] = TMath::Sqrt(1. + a + b );
672 //printf ("rIndexQuarz: %e\n",rIndexQuarz[i]);
675 //Refraction index for opaque quarz, methane and grid
676 Float_t rIndexOpaqueQuarz[26];
677 Float_t rIndexMethane[26];
678 Float_t rIndexGrid[26];
681 rIndexOpaqueQuarz[i]=1;
682 rIndexMethane[i]=1.000444;
684 //printf ("rIndexOpaqueQuarz , etc: %e, %e, %e\n",rIndexOpaqueQuarz[i], rIndexMethane[i], rIndexGrid[i]=1);
687 //Absorption index for freon
688 Float_t abscoFreon[26] = {179.0987, 179.0987, 179.0987, 179.0987, 179.0987, 179.0987, 179.0987, 179.0987,
689 179.0987, 142.9206, 56.64957, 25.58622, 13.95293, 12.03905, 10.42953, 8.804196,
690 7.069031, 4.461292, 2.028366, 1.293013, .577267, .40746, .334964, 0., 0., 0.};
692 //Absorption index for quarz
693 /*Float_t Qzt [21] = {.0,.0,.005,.04,.35,.647,.769,.808,.829,.844,.853,.858,.869,.887,.903,.902,.902,
694 .906,.907,.907,.907};
695 Float_t Wavl2[] = {150.,155.,160.0,165.0,170.0,175.0,180.0,185.0,190.0,195.0,200.0,205.0,210.0,
696 215.0,220.0,225.0,230.0,235.0,240.0,245.0,250.0};
697 Float_t abscoQuarz[31];
698 for (Int_t i=0;i<31;i++)
700 Float_t Xlam = 1237.79 / (ppckov[i]*1e9);
701 if (Xlam <= 160) abscoQuarz[i] = 0;
702 if (Xlam > 250) abscoQuarz[i] = 1;
705 for (Int_t j=0;j<21;j++)
707 //printf ("Passed\n");
708 if (Xlam > Wavl2[j] && Xlam < Wavl2[j+1])
710 Float_t Dabs = (Qzt[j+1] - Qzt[j])/(Wavl2[j+1] - Wavl2[j]);
711 Float_t Abso = Qzt[j] + Dabs*(Xlam - Wavl2[j]);
712 abscoQuarz[i] = -5.0/(TMath::Log(Abso));
716 printf ("abscoQuarz: %e abscoFreon: %e for energy: %e\n",abscoQuarz[i],abscoFreon[i],ppckov[i]);
719 /*Float_t abscoQuarz[31] = {49.64211, 48.41296, 47.46989, 46.50492, 45.13682, 44.47883, 43.1929 , 41.30922, 40.5943 ,
720 39.82956, 38.98623, 38.6247 , 38.43448, 37.41084, 36.22575, 33.74852, 30.73901, 24.25086,
721 17.94531, 11.88753, 5.99128, 3.83503, 2.36661, 1.53155, 1.30582, 1.08574, .8779708,
722 .675275, 0., 0., 0.};
724 for (Int_t i=0;i<31;i++)
726 abscoQuarz[i] = abscoQuarz[i]/10;
729 Float_t abscoQuarz [26] = {105.8, 65.52, 48.58, 42.85, 35.79, 31.262, 28.598, 27.527, 25.007, 22.815, 21.004,
730 19.266, 17.525, 15.878, 14.177, 11.719, 9.282, 6.62, 4.0925, 2.601, 1.149, .667, .3627,
731 .192, .1497, .10857};
733 //Absorption index for methane
734 Float_t abscoMethane[26];
737 abscoMethane[i]=AbsoCH4(ppckov[i]*1e9);
738 //printf("abscoMethane: %e for energy: %e\n", abscoMethane[i],ppckov[i]*1e9);
741 //Absorption index for opaque quarz, csi and grid, efficiency for all and grid
742 Float_t abscoOpaqueQuarz[26];
743 Float_t abscoCsI[26];
744 Float_t abscoGrid[26];
745 Float_t efficAll[26];
746 Float_t efficGrid[26];
749 abscoOpaqueQuarz[i]=1e-5;
754 //printf ("All must be 1: %e, %e, %e, %e, %e\n",abscoOpaqueQuarz[i],abscoCsI[i],abscoGrid[i],efficAll[i],efficGrid[i]);
759 Float_t efficCsI[26] = {0.000199999995, 0.000600000028, 0.000699999975, 0.00499999989, 0.00749999983, 0.010125,
760 0.0242999997, 0.0405000001, 0.0688500032, 0.105299994, 0.121500008, 0.141749993, 0.157949999,
761 0.162, 0.166050002, 0.167669997, 0.174299985, 0.176789999, 0.179279998, 0.182599992, 0.18592,
762 0.187579989, 0.189239994, 0.190899998, 0.207499996, 0.215799987};
766 //FRESNEL LOSS CORRECTION FOR PERPENDICULAR INCIDENCE AND
767 //UNPOLARIZED PHOTONS
771 efficCsI[i] = efficCsI[i]/(1.-Fresnel(ppckov[i]*1e9,1.,0));
772 //printf ("Fresnel result: %e for energy: %e\n",Fresnel(ppckov[i]*1e9,1.,0),ppckov[i]*1e9);
775 /*******************************************End of rich_media.f***************************************/
782 Float_t afre[2], agri, amet[2], aqua[2], ahon, zfre[2], zgri, zhon,
786 Int_t nlmatmet, nlmatqua;
787 Float_t wmatquao[2], rIndexFreon[26];
788 Float_t aquao[2], epsil, stmin, zquao[2];
790 Float_t radlal, densal, tmaxfd, deemax, stemax;
791 Float_t aal, zal, radlgri, densfre, radlhon, densgri, denshon,densqua, densmet, wmatfre[2], wmatmet[2], wmatqua[2];
793 Int_t *idtmed = fIdtmed->GetArray()-999;
795 TGeant3 *geant3 = (TGeant3*) gMC;
797 // --- Photon energy (GeV)
798 // --- Refraction indexes
799 for (i = 0; i < 26; ++i) {
800 rIndexFreon[i] = ppckov[i] * .0172 * 1e9 + 1.177;
801 //rIndexFreon[i] = 1;
802 //printf ("rIndexFreon: %e \n efficCsI: %e for energy: %e\n",rIndexFreon[i], efficCsI[i], ppckov[i]);
805 // --- Detection efficiencies (quantum efficiency for CsI)
806 // --- Define parameters for honeycomb.
807 // Used carbon of equivalent rad. lenght
814 // --- Parameters to include in GSMIXT, relative to Quarz (SiO2)
825 // --- Parameters to include in GSMIXT, relative to opaque Quarz (SiO2)
836 // --- Parameters to include in GSMIXT, relative to Freon (C6F14)
847 // --- Parameters to include in GSMIXT, relative to methane (CH4)
858 // --- Parameters to include in GSMIXT, relative to anode grid (Cu)
865 // --- Parameters to include in GSMATE related to aluminium sheet
872 AliMaterial(1, "Air $", 14.61, 7.3, .001205, 30420., 67500);
873 AliMaterial(6, "HON", ahon, zhon, denshon, radlhon, 0);
874 AliMaterial(16, "CSI", ahon, zhon, denshon, radlhon, 0);
875 AliMixture(20, "QUA", aqua, zqua, densqua, nlmatqua, wmatqua);
876 AliMixture(21, "QUAO", aquao, zquao, densquao, nlmatquao, wmatquao);
877 AliMixture(30, "FRE", afre, zfre, densfre, nlmatfre, wmatfre);
878 AliMixture(40, "MET", amet, zmet, densmet, nlmatmet, wmatmet);
879 AliMixture(41, "METG", amet, zmet, densmet, nlmatmet, wmatmet);
880 AliMaterial(11, "GRI", agri, zgri, densgri, radlgri, 0);
881 AliMaterial(50, "ALUM", aal, zal, densal, radlal, 0);
889 AliMedium(1, "DEFAULT MEDIUM AIR$", 1, 0, isxfld, sxmgmx, tmaxfd, stemax, deemax, epsil, stmin);
890 AliMedium(2, "HONEYCOMB$", 6, 0, isxfld, sxmgmx, tmaxfd, stemax, deemax, epsil, stmin);
891 AliMedium(3, "QUARZO$", 20, 1, isxfld, sxmgmx, tmaxfd, stemax, deemax, epsil, stmin);
892 AliMedium(4, "FREON$", 30, 1, isxfld, sxmgmx, tmaxfd, stemax, deemax, epsil, stmin);
893 AliMedium(5, "METANO$", 40, 1, isxfld, sxmgmx, tmaxfd, stemax, deemax, epsil, stmin);
894 AliMedium(6, "CSI$", 16, 1, isxfld, sxmgmx,tmaxfd, stemax, deemax, epsil, stmin);
895 AliMedium(7, "GRIGLIA$", 11, 0, isxfld, sxmgmx, tmaxfd, stemax, deemax, epsil, stmin);
896 AliMedium(8, "QUARZOO$", 21, 1, isxfld, sxmgmx, tmaxfd, stemax, deemax, epsil, stmin);
897 AliMedium(9, "GAP$", 41, 1, isxfld, sxmgmx,tmaxfd, .1, -deemax, epsil, -stmin);
898 AliMedium(10, "ALUMINUM$", 50, 1, isxfld, sxmgmx, tmaxfd, stemax, deemax, epsil, stmin);
901 geant3->Gsckov(idtmed[1000], 26, ppckov, abscoMethane, efficAll, rIndexMethane);
902 geant3->Gsckov(idtmed[1001], 26, ppckov, abscoMethane, efficAll, rIndexMethane);
903 geant3->Gsckov(idtmed[1002], 26, ppckov, abscoQuarz, efficAll,rIndexQuarz);
904 geant3->Gsckov(idtmed[1003], 26, ppckov, abscoFreon, efficAll,rIndexFreon);
905 geant3->Gsckov(idtmed[1004], 26, ppckov, abscoMethane, efficAll, rIndexMethane);
906 geant3->Gsckov(idtmed[1005], 26, ppckov, abscoCsI, efficCsI, rIndexMethane);
907 geant3->Gsckov(idtmed[1006], 26, ppckov, abscoGrid, efficGrid, rIndexGrid);
908 geant3->Gsckov(idtmed[1007], 26, ppckov, abscoOpaqueQuarz, efficAll, rIndexOpaqueQuarz);
909 geant3->Gsckov(idtmed[1008], 26, ppckov, abscoMethane, efficAll, rIndexMethane);
910 geant3->Gsckov(idtmed[1009], 26, ppckov, abscoGrid, efficGrid, rIndexGrid);
913 //___________________________________________
915 Float_t AliRICH::Fresnel(Float_t ene,Float_t pdoti, Bool_t pola)
918 //ENE(EV), PDOTI=COS(INC.ANG.), PDOTR=COS(POL.PLANE ROT.ANG.)
920 Float_t en[36] = {5.0,5.1,5.2,5.3,5.4,5.5,5.6,5.7,5.8,5.9,6.0,6.1,6.2,
921 6.3,6.4,6.5,6.6,6.7,6.8,6.9,7.0,7.1,7.2,7.3,7.4,7.5,7.6,7.7,
922 7.8,7.9,8.0,8.1,8.2,8.3,8.4,8.5};
925 Float_t csin[36] = {2.14,2.21,2.33,2.48,2.76,2.97,2.99,2.59,2.81,3.05,
926 2.86,2.53,2.55,2.66,2.79,2.96,3.18,3.05,2.84,2.81,2.38,2.11,
927 2.01,2.13,2.39,2.73,3.08,3.15,2.95,2.73,2.56,2.41,2.12,1.95,
930 Float_t csik[36] = {0.,0.,0.,0.,0.,0.196,0.408,0.208,0.118,0.49,0.784,0.543,
931 0.424,0.404,0.371,0.514,0.922,1.102,1.139,1.376,1.461,1.253,0.878,
932 0.69,0.612,0.649,0.824,1.347,1.571,1.678,1.763,1.857,1.824,1.824,
935 Int_t j=Int_t(xe*10)-49;
936 Float_t cn=csin[j]+((csin[j+1]-csin[j])/0.1)*(xe-en[j]);
937 Float_t ck=csik[j]+((csik[j+1]-csik[j])/0.1)*(xe-en[j]);
939 //FORMULAE FROM HANDBOOK OF OPTICS, 33.23 OR
940 //W.R. HUNTER, J.O.S.A. 54 (1964),15 , J.O.S.A. 55(1965),1197
942 Float_t sinin=TMath::Sqrt(1-pdoti*pdoti);
943 Float_t tanin=sinin/pdoti;
945 Float_t c1=cn*cn-ck*ck-sinin*sinin;
946 Float_t c2=4*cn*cn*ck*ck;
947 Float_t aO=TMath::Sqrt(0.5*(TMath::Sqrt(c1*c1+c2)+c1));
948 Float_t b2=0.5*(TMath::Sqrt(c1*c1+c2)-c1);
950 Float_t rs=((aO-pdoti)*(aO-pdoti)+b2)/((aO+pdoti)*(aO+pdoti)+b2);
951 Float_t rp=rs*((aO-sinin*tanin)*(aO-sinin*tanin)+b2)/((aO+sinin*tanin)*(aO+sinin*tanin)+b2);
954 //CORRECTION FACTOR FOR SURFACE ROUGHNESS
955 //B.J. STAGG APPLIED OPTICS, 30(1991),4113
958 Float_t lamb=1240/ene;
961 Float_t rO=TMath::Exp(-(4*TMath::Pi()*pdoti*sigraf/lamb)*(4*TMath::Pi()*pdoti*sigraf/lamb));
965 Float_t pdotr=0.8; //DEGREE OF POLARIZATION : 1->P , -1->S
966 fresn=0.5*(rp*(1+pdotr)+rs*(1-pdotr));
975 //__________________________________________
976 Float_t AliRICH::AbsoCH4(Float_t x)
979 //KLOSCH,SCH4(9),WL(9),EM(9),ALENGTH(31)
980 Float_t sch4[9] = {.12,.16,.23,.38,.86,2.8,7.9,28.,80.}; //MB X 10^22
981 //Float_t wl[9] = {153.,152.,151.,150.,149.,148.,147.,146.,145};
982 Float_t em[9] = {8.1,8.158,8.212,8.267,8.322,8.378,8.435,8.493,8.55};
983 const Float_t kLosch=2.686763E19; // LOSCHMIDT NUMBER IN CM-3
984 const Float_t kIgas1=100, kIgas2=0, kOxy=10., kWater=5., kPressure=750.,kTemperature=283.;
985 Float_t pn=kPressure/760.;
986 Float_t tn=kTemperature/273.16;
989 // ------- METHANE CROSS SECTION -----------------
990 // ASTROPH. J. 214, L47 (1978)
996 if(x>=7.75 && x<=8.1)
998 Float_t c0=-1.655279e-1;
999 Float_t c1=6.307392e-2;
1000 Float_t c2=-8.011441e-3;
1001 Float_t c3=3.392126e-4;
1002 sm=(c0+c1*x+c2*x*x+c3*x*x*x)*1.e-18;
1008 while (x<=em[j] && x>=em[j+1])
1011 Float_t a=(sch4[j+1]-sch4[j])/(em[j+1]-em[j]);
1012 sm=(sch4[j]+a*(x-em[j]))*1e-22;
1016 Float_t dm=(kIgas1/100.)*(1.-((kOxy+kWater)/1.e6))*kLosch*pn/tn;
1017 Float_t abslm=1./sm/dm;
1019 // ------- ISOBUTHANE CROSS SECTION --------------
1020 // i-C4H10 (ai) abs. length from curves in
1021 // Lu-McDonald paper for BARI RICH workshop .
1022 // -----------------------------------------------------------
1031 if(x>=7.25 && x<7.375)
1037 Float_t si = 1./(ai*kLosch*273.16/293.); // ISOB. CRO.SEC.IN CM2
1038 Float_t di=(kIgas2/100.)*(1.-((kOxy+kWater)/1.e6))*kLosch*pn/tn;
1043 // ---------------------------------------------------------
1045 // transmission of O2
1047 // y= path in cm, x=energy in eV
1048 // so= cross section for UV absorption in cm2
1049 // do= O2 molecular density in cm-3
1050 // ---------------------------------------------------------
1058 so=3.392709e-13 * TMath::Exp(2.864104 *x);
1064 so=2.910039e-34 * TMath::Exp(10.3337*x);
1071 Float_t a0=-73770.76;
1072 Float_t a1=46190.69;
1073 Float_t a2=-11475.44;
1074 Float_t a3=1412.611;
1075 Float_t a4=-86.07027;
1076 Float_t a5=2.074234;
1077 so= a0+(a1*x)+(a2*x*x)+(a3*x*x*x)+(a4*x*x*x*x)+(a5*x*x*x*x*x);
1081 Float_t dox=(kOxy/1e6)*kLosch*pn/tn;
1086 // ---------------------------------------------------------
1088 // transmission of H2O
1090 // y= path in cm, x=energy in eV
1091 // sw= cross section for UV absorption in cm2
1092 // dw= H2O molecular density in cm-3
1093 // ---------------------------------------------------------
1097 Float_t b0=29231.65;
1098 Float_t b1=-15807.74;
1099 Float_t b2=3192.926;
1100 Float_t b3=-285.4809;
1101 Float_t b4=9.533944;
1105 Float_t sw= b0+(b1*x)+(b2*x*x)+(b3*x*x*x)+(b4*x*x*x*x);
1107 Float_t dw=(kWater/1e6)*kLosch*pn/tn;
1113 // ---------------------------------------------------------
1115 Float_t alength=1./(1./abslm+1./absli+1./abslo+1./abslw);
1121 //___________________________________________
1122 Int_t AliRICH::DistancetoPrimitive(Int_t , Int_t )
1130 //___________________________________________
1131 void AliRICH::MakeBranch(Option_t* option)
1133 // Create Tree branches for the RICH.
1135 const Int_t kBufferSize = 4000;
1136 char branchname[20];
1139 AliDetector::MakeBranch(option);
1140 sprintf(branchname,"%sCerenkov",GetName());
1141 if (fCerenkovs && gAlice->TreeH()) {
1142 gAlice->TreeH()->Branch(branchname,&fCerenkovs, kBufferSize);
1143 printf("Making Branch %s for Cerenkov Hits\n",branchname);
1146 sprintf(branchname,"%sPadHits",GetName());
1147 if (fPadHits && gAlice->TreeH()) {
1148 gAlice->TreeH()->Branch(branchname,&fPadHits, kBufferSize);
1149 printf("Making Branch %s for PadHits\n",branchname);
1152 // one branch for digits per chamber
1155 for (i=0; i<kNCH ;i++) {
1156 sprintf(branchname,"%sDigits%d",GetName(),i+1);
1158 if (fDchambers && gAlice->TreeD()) {
1159 gAlice->TreeD()->Branch(branchname,&((*fDchambers)[i]), kBufferSize);
1160 printf("Making Branch %s for digits in chamber %d\n",branchname,i+1);
1164 // one branch for raw clusters per chamber
1165 for (i=0; i<kNCH ;i++) {
1166 sprintf(branchname,"%sRawClusters%d",GetName(),i+1);
1168 if (fRawClusters && gAlice->TreeR()) {
1169 gAlice->TreeR()->Branch(branchname,&((*fRawClusters)[i]), kBufferSize);
1170 printf("Making Branch %s for raw clusters in chamber %d\n",branchname,i+1);
1174 // one branch for rec hits per chamber
1175 for (i=0; i<kNCH ;i++) {
1176 sprintf(branchname,"%sRecHits%d",GetName(),i+1);
1178 if (fRecHits && gAlice->TreeR()) {
1179 gAlice->TreeR()->Branch(branchname,&((*fRecHits)[i]), kBufferSize);
1180 printf("Making Branch %s for rec. hits in chamber %d\n",branchname,i+1);
1185 //___________________________________________
1186 void AliRICH::SetTreeAddress()
1188 // Set branch address for the Hits and Digits Tree.
1189 char branchname[20];
1192 AliDetector::SetTreeAddress();
1195 TTree *treeH = gAlice->TreeH();
1196 TTree *treeD = gAlice->TreeD();
1197 TTree *treeR = gAlice->TreeR();
1201 branch = treeH->GetBranch("RICHPadHits");
1202 if (branch) branch->SetAddress(&fPadHits);
1205 branch = treeH->GetBranch("RICHCerenkov");
1206 if (branch) branch->SetAddress(&fCerenkovs);
1211 for (int i=0; i<kNCH; i++) {
1212 sprintf(branchname,"%sDigits%d",GetName(),i+1);
1214 branch = treeD->GetBranch(branchname);
1215 if (branch) branch->SetAddress(&((*fDchambers)[i]));
1220 for (i=0; i<kNCH; i++) {
1221 sprintf(branchname,"%sRawClusters%d",GetName(),i+1);
1223 branch = treeR->GetBranch(branchname);
1224 if (branch) branch->SetAddress(&((*fRawClusters)[i]));
1228 for (i=0; i<kNCH; i++) {
1229 sprintf(branchname,"%sRecHits%d",GetName(),i+1);
1231 branch = treeR->GetBranch(branchname);
1232 if (branch) branch->SetAddress(&((*fRecHits)[i]));
1238 //___________________________________________
1239 void AliRICH::ResetHits()
1241 // Reset number of clusters and the cluster array for this detector
1242 AliDetector::ResetHits();
1245 if (fPadHits) fPadHits->Clear();
1246 if (fCerenkovs) fCerenkovs->Clear();
1250 //____________________________________________
1251 void AliRICH::ResetDigits()
1254 // Reset number of digits and the digits array for this detector
1256 for ( int i=0;i<kNCH;i++ ) {
1257 if ((*fDchambers)[i]) (*fDchambers)[i]->Clear();
1258 if (fNdch) fNdch[i]=0;
1262 //____________________________________________
1263 void AliRICH::ResetRawClusters()
1266 // Reset number of raw clusters and the raw clust array for this detector
1268 for ( int i=0;i<kNCH;i++ ) {
1269 if ((*fRawClusters)[i]) ((TClonesArray*)(*fRawClusters)[i])->Clear();
1270 if (fNrawch) fNrawch[i]=0;
1274 //____________________________________________
1275 void AliRICH::ResetRecHits()
1278 // Reset number of raw clusters and the raw clust array for this detector
1281 for ( int i=0;i<kNCH;i++ ) {
1282 if ((*fRecHits)[i]) ((TClonesArray*)(*fRecHits)[i])->Clear();
1283 if (fNrechits) fNrechits[i]=0;
1287 //___________________________________________
1288 void AliRICH::SetGeometryModel(Int_t id, AliRICHGeometry *geometry)
1292 // Setter for the RICH geometry model
1296 ((AliRICHChamber*) (*fChambers)[id])->GeometryModel(geometry);
1299 //___________________________________________
1300 void AliRICH::SetSegmentationModel(Int_t id, AliRICHSegmentation *segmentation)
1304 // Setter for the RICH segmentation model
1307 ((AliRICHChamber*) (*fChambers)[id])->SegmentationModel(segmentation);
1310 //___________________________________________
1311 void AliRICH::SetResponseModel(Int_t id, AliRICHResponse *response)
1315 // Setter for the RICH response model
1318 ((AliRICHChamber*) (*fChambers)[id])->ResponseModel(response);
1321 void AliRICH::SetReconstructionModel(Int_t id, AliRICHClusterFinder *reconst)
1325 // Setter for the RICH reconstruction model (clusters)
1328 ((AliRICHChamber*) (*fChambers)[id])->ReconstructionModel(reconst);
1331 void AliRICH::SetNsec(Int_t id, Int_t nsec)
1335 // Sets the number of padplanes
1338 ((AliRICHChamber*) (*fChambers)[id])->SetNsec(nsec);
1342 //___________________________________________
1343 void AliRICH::StepManager()
1346 // Full Step Manager
1350 static Int_t vol[2];
1352 static Float_t hits[18];
1353 static Float_t ckovData[19];
1354 TLorentzVector position;
1355 TLorentzVector momentum;
1358 Float_t localPos[3];
1359 Float_t localMom[4];
1360 Float_t localTheta,localPhi;
1362 Float_t destep, step;
1365 Float_t coscerenkov;
1366 static Float_t eloss, xhit, yhit, tlength;
1367 const Float_t kBig=1.e10;
1369 TClonesArray &lhits = *fHits;
1370 TGeant3 *geant3 = (TGeant3*) gMC;
1371 TParticle *current = (TParticle*)(*gAlice->Particles())[gAlice->CurrentTrack()];
1373 //if (current->Energy()>1)
1376 // Only gas gap inside chamber
1377 // Tag chambers and record hits when track enters
1380 id=gMC->CurrentVolID(copy);
1381 Float_t cherenkovLoss=0;
1382 //gAlice->KeepTrack(gAlice->CurrentTrack());
1384 gMC->TrackPosition(position);
1388 bzero((char *)ckovData,sizeof(ckovData)*19);
1389 ckovData[1] = pos[0]; // X-position for hit
1390 ckovData[2] = pos[1]; // Y-position for hit
1391 ckovData[3] = pos[2]; // Z-position for hit
1392 //ckovData[11] = gAlice->CurrentTrack();
1394 //AliRICH *RICH = (AliRICH *) gAlice->GetDetector("RICH");
1396 /********************Store production parameters for Cerenkov photons************************/
1397 //is it a Cerenkov photon?
1398 if (gMC->TrackPid() == 50000050) {
1400 //if (gMC->VolId("GAP ")==gMC->CurrentVolID(copy))
1402 Float_t ckovEnergy = current->Energy();
1403 //energy interval for tracking
1404 if (ckovEnergy > 5.6e-09 && ckovEnergy < 7.8e-09 )
1405 //if (ckovEnergy > 0)
1407 if (gMC->IsTrackEntering()){ //is track entering?
1408 if (gMC->VolId("FRE1")==gMC->CurrentVolID(copy) || gMC->VolId("FRE2")==gMC->CurrentVolID(copy))
1410 if (geant3->Gctrak()->nstep<1){ //is it the first step?
1411 Int_t mother = current->GetFirstMother();
1413 //printf("Second Mother:%d\n",current->GetSecondMother());
1415 ckovData[10] = mother;
1416 ckovData[11] = gAlice->CurrentTrack();
1417 ckovData[12] = 1; //Media where photon was produced 1->Freon, 2->Quarz
1420 //printf("Index: %d\n",fCkovNumber);
1421 } //first step question
1424 if (geant3->Gctrak()->nstep<1){ //is it first step?
1425 if (gMC->VolId("QUAR")==gMC->CurrentVolID(copy)) //is it in quarz?
1429 } //first step question
1431 //printf("Before %d\n",fFreonProd);
1432 } //track entering question
1434 if (ckovData[12] == 1) //was it produced in Freon?
1435 //if (fFreonProd == 1)
1437 if (gMC->IsTrackEntering()){ //is track entering?
1439 if (gMC->VolId("META")==gMC->CurrentVolID(copy)) //is it in gap?
1441 //printf("Got in\n");
1442 gMC->TrackMomentum(momentum);
1447 // Z-position for hit
1450 /**************** Photons lost in second grid have to be calculated by hand************/
1452 Float_t cophi = TMath::Cos(TMath::ATan2(mom[0], mom[1]));
1453 Float_t t = (1. - .025 / cophi) * (1. - .05 / cophi);
1455 //printf("grid calculation:%f\n",t);
1457 //geant3->StopTrack();
1459 AddCerenkov(gAlice->CurrentTrack(),vol,ckovData);
1460 //printf("Lost one in grid\n");
1462 /**********************************************************************************/
1465 if (gMC->VolId("CSI ")==gMC->CurrentVolID(copy)) //is it in csi?
1467 gMC->TrackMomentum(momentum);
1473 /********* Photons lost by Fresnel reflection have to be calculated by hand********/
1474 /***********************Cerenkov phtons (always polarised)*************************/
1476 Float_t cophi = TMath::Cos(TMath::ATan2(mom[0], mom[1]));
1477 Float_t t = Fresnel(ckovEnergy*1e9,cophi,1);
1480 //geant3->StopTrack();
1482 AddCerenkov(gAlice->CurrentTrack(),vol,ckovData);
1483 //printf("Lost by Fresnel\n");
1485 /**********************************************************************************/
1490 /********************Evaluation of losses************************/
1491 /******************still in the old fashion**********************/
1493 Int_t i1 = geant3->Gctrak()->nmec; //number of physics mechanisms acting on the particle
1494 for (Int_t i = 0; i < i1; ++i) {
1496 if (geant3->Gctrak()->lmec[i] == 106) { //was it reflected
1498 if (gMC->VolId("FRE1")==gMC->CurrentVolID(copy) || gMC->VolId("FRE2")==gMC->CurrentVolID(copy))
1500 if (gMC->CurrentVolID(copy) == gMC->VolId("QUAR"))
1502 //geant3->StopTrack();
1503 AddCerenkov(gAlice->CurrentTrack(),vol,ckovData);
1504 } //reflection question
1508 else if (geant3->Gctrak()->lmec[i] == 101) { //was it absorbed?
1510 if (gMC->VolId("FRE1")==gMC->CurrentVolID(copy) || gMC->VolId("FRE2")==gMC->CurrentVolID(copy))
1512 if (gMC->CurrentVolID(copy) == gMC->VolId("QUAR"))
1514 if (gMC->CurrentVolID(copy) == gMC->VolId("META"))
1516 if (gMC->CurrentVolID(copy) == gMC->VolId("GAP "))
1519 if (gMC->CurrentVolID(copy) == gMC->VolId("SRIC"))
1523 if (gMC->CurrentVolID(copy) == gMC->VolId("CSI ")) {
1526 //geant3->StopTrack();
1527 AddCerenkov(gAlice->CurrentTrack(),vol,ckovData);
1528 //printf("Added cerenkov %d\n",fCkovNumber);
1529 } //absorption question
1532 // Photon goes out of tracking scope
1533 else if (geant3->Gctrak()->lmec[i] == 30) { //is it below energy treshold?
1535 //geant3->StopTrack();
1536 AddCerenkov(gAlice->CurrentTrack(),vol,ckovData);
1537 } // energy treshold question
1538 } //number of mechanisms cycle
1539 /**********************End of evaluation************************/
1540 } //freon production question
1541 } //energy interval question
1542 //}//inside the proximity gap question
1543 } //cerenkov photon question
1545 /**************************************End of Production Parameters Storing*********************/
1548 /*******************************Treat photons that hit the CsI (Ckovs and Feedbacks)************/
1550 if (gMC->TrackPid() == 50000050 || gMC->TrackPid() == 50000051) {
1551 //printf("Cerenkov\n");
1552 if (gMC->VolId("CSI ")==gMC->CurrentVolID(copy))
1555 if (gMC->Edep() > 0.){
1556 gMC->TrackPosition(position);
1557 gMC->TrackMomentum(momentum);
1565 Double_t tc = mom[0]*mom[0]+mom[1]*mom[1];
1566 Double_t rt = TMath::Sqrt(tc);
1567 theta = Float_t(TMath::ATan2(rt,Double_t(mom[2])))*kRaddeg;
1568 phi = Float_t(TMath::ATan2(Double_t(mom[1]),Double_t(mom[0])))*kRaddeg;
1569 gMC->Gmtod(pos,localPos,1);
1570 gMC->Gmtod(mom,localMom,2);
1572 gMC->CurrentVolOffID(2,copy);
1576 //Int_t sector=((AliRICHChamber*) (*fChambers)[idvol])
1577 //->Sector(localPos[0], localPos[2]);
1578 //printf("Sector:%d\n",sector);
1580 /*if (gMC->TrackPid() == 50000051){
1582 printf("Feedbacks:%d\n",fFeedbacks);
1585 ((AliRICHChamber*) (*fChambers)[idvol])
1586 ->SigGenInit(localPos[0], localPos[2], localPos[1]);
1588 ckovData[0] = gMC->TrackPid(); // particle type
1589 ckovData[1] = pos[0]; // X-position for hit
1590 ckovData[2] = pos[1]; // Y-position for hit
1591 ckovData[3] = pos[2]; // Z-position for hit
1592 ckovData[4] = theta; // theta angle of incidence
1593 ckovData[5] = phi; // phi angle of incidence
1594 ckovData[8] = (Float_t) fNPadHits; // first padhit
1595 ckovData[9] = -1; // last pad hit
1596 ckovData[13] = 4; // photon was detected
1597 ckovData[14] = mom[0];
1598 ckovData[15] = mom[1];
1599 ckovData[16] = mom[2];
1601 destep = gMC->Edep();
1602 gMC->SetMaxStep(kBig);
1603 cherenkovLoss += destep;
1604 ckovData[7]=cherenkovLoss;
1606 nPads = MakePadHits(localPos[0],localPos[2],cherenkovLoss,idvol,kCerenkov);
1607 if (fNPadHits > (Int_t)ckovData[8]) {
1608 ckovData[8]= ckovData[8]+1;
1609 ckovData[9]= (Float_t) fNPadHits;
1612 ckovData[17] = nPads;
1613 //printf("nPads:%d",nPads);
1615 //TClonesArray *Hits = RICH->Hits();
1616 AliRICHHit *mipHit = (AliRICHHit*) (fHits->UncheckedAt(0));
1619 mom[0] = current->Px();
1620 mom[1] = current->Py();
1621 mom[2] = current->Pz();
1622 Float_t mipPx = mipHit->fMomX;
1623 Float_t mipPy = mipHit->fMomY;
1624 Float_t mipPz = mipHit->fMomZ;
1626 Float_t r = mom[0]*mom[0] + mom[1]*mom[1] + mom[2]*mom[2];
1627 Float_t rt = TMath::Sqrt(r);
1628 Float_t mipR = mipPx*mipPx + mipPy*mipPy + mipPz*mipPz;
1629 Float_t mipRt = TMath::Sqrt(mipR);
1632 coscerenkov = (mom[0]*mipPx + mom[1]*mipPy + mom[2]*mipPz)/(rt*mipRt);
1638 Float_t cherenkov = TMath::ACos(coscerenkov);
1639 ckovData[18]=cherenkov;
1643 AddHit(gAlice->CurrentTrack(),vol,ckovData);
1644 AddCerenkov(gAlice->CurrentTrack(),vol,ckovData);
1651 /***********************************************End of photon hits*********************************************/
1654 /**********************************************Charged particles treatment*************************************/
1656 else if (gMC->TrackCharge())
1660 /*if (gMC->IsTrackEntering())
1662 hits[13]=20;//is track entering?
1664 if (gMC->VolId("FRE1")==gMC->CurrentVolID(copy) || gMC->VolId("FRE2")==gMC->CurrentVolID(copy))
1669 if (gMC->VolId("GAP ")== gMC->CurrentVolID(copy)) {
1670 // Get current particle id (ipart), track position (pos) and momentum (mom)
1672 gMC->CurrentVolOffID(3,copy);
1676 //Int_t sector=((AliRICHChamber*) (*fChambers)[idvol])
1677 //->Sector(localPos[0], localPos[2]);
1678 //printf("Sector:%d\n",sector);
1680 gMC->TrackPosition(position);
1681 gMC->TrackMomentum(momentum);
1689 gMC->Gmtod(pos,localPos,1);
1690 gMC->Gmtod(mom,localMom,2);
1692 ipart = gMC->TrackPid();
1694 // momentum loss and steplength in last step
1695 destep = gMC->Edep();
1696 step = gMC->TrackStep();
1699 // record hits when track enters ...
1700 if( gMC->IsTrackEntering()) {
1701 // gMC->SetMaxStep(fMaxStepGas);
1702 Double_t tc = mom[0]*mom[0]+mom[1]*mom[1];
1703 Double_t rt = TMath::Sqrt(tc);
1704 theta = Float_t(TMath::ATan2(rt,Double_t(mom[2])))*kRaddeg;
1705 phi = Float_t(TMath::ATan2(Double_t(mom[1]),Double_t(mom[0])))*kRaddeg;
1708 Double_t localTc = localMom[0]*localMom[0]+localMom[2]*localMom[2];
1709 Double_t localRt = TMath::Sqrt(localTc);
1710 localTheta = Float_t(TMath::ATan2(localRt,Double_t(localMom[1])))*kRaddeg;
1711 localPhi = Float_t(TMath::ATan2(Double_t(localMom[2]),Double_t(localMom[0])))*kRaddeg;
1713 hits[0] = Float_t(ipart); // particle type
1714 hits[1] = localPos[0]; // X-position for hit
1715 hits[2] = localPos[1]; // Y-position for hit
1716 hits[3] = localPos[2]; // Z-position for hit
1717 hits[4] = localTheta; // theta angle of incidence
1718 hits[5] = localPhi; // phi angle of incidence
1719 hits[8] = (Float_t) fNPadHits; // first padhit
1720 hits[9] = -1; // last pad hit
1721 hits[13] = fFreonProd; // did id hit the freon?
1730 Chamber(idvol).LocaltoGlobal(localPos,hits+1);
1733 //To make chamber coordinates x-y had to pass localPos[0], localPos[2]
1736 // Only if not trigger chamber
1739 // Initialize hit position (cursor) in the segmentation model
1740 ((AliRICHChamber*) (*fChambers)[idvol])
1741 ->SigGenInit(localPos[0], localPos[2], localPos[1]);
1746 // Calculate the charge induced on a pad (disintegration) in case
1748 // Mip left chamber ...
1749 if( gMC->IsTrackExiting() || gMC->IsTrackStop() || gMC->IsTrackDisappeared()){
1750 gMC->SetMaxStep(kBig);
1755 // Only if not trigger chamber
1759 if(gMC->TrackPid() == kNeutron)
1760 printf("\n\n\n\n\n Neutron Making Pad Hit!!! \n\n\n\n");
1761 nPads = MakePadHits(xhit,yhit,eloss,idvol,kMip);
1763 //printf("nPads:%d",nPads);
1769 if (fNPadHits > (Int_t)hits[8]) {
1771 hits[9]= (Float_t) fNPadHits;
1775 new(lhits[fNhits++]) AliRICHHit(fIshunt,gAlice->CurrentTrack(),vol,hits);
1778 // Check additional signal generation conditions
1779 // defined by the segmentation
1780 // model (boundary crossing conditions)
1782 (((AliRICHChamber*) (*fChambers)[idvol])
1783 ->SigGenCond(localPos[0], localPos[2], localPos[1]))
1785 ((AliRICHChamber*) (*fChambers)[idvol])
1786 ->SigGenInit(localPos[0], localPos[2], localPos[1]);
1789 if(gMC->TrackPid() == kNeutron)
1790 printf("\n\n\n\n\n Neutron Making Pad Hit!!! \n\n\n\n");
1791 nPads = MakePadHits(xhit,yhit,eloss,idvol,kMip);
1793 //printf("Npads:%d",NPads);
1800 // nothing special happened, add up energy loss
1807 /*************************************************End of MIP treatment**************************************/
1811 void AliRICH::FindClusters(Int_t nev,Int_t lastEntry)
1815 // Loop on chambers and on cathode planes
1817 for (Int_t icat=1;icat<2;icat++) {
1818 gAlice->ResetDigits();
1819 gAlice->TreeD()->GetEvent(1); // spurious +1 ...
1820 for (Int_t ich=0;ich<kNCH;ich++) {
1821 AliRICHChamber* iChamber=(AliRICHChamber*) (*fChambers)[ich];
1822 TClonesArray *pRICHdigits = this->DigitsAddress(ich);
1823 if (pRICHdigits == 0)
1826 // Get ready the current chamber stuff
1828 AliRICHResponse* response = iChamber->GetResponseModel();
1829 AliRICHSegmentation* seg = iChamber->GetSegmentationModel();
1830 AliRICHClusterFinder* rec = iChamber->GetReconstructionModel();
1832 rec->SetSegmentation(seg);
1833 rec->SetResponse(response);
1834 rec->SetDigits(pRICHdigits);
1835 rec->SetChamber(ich);
1836 if (nev==0) rec->CalibrateCOG();
1837 rec->FindRawClusters();
1840 fRch=RawClustAddress(ich);
1844 gAlice->TreeR()->Fill();
1846 for (int i=0;i<kNCH;i++) {
1847 fRch=RawClustAddress(i);
1848 int nraw=fRch->GetEntriesFast();
1849 printf ("Chamber %d, raw clusters %d\n",i,nraw);
1857 sprintf(hname,"TreeR%d",nev);
1858 gAlice->TreeR()->Write(hname,kOverwrite,0);
1859 gAlice->TreeR()->Reset();
1861 //gObjectTable->Print();
1865 //______________________________________________________________________________
1866 void AliRICH::Streamer(TBuffer &R__b)
1868 // Stream an object of class AliRICH.
1869 AliRICHChamber *iChamber;
1870 AliRICHSegmentation *segmentation;
1871 AliRICHResponse *response;
1872 TClonesArray *digitsaddress;
1873 TClonesArray *rawcladdress;
1874 TClonesArray *rechitaddress;
1876 if (R__b.IsReading()) {
1877 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
1878 AliDetector::Streamer(R__b);
1880 R__b >> fPadHits; // diff
1881 R__b >> fNcerenkovs;
1882 R__b >> fCerenkovs; // diff
1884 R__b >> fRawClusters;
1885 R__b >> fRecHits; //diff
1886 R__b >> fDebugLevel; //diff
1887 R__b.ReadArray(fNdch);
1888 R__b.ReadArray(fNrawch);
1889 R__b.ReadArray(fNrechits);
1892 // Stream chamber related information
1893 for (Int_t i =0; i<kNCH; i++) {
1894 iChamber=(AliRICHChamber*) (*fChambers)[i];
1895 iChamber->Streamer(R__b);
1896 segmentation=iChamber->GetSegmentationModel();
1897 segmentation->Streamer(R__b);
1898 response=iChamber->GetResponseModel();
1899 response->Streamer(R__b);
1900 rawcladdress=(TClonesArray*) (*fRawClusters)[i];
1901 rawcladdress->Streamer(R__b);
1902 rechitaddress=(TClonesArray*) (*fRecHits)[i];
1903 rechitaddress->Streamer(R__b);
1904 digitsaddress=(TClonesArray*) (*fDchambers)[i];
1905 digitsaddress->Streamer(R__b);
1907 R__b >> fDebugLevel;
1908 R__b >> fCkovNumber;
1915 R__b >> fLostAquarz;
1923 R__b >> fLostFresnel;
1926 R__b.WriteVersion(AliRICH::IsA());
1927 AliDetector::Streamer(R__b);
1929 R__b << fPadHits; // diff
1930 R__b << fNcerenkovs;
1931 R__b << fCerenkovs; // diff
1933 R__b << fRawClusters;
1934 R__b << fRecHits; //diff
1935 R__b << fDebugLevel; //diff
1936 R__b.WriteArray(fNdch, kNCH);
1937 R__b.WriteArray(fNrawch, kNCH);
1938 R__b.WriteArray(fNrechits, kNCH);
1941 // Stream chamber related information
1942 for (Int_t i =0; i<kNCH; i++) {
1943 iChamber=(AliRICHChamber*) (*fChambers)[i];
1944 iChamber->Streamer(R__b);
1945 segmentation=iChamber->GetSegmentationModel();
1946 segmentation->Streamer(R__b);
1947 response=iChamber->GetResponseModel();
1948 response->Streamer(R__b);
1949 rawcladdress=(TClonesArray*) (*fRawClusters)[i];
1950 rawcladdress->Streamer(R__b);
1951 rechitaddress=(TClonesArray*) (*fRecHits)[i];
1952 rechitaddress->Streamer(R__b);
1953 digitsaddress=(TClonesArray*) (*fDchambers)[i];
1954 digitsaddress->Streamer(R__b);
1956 R__b << fDebugLevel;
1957 R__b << fCkovNumber;
1964 R__b << fLostAquarz;
1972 R__b << fLostFresnel;
1975 AliRICHPadHit* AliRICH::FirstPad(AliRICHHit* hit,TClonesArray *clusters )
1978 // Initialise the pad iterator
1979 // Return the address of the first padhit for hit
1980 TClonesArray *theClusters = clusters;
1981 Int_t nclust = theClusters->GetEntriesFast();
1982 if (nclust && hit->fPHlast > 0) {
1983 sMaxIterPad=Int_t(hit->fPHlast);
1984 sCurIterPad=Int_t(hit->fPHfirst);
1985 return (AliRICHPadHit*) clusters->UncheckedAt(sCurIterPad-1);
1992 AliRICHPadHit* AliRICH::NextPad(TClonesArray *clusters)
1995 // Iterates over pads
1998 if (sCurIterPad <= sMaxIterPad) {
1999 return (AliRICHPadHit*) clusters->UncheckedAt(sCurIterPad-1);
2006 void AliRICH::Digitise(Int_t nev, Int_t flag, Option_t *option,Text_t *filename)
2008 // keep galice.root for signal and name differently the file for
2009 // background when add! otherwise the track info for signal will be lost !
2011 static Bool_t first=kTRUE;
2012 static TFile *pFile;
2013 char *addBackground = strstr(option,"Add");
2015 FILE* points; //these will be the digits...
2017 points=fopen("points.dat","w");
2019 AliRICHChamber* iChamber;
2020 AliRICHSegmentation* segmentation;
2025 TObjArray *list=new TObjArray;
2026 static TClonesArray *pAddress=0;
2027 if(!pAddress) pAddress=new TClonesArray("TVector",1000);
2030 AliRICH *pRICH = (AliRICH *) gAlice->GetDetector("RICH");
2031 AliRICHHitMap* pHitMap[10];
2033 for (i=0; i<10; i++) {pHitMap[i]=0;}
2034 if (addBackground ) {
2037 cout<<"filename"<<fFileName<<endl;
2038 pFile=new TFile(fFileName);
2039 cout<<"I have opened "<<fFileName<<" file "<<endl;
2040 fHits2 = new TClonesArray("AliRICHHit",1000 );
2041 fClusters2 = new TClonesArray("AliRICHPadHit",10000);
2045 // Get Hits Tree header from file
2046 if(fHits2) fHits2->Clear();
2047 if(fClusters2) fClusters2->Clear();
2048 if(TrH1) delete TrH1;
2052 sprintf(treeName,"TreeH%d",nev);
2053 TrH1 = (TTree*)gDirectory->Get(treeName);
2055 printf("ERROR: cannot find Hits Tree for event:%d\n",nev);
2057 // Set branch addresses
2059 char branchname[20];
2060 sprintf(branchname,"%s",GetName());
2061 if (TrH1 && fHits2) {
2062 branch = TrH1->GetBranch(branchname);
2063 if (branch) branch->SetAddress(&fHits2);
2065 if (TrH1 && fClusters2) {
2066 branch = TrH1->GetBranch("RICHCluster");
2067 if (branch) branch->SetAddress(&fClusters2);
2074 for (i =0; i<kNCH; i++) {
2075 iChamber=(AliRICHChamber*) (*fChambers)[i];
2076 segmentation=iChamber->GetSegmentationModel(1);
2077 pHitMap[i] = new AliRICHHitMapA1(segmentation, list);
2083 TTree *treeH = gAlice->TreeH();
2084 Int_t ntracks =(Int_t) treeH->GetEntries();
2085 for (Int_t track=0; track<ntracks; track++) {
2086 gAlice->ResetHits();
2087 treeH->GetEvent(track);
2090 for(AliRICHHit* mHit=(AliRICHHit*)pRICH->FirstHit(-1);
2092 mHit=(AliRICHHit*)pRICH->NextHit())
2097 Int_t nch = mHit->fChamber-1; // chamber number
2098 if (nch >kNCH) continue;
2099 iChamber = &(pRICH->Chamber(nch));
2101 TParticle *current = (TParticle*)(*gAlice->Particles())[track];
2103 Int_t particle = current->GetPdgCode();
2105 //printf("Flag:%d\n",flag);
2106 //printf("Track:%d\n",track);
2107 //printf("Particle:%d\n",particle);
2113 if(TMath::Abs(particle) == 211 || TMath::Abs(particle) == 111)
2117 if(TMath::Abs(particle)==321 || TMath::Abs(particle)==130 || TMath::Abs(particle)==310
2118 || TMath::Abs(particle)==311)
2121 if (flag == 3 && TMath::Abs(particle)==2212)
2124 if (flag == 4 && TMath::Abs(particle)==13)
2127 if (flag == 5 && TMath::Abs(particle)==11)
2130 if (flag == 6 && TMath::Abs(particle)==2112)
2134 //printf ("Particle: %d, Flag: %d, Digitse: %d\n",particle,flag,digitse);
2141 // Loop over pad hits
2142 for (AliRICHPadHit* mPad=
2143 (AliRICHPadHit*)pRICH->FirstPad(mHit,fPadHits);
2145 mPad=(AliRICHPadHit*)pRICH->NextPad(fPadHits))
2147 Int_t cathode = mPad->fCathode; // cathode number
2148 Int_t ipx = mPad->fPadX; // pad number on X
2149 Int_t ipy = mPad->fPadY; // pad number on Y
2150 Int_t iqpad = mPad->fQpad; // charge per pad
2153 //printf("X:%d, Y:%d, Q:%d\n",ipx,ipy,iqpad);
2156 segmentation=iChamber->GetSegmentationModel(cathode);
2157 segmentation->GetPadCxy(ipx,ipy,thex,they);
2158 new((*pAddress)[countadr++]) TVector(2);
2159 TVector &trinfo=*((TVector*) (*pAddress)[countadr-1]);
2160 trinfo(0)=(Float_t)track;
2161 trinfo(1)=(Float_t)iqpad;
2167 AliRICHTransientDigit* pdigit;
2168 // build the list of fired pads and update the info
2169 if (!pHitMap[nch]->TestHit(ipx, ipy)) {
2170 list->AddAtAndExpand(new AliRICHTransientDigit(nch,digits),counter);
2171 pHitMap[nch]->SetHit(ipx, ipy, counter);
2173 pdigit=(AliRICHTransientDigit*)list->At(list->GetLast());
2175 TObjArray *trlist=(TObjArray*)pdigit->TrackList();
2176 trlist->Add(&trinfo);
2178 pdigit=(AliRICHTransientDigit*) pHitMap[nch]->GetHit(ipx, ipy);
2180 (*pdigit).fSignal+=iqpad;
2181 // update list of tracks
2182 TObjArray* trlist=(TObjArray*)pdigit->TrackList();
2183 Int_t lastEntry=trlist->GetLast();
2184 TVector *ptrkP=(TVector*)trlist->At(lastEntry);
2185 TVector &ptrk=*ptrkP;
2186 Int_t lastTrack=Int_t(ptrk(0));
2187 Int_t lastCharge=Int_t(ptrk(1));
2188 if (lastTrack==track) {
2190 trlist->RemoveAt(lastEntry);
2191 trinfo(0)=lastTrack;
2192 trinfo(1)=lastCharge;
2193 trlist->AddAt(&trinfo,lastEntry);
2195 trlist->Add(&trinfo);
2197 // check the track list
2198 Int_t nptracks=trlist->GetEntriesFast();
2200 printf("Attention - tracks: %d (>2)\n",nptracks);
2201 //printf("cat,nch,ix,iy %d %d %d %d \n",icat+1,nch,ipx,ipy);
2202 for (Int_t tr=0;tr<nptracks;tr++) {
2203 TVector *pptrkP=(TVector*)trlist->At(tr);
2204 TVector &pptrk=*pptrkP;
2205 trk[tr]=Int_t(pptrk(0));
2206 chtrk[tr]=Int_t(pptrk(1));
2208 } // end if nptracks
2210 } //end loop over clusters
2211 }// track type condition
2215 // open the file with background
2217 if (addBackground ) {
2218 ntracks =(Int_t)TrH1->GetEntries();
2219 //printf("background - icat,ntracks1 %d %d\n",icat,ntracks);
2220 //printf("background - Start loop over tracks \n");
2224 for (Int_t trak=0; trak<ntracks; trak++) {
2225 if (fHits2) fHits2->Clear();
2226 if (fClusters2) fClusters2->Clear();
2227 TrH1->GetEvent(trak);
2231 for(int j=0;j<fHits2->GetEntriesFast();++j)
2233 mHit=(AliRICHHit*) (*fHits2)[j];
2234 Int_t nch = mHit->fChamber-1; // chamber number
2235 if (nch >6) continue;
2236 iChamber = &(pRICH->Chamber(nch));
2237 Int_t rmin = (Int_t)iChamber->RInner();
2238 Int_t rmax = (Int_t)iChamber->ROuter();
2240 // Loop over pad hits
2241 for (AliRICHPadHit* mPad=
2242 (AliRICHPadHit*)pRICH->FirstPad(mHit,fClusters2);
2244 mPad=(AliRICHPadHit*)pRICH->NextPad(fClusters2))
2246 Int_t cathode = mPad->fCathode; // cathode number
2247 Int_t ipx = mPad->fPadX; // pad number on X
2248 Int_t ipy = mPad->fPadY; // pad number on Y
2249 Int_t iqpad = mPad->fQpad; // charge per pad
2252 segmentation=iChamber->GetSegmentationModel(cathode);
2253 segmentation->GetPadCxy(ipx,ipy,thex,they);
2254 Float_t rpad=TMath::Sqrt(thex*thex+they*they);
2255 if (rpad < rmin || iqpad ==0 || rpad > rmax) continue;
2256 new((*pAddress)[countadr++]) TVector(2);
2257 TVector &trinfo=*((TVector*) (*pAddress)[countadr-1]);
2258 trinfo(0)=-1; // tag background
2263 if (trak <4 && nch==0)
2264 printf("bgr - pHitMap[nch]->TestHit(ipx, ipy),trak %d %d\n",
2265 pHitMap[nch]->TestHit(ipx, ipy),trak);
2266 AliRICHTransientDigit* pdigit;
2267 // build the list of fired pads and update the info
2268 if (!pHitMap[nch]->TestHit(ipx, ipy)) {
2269 list->AddAtAndExpand(new AliRICHTransientDigit(nch,digits),counter);
2271 pHitMap[nch]->SetHit(ipx, ipy, counter);
2273 printf("bgr new elem in list - counter %d\n",counter);
2275 pdigit=(AliRICHTransientDigit*)list->At(list->GetLast());
2277 TObjArray *trlist=(TObjArray*)pdigit->TrackList();
2278 trlist->Add(&trinfo);
2280 pdigit=(AliRICHTransientDigit*) pHitMap[nch]->GetHit(ipx, ipy);
2282 (*pdigit).fSignal+=iqpad;
2283 // update list of tracks
2284 TObjArray* trlist=(TObjArray*)pdigit->TrackList();
2285 Int_t lastEntry=trlist->GetLast();
2286 TVector *ptrkP=(TVector*)trlist->At(lastEntry);
2287 TVector &ptrk=*ptrkP;
2288 Int_t lastTrack=Int_t(ptrk(0));
2289 if (lastTrack==-1) {
2292 trlist->Add(&trinfo);
2294 // check the track list
2295 Int_t nptracks=trlist->GetEntriesFast();
2297 for (Int_t tr=0;tr<nptracks;tr++) {
2298 TVector *pptrkP=(TVector*)trlist->At(tr);
2299 TVector &pptrk=*pptrkP;
2300 trk[tr]=Int_t(pptrk(0));
2301 chtrk[tr]=Int_t(pptrk(1));
2303 } // end if nptracks
2305 } //end loop over clusters
2308 TTree *fAli=gAlice->TreeK();
2309 if (fAli) pFile =fAli->GetCurrentFile();
2315 //cout<<"Start filling digits \n "<<endl;
2316 Int_t nentries=list->GetEntriesFast();
2317 //printf(" \n \n nentries %d \n",nentries);
2319 // start filling the digits
2321 for (Int_t nent=0;nent<nentries;nent++) {
2322 AliRICHTransientDigit *address=(AliRICHTransientDigit*)list->At(nent);
2323 if (address==0) continue;
2325 Int_t ich=address->fChamber;
2326 Int_t q=address->fSignal;
2327 iChamber=(AliRICHChamber*) (*fChambers)[ich];
2328 AliRICHResponse * response=iChamber->GetResponseModel();
2329 Int_t adcmax= (Int_t) response->MaxAdc();
2332 // add white noise and do zero-suppression and signal truncation (new electronics,old electronics gaus 1.2,0.2)
2333 //printf("Treshold: %d\n",iChamber->fTresh->GetHitIndex(address->fPadX,address->fPadY));
2334 Int_t pedestal = iChamber->fTresh->GetHitIndex(address->fPadX,address->fPadY);
2336 //printf("Pedestal:%d\n",pedestal);
2338 Float_t treshold = (pedestal + 4*1.7);
2340 Float_t meanNoise = gRandom->Gaus(1.7, 0.25);
2341 Float_t noise = gRandom->Gaus(0, meanNoise);
2342 q+=(Int_t)(noise + pedestal);
2343 //q+=(Int_t)(noise);
2344 // magic number to be parametrised !!!
2351 if ( q >= adcmax) q=adcmax;
2352 digits[0]=address->fPadX;
2353 digits[1]=address->fPadY;
2356 TObjArray* trlist=(TObjArray*)address->TrackList();
2357 Int_t nptracks=trlist->GetEntriesFast();
2359 // this was changed to accomodate the real number of tracks
2360 if (nptracks > 10) {
2361 cout<<"Attention - tracks > 10 "<<nptracks<<endl;
2365 printf("Attention - tracks > 2 %d \n",nptracks);
2366 //printf("cat,ich,ix,iy,q %d %d %d %d %d \n",
2367 //icat,ich,digits[0],digits[1],q);
2369 for (Int_t tr=0;tr<nptracks;tr++) {
2370 TVector *ppP=(TVector*)trlist->At(tr);
2372 tracks[tr]=Int_t(pp(0));
2373 charges[tr]=Int_t(pp(1));
2374 } //end loop over list of tracks for one pad
2375 if (nptracks < 10 ) {
2376 for (Int_t t=nptracks; t<10; t++) {
2383 fprintf(points,"%4d, %4d, %4d\n",digits[0],digits[1],digits[2]);
2386 pRICH->AddDigits(ich,tracks,charges,digits);
2388 gAlice->TreeD()->Fill();
2391 for(Int_t ii=0;ii<kNCH;++ii) {
2399 //TTree *TD=gAlice->TreeD();
2400 //Stat_t ndig=TD->GetEntries();
2401 //cout<<"number of digits "<<ndig<<endl;
2403 for (int k=0;k<kNCH;k++) {
2404 fDch= pRICH->DigitsAddress(k);
2405 int ndigit=fDch->GetEntriesFast();
2406 printf ("Chamber %d digits %d \n",k,ndigit);
2408 pRICH->ResetDigits();
2410 sprintf(hname,"TreeD%d",nev);
2411 gAlice->TreeD()->Write(hname,kOverwrite,0);
2414 // gAlice->TreeD()->Reset();
2417 // gObjectTable->Print();
2420 AliRICH& AliRICH::operator=(const AliRICH& rhs)
2422 // Assignment operator
2427 Int_t AliRICH::MakePadHits(Float_t xhit,Float_t yhit,Float_t eloss, Int_t idvol, ResponseType res)
2430 // Calls the charge disintegration method of the current chamber and adds
2431 // the simulated cluster to the root treee
2434 Float_t newclust[6][500];
2438 // Integrated pulse height on chamber
2442 ((AliRICHChamber*) (*fChambers)[idvol])->DisIntegration(eloss, xhit, yhit, nnew, newclust, res);
2447 for (Int_t i=0; i<nnew; i++) {
2448 if (Int_t(newclust[3][i]) > 0) {
2451 clhits[1] = Int_t(newclust[5][i]);
2453 clhits[2] = Int_t(newclust[0][i]);
2455 clhits[3] = Int_t(newclust[1][i]);
2457 clhits[4] = Int_t(newclust[2][i]);
2459 clhits[5] = Int_t(newclust[3][i]);
2460 // Pad: chamber sector
2461 clhits[6] = Int_t(newclust[4][i]);