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 **************************************************************************/
15 //-----------------------------------------------------//
17 // Source File : PMDDigitizer.cxx, Version 00 //
19 // Date : September 20 2002 //
21 //-----------------------------------------------------//
23 #include <Riostream.h>
27 #include <TGeometry.h>
28 #include <TObjArray.h>
29 #include <TClonesArray.h>
32 #include <TParticle.h>
38 #include "AliDetector.h"
39 #include "AliRunLoader.h"
40 #include "AliLoader.h"
41 #include "AliConfig.h"
43 #include "AliRunDigitizer.h"
44 #include "AliDigitizer.h"
45 #include "AliHeader.h"
46 #include "AliCDBManager.h"
47 #include "AliCDBStorage.h"
48 #include "AliCDBEntry.h"
52 #include "AliPMDhit.h"
53 #include "AliPMDcell.h"
54 #include "AliPMDsdigit.h"
55 #include "AliPMDdigit.h"
56 #include "AliPMDCalibData.h"
57 #include "AliPMDPedestal.h"
58 #include "AliPMDDigitizer.h"
61 ClassImp(AliPMDDigitizer)
63 AliPMDDigitizer::AliPMDDigitizer() :
68 fCalibGain(GetCalibGain()),
69 fCalibPed(GetCalibPed()),
77 fZPos(361.5) // in units of cm, default position of PMD
79 // Default Constructor
81 for (Int_t i = 0; i < fgkTotUM; i++)
83 for (Int_t j = 0; j < fgkRow; j++)
85 for (Int_t k = 0; k < fgkCol; k++)
89 fCPVCounter[i][j][k] = 0;
90 fPRECounter[i][j][k] = 0;
91 fCPVTrackNo[i][j][k] = -1;
92 fPRETrackNo[i][j][k] = -1;
93 fCPVTrackPid[i][j][k] = -1;
94 fPRETrackPid[i][j][k] = -1;
100 //____________________________________________________________________________
101 AliPMDDigitizer::AliPMDDigitizer(const AliPMDDigitizer& digitizer):
102 AliDigitizer(digitizer),
107 fCalibGain(GetCalibGain()),
108 fCalibPed(GetCalibPed()),
116 fZPos(361.5) // in units of cm, default position of PMD
119 AliError("Copy constructor not allowed ");
122 //____________________________________________________________________________
123 AliPMDDigitizer & AliPMDDigitizer::operator=(const AliPMDDigitizer& /*digitizer*/)
125 // Assignment operator
126 AliError("Assignement operator not allowed ");
130 //____________________________________________________________________________
131 AliPMDDigitizer::AliPMDDigitizer(AliRunDigitizer* manager):
132 AliDigitizer(manager),
137 fCalibGain(GetCalibGain()),
138 fCalibPed(GetCalibPed()),
139 fSDigits(new TClonesArray("AliPMDsdigit", 1000)),
140 fDigits(new TClonesArray("AliPMDdigit", 1000)),
146 fZPos(361.5)// in units of cm, This is the default position of PMD
148 // ctor which should be used
151 for (Int_t i = 0; i < fgkTotUM; i++)
153 for (Int_t j = 0; j < fgkRow; j++)
155 for (Int_t k = 0; k < fgkCol; k++)
159 fCPVCounter[i][j][k] = 0;
160 fPRECounter[i][j][k] = 0;
161 fCPVTrackNo[i][j][k] = -1;
162 fPRETrackNo[i][j][k] = -1;
163 fCPVTrackPid[i][j][k] = -1;
164 fPRETrackPid[i][j][k] = -1;
170 //____________________________________________________________________________
171 AliPMDDigitizer::~AliPMDDigitizer()
173 // Default Destructor
191 //____________________________________________________________________________
192 void AliPMDDigitizer::OpengAliceFile(const char *file, Option_t *option)
194 // Loads galice.root file and corresponding header, kinematics
195 // hits and sdigits or digits depending on the option
198 TString evfoldname = AliConfig::GetDefaultEventFolderName();
199 fRunLoader = AliRunLoader::GetRunLoader(evfoldname);
201 fRunLoader = AliRunLoader::Open(file,AliConfig::GetDefaultEventFolderName(), "UPDATE");
205 AliError(Form("Can not open session for file %s.",file));
208 const char *cHS = strstr(option,"HS");
209 const char *cHD = strstr(option,"HD");
210 const char *cSD = strstr(option,"SD");
214 if (!fRunLoader->GetAliRun()) fRunLoader->LoadgAlice();
215 if (!fRunLoader->TreeE()) fRunLoader->LoadHeader();
216 if (!fRunLoader->TreeK()) fRunLoader->LoadKinematics();
218 gAlice = fRunLoader->GetAliRun();
222 AliDebug(1,"Alirun object found");
226 AliError("Could not found Alirun object");
229 fPMD = (AliPMD*)gAlice->GetDetector("PMD");
232 fPMDLoader = fRunLoader->GetLoader("PMDLoader");
233 if (fPMDLoader == 0x0)
235 AliError("Can not find PMDLoader");
241 fPMDLoader->LoadHits("READ");
242 fPMDLoader->LoadSDigits("recreate");
246 fPMDLoader->LoadHits("READ");
247 fPMDLoader->LoadDigits("recreate");
251 fPMDLoader->LoadSDigits("READ");
252 fPMDLoader->LoadDigits("recreate");
255 //____________________________________________________________________________
256 void AliPMDDigitizer::Hits2SDigits(Int_t ievt)
258 // This reads the PMD Hits tree and assigns the right track number
259 // to a cell and stores in the summable digits tree
262 const Int_t kPi0 = 111;
263 const Int_t kGamma = 22;
271 Float_t xPos, yPos, zPos;
272 Int_t xpad = -1, ypad = -1;
274 Float_t vx = -999.0, vy = -999.0, vz = -999.0;
277 if (!fSDigits) fSDigits = new TClonesArray("AliPMDsdigit", 1000);
280 AliDebug(1,Form("Event Number = %d",ievt));
281 Int_t nparticles = fRunLoader->GetHeader()->GetNtrack();
282 AliDebug(1,Form("Number of Particles = %d",nparticles));
286 fRunLoader->GetEvent(ievt);
287 // ------------------------------------------------------- //
288 // Pointer to specific detector hits.
289 // Get pointers to Alice detectors and Hits containers
291 TTree* treeH = fPMDLoader->TreeH();
293 Int_t ntracks = (Int_t) treeH->GetEntries();
294 AliDebug(1,Form("Number of Tracks in the TreeH = %d", ntracks));
295 TTree* treeS = fPMDLoader->TreeS();
298 fPMDLoader->MakeTree("S");
299 treeS = fPMDLoader->TreeS();
301 Int_t bufsize = 16000;
302 treeS->Branch("PMDSDigit", &fSDigits, bufsize);
304 TClonesArray* hits = 0;
305 if (fPMD) hits = fPMD->Hits();
307 // Start loop on tracks in the hits containers
309 for (Int_t track=0; track<ntracks;track++)
311 gAlice->GetMCApp()->ResetHits();
312 treeH->GetEvent(track);
315 npmd = hits->GetEntriesFast();
316 for (Int_t ipmd = 0; ipmd < npmd; ipmd++)
318 fPMDHit = (AliPMDhit*) hits->UncheckedAt(ipmd);
319 trackno = fPMDHit->GetTrack();
320 // get kinematics of the particles
322 TParticle* mparticle = gAlice->GetMCApp()->Particle(trackno);
323 trackpid = mparticle->GetPdgCode();
324 Int_t ks = mparticle->GetStatusCode();
326 Int_t tracknoOld=0, trackpidOld=0, statusOld = 0;
328 if (mparticle->GetFirstMother() == -1)
330 tracknoOld = trackno;
331 trackpidOld = trackpid;
336 Int_t trnotemp = trackno; // Modified on 25th Nov 2009
337 if(ks==1||(imo = mparticle->GetFirstMother())<0 ){
338 vx = mparticle->Vx();
339 vy = mparticle->Vy();
340 vz = mparticle->Vz();
342 if(trackpid==kGamma||trackpid==11||trackpid==-11||
343 trackpid==kPi0)igstatus=1;
347 while(((imo = mparticle->GetFirstMother()) >= 0)&&
348 (ks = mparticle->GetStatusCode() <1) )
350 mparticle = gAlice->GetMCApp()->Particle(imo);
351 trackpid = mparticle->GetPdgCode();
352 ks = mparticle->GetStatusCode();
353 vx = mparticle->Vx();
354 vy = mparticle->Vy();
355 vz = mparticle->Vz();
357 // Modified on 25th Nov 2009
368 // end of modification on 25th Nov 2009
371 if(trackpid==kGamma||trackpid==11||trackpid==-11||
372 trackpid==kPi0)igstatus=1;
375 trackpid=trackpidOld;
378 //-----------------end of modification----------------
379 Float_t ptime = fPMDHit->GetTime()*1e6; // time in microsec
380 if (ptime < 0. || ptime > 1.2) continue;
386 edep = fPMDHit->GetEnergy();
387 Int_t vol1 = fPMDHit->GetVolume(1); // Column
388 Int_t vol2 = fPMDHit->GetVolume(2); // Row
389 Int_t vol7 = fPMDHit->GetVolume(4); // Serial Module No
392 // -----------------------------------------//
393 // In new geometry after adding electronics //
394 // For Super Module 1 & 2 //
395 // nrow = 48, ncol = 96 //
396 // For Super Module 3 & 4 //
397 // nrow = 96, ncol = 48 //
398 // -----------------------------------------//
406 smnumber = vol7 - 24;
408 Int_t vol8 = smnumber/6 + 1; // fake supermodule
410 if (vol8 == 1 || vol8 == 2)
415 else if (vol8 == 3 || vol8 == 4)
421 AliDebug(2,Form("Zposition = %f Edeposition = %f",zPos,edep));
434 Int_t smn = smnumber;
435 Int_t ixx = xpad - 1;
436 Int_t iyy = ypad - 1;
439 fPRE[smn][ixx][iyy] += edep;
440 fPRECounter[smn][ixx][iyy]++;
442 AliPMDcell* cell = new AliPMDcell(mtrackno,smn,ixx,iyy,edep);
447 fCPV[smn][ixx][iyy] += edep;
448 fCPVCounter[smn][ixx][iyy]++;
449 AliPMDcell* cpvcell = new AliPMDcell(mtrackno,smn,ixx,iyy,edep);
450 fCPVCell.Add(cpvcell);
454 } // Track Loop ended
455 TrackAssignment2CPVCell();
456 TrackAssignment2Cell();
464 for (Int_t idet = 0; idet < 2; idet++)
466 for (Int_t ism = 0; ism < fgkTotUM; ism++)
468 for (Int_t jrow = 0; jrow < fgkRow; jrow++)
470 for (Int_t kcol = 0; kcol < fgkCol; kcol++)
474 deltaE = fPRE[ism][jrow][kcol];
475 trno = fPRETrackNo[ism][jrow][kcol];
480 deltaE = fCPV[ism][jrow][kcol];
481 trno = fCPVTrackNo[ism][jrow][kcol];
487 TParticle *mparticle = gAlice->GetMCApp()->Particle(trno);
488 trpid = mparticle->GetPdgCode();
489 AddSDigit(trno,trpid,detno,ism,jrow,kcol,deltaE);
497 fPMDLoader->WriteSDigits("OVERWRITE");
500 //____________________________________________________________________________
502 void AliPMDDigitizer::Hits2Digits(Int_t ievt)
504 // This reads the PMD Hits tree and assigns the right track number
505 // to a cell and stores in the digits tree
507 const Int_t kPi0 = 111;
508 const Int_t kGamma = 22;
516 Float_t xPos, yPos, zPos;
517 Int_t xpad = -1, ypad = -1;
519 Float_t vx = -999.0, vy = -999.0, vz = -999.0;
521 if (!fDigits) fDigits = new TClonesArray("AliPMDdigit", 1000);
524 AliDebug(1,Form("Event Number = %d",ievt));
525 Int_t nparticles = fRunLoader->GetHeader()->GetNtrack();
526 AliDebug(1,Form("Number of Particles = %d", nparticles));
528 fRunLoader->GetEvent(ievt);
529 // ------------------------------------------------------- //
530 // Pointer to specific detector hits.
531 // Get pointers to Alice detectors and Hits containers
533 fPMD = (AliPMD*)gAlice->GetDetector("PMD");
534 fPMDLoader = fRunLoader->GetLoader("PMDLoader");
536 if (fPMDLoader == 0x0)
538 AliError("Can not find PMD or PMDLoader");
540 TTree* treeH = fPMDLoader->TreeH();
541 Int_t ntracks = (Int_t) treeH->GetEntries();
542 AliDebug(1,Form("Number of Tracks in the TreeH = %d", ntracks));
543 fPMDLoader->LoadDigits("recreate");
544 TTree* treeD = fPMDLoader->TreeD();
547 fPMDLoader->MakeTree("D");
548 treeD = fPMDLoader->TreeD();
550 Int_t bufsize = 16000;
551 treeD->Branch("PMDDigit", &fDigits, bufsize);
553 TClonesArray* hits = 0;
554 if (fPMD) hits = fPMD->Hits();
556 // Start loop on tracks in the hits containers
558 for (Int_t track=0; track<ntracks;track++)
560 gAlice->GetMCApp()->ResetHits();
561 treeH->GetEvent(track);
565 npmd = hits->GetEntriesFast();
566 for (Int_t ipmd = 0; ipmd < npmd; ipmd++)
568 fPMDHit = (AliPMDhit*) hits->UncheckedAt(ipmd);
569 trackno = fPMDHit->GetTrack();
571 // get kinematics of the particles
573 TParticle* mparticle = gAlice->GetMCApp()->Particle(trackno);
574 trackpid = mparticle->GetPdgCode();
575 Int_t ks = mparticle->GetStatusCode();
577 Int_t tracknoOld=0, trackpidOld=0, statusOld = 0;
578 if (mparticle->GetFirstMother() == -1)
580 tracknoOld = trackno;
581 trackpidOld = trackpid;
587 Int_t trnotemp = trackno; // modified on 25th Nov 2009
588 if(ks==1||(imo = mparticle->GetFirstMother())<0 ){
589 vx = mparticle->Vx();
590 vy = mparticle->Vy();
591 vz = mparticle->Vz();
593 if(trackpid==kGamma||trackpid==11||trackpid==-11||trackpid==kPi0)
598 while(((imo = mparticle->GetFirstMother()) >= 0)&&
599 (ks = mparticle->GetStatusCode() <1) )
601 mparticle = gAlice->GetMCApp()->Particle(imo);
602 trackpid = mparticle->GetPdgCode();
603 ks = mparticle->GetStatusCode();
604 vx = mparticle->Vx();
605 vy = mparticle->Vy();
606 vz = mparticle->Vz();
608 // Modified on 25th Nov 2009
621 if(trackpid==kGamma||trackpid==11||trackpid==-11||trackpid==kPi0)
625 trackpid=trackpidOld;
628 Float_t ptime = fPMDHit->GetTime()*1e6;
629 if (ptime < 0. || ptime > 1.2) continue;
634 edep = fPMDHit->GetEnergy();
635 Int_t vol1 = fPMDHit->GetVolume(1); // Column
636 Int_t vol2 = fPMDHit->GetVolume(2); // Row
637 Int_t vol7 = fPMDHit->GetVolume(4); // Serial Module No
639 // -----------------------------------------//
640 // In new geometry after adding electronics //
641 // For Super Module 1 & 2 //
642 // nrow = 48, ncol = 96 //
643 // For Super Module 3 & 4 //
644 // nrow = 96, ncol = 48 //
645 // -----------------------------------------//
653 smnumber = vol7 - 24;
655 Int_t vol8 = smnumber/6 + 1; // fake supermodule
657 if (vol8 == 1 || vol8 == 2)
662 else if (vol8 == 3 || vol8 == 4)
668 AliDebug(2,Form("ZPosition = %f Edeposition = %f",zPos,edep));
685 else if (zPos > fZPos)
691 Int_t smn = smnumber;
692 Int_t ixx = xpad - 1;
693 Int_t iyy = ypad - 1;
696 fPRE[smn][ixx][iyy] += edep;
697 fPRECounter[smn][ixx][iyy]++;
699 AliPMDcell* cell = new AliPMDcell(mtrackno,smn,ixx,iyy,edep);
705 fCPV[smn][ixx][iyy] += edep;
706 fCPVCounter[smn][ixx][iyy]++;
707 AliPMDcell* cpvcell = new AliPMDcell(mtrackno,smn,ixx,iyy,edep);
708 fCPVCell.Add(cpvcell);
712 } // Track Loop ended
713 TrackAssignment2CPVCell();
714 TrackAssignment2Cell();
724 for (Int_t idet = 0; idet < 2; idet++)
726 for (Int_t ism = 0; ism < fgkTotUM; ism++)
728 for (Int_t jrow = 0; jrow < fgkRow; jrow++)
730 for (Int_t kcol = 0; kcol < fgkCol; kcol++)
734 deltaE = fPRE[ism][jrow][kcol];
735 trno = fPRETrackNo[ism][jrow][kcol];
740 deltaE = fCPV[ism][jrow][kcol];
741 trno = fCPVTrackNo[ism][jrow][kcol];
748 // To decalibrate the adc values
750 gain1 = Gain(idet,ism,jrow,kcol);
753 Int_t adcDecalib = (Int_t)(adc/gain1);
754 adc = (Float_t) adcDecalib;
761 // Pedestal Decalibration
763 fCalibPed->GetPedMeanRms(idet,ism,jrow,kcol);
764 Int_t pedrms1 = (Int_t) pedmeanrms%100;
765 Float_t pedrms = (Float_t)pedrms1/10.;
767 (Float_t) (pedmeanrms - pedrms1)/1000.0;
770 adc += (pedmean + 3.0*pedrms);
772 = gAlice->GetMCApp()->Particle(trno);
773 trpid = mparticle->GetPdgCode();
775 AddDigit(trno,trpid,detno,ism,jrow,kcol,adc);
782 } // supermodule loop
785 fPMDLoader->WriteDigits("OVERWRITE");
789 //____________________________________________________________________________
792 void AliPMDDigitizer::SDigits2Digits(Int_t ievt)
794 // This reads the PMD sdigits tree and converts energy deposition
795 // in a cell to ADC and stores in the digits tree
798 fRunLoader->GetEvent(ievt);
800 TTree* treeS = fPMDLoader->TreeS();
801 AliPMDsdigit *pmdsdigit;
802 TBranch *branch = treeS->GetBranch("PMDSDigit");
805 AliError("PMD Sdigit branch does not exist");
808 if (!fSDigits) fSDigits = new TClonesArray("AliPMDsdigit", 1000);
809 branch->SetAddress(&fSDigits);
811 TTree* treeD = fPMDLoader->TreeD();
814 fPMDLoader->MakeTree("D");
815 treeD = fPMDLoader->TreeD();
817 Int_t bufsize = 16000;
818 if (!fDigits) fDigits = new TClonesArray("AliPMDdigit", 1000);
819 treeD->Branch("PMDDigit", &fDigits, bufsize);
821 Int_t trno, trpid, det, smn;
825 Int_t nmodules = (Int_t) treeS->GetEntries();
826 AliDebug(1,Form("Number of modules = %d",nmodules));
828 for (Int_t imodule = 0; imodule < nmodules; imodule++)
830 treeS->GetEntry(imodule);
831 Int_t nentries = fSDigits->GetLast();
832 AliDebug(2,Form("Number of entries per module = %d",nentries+1));
833 for (Int_t ient = 0; ient < nentries+1; ient++)
835 pmdsdigit = (AliPMDsdigit*)fSDigits->UncheckedAt(ient);
836 trno = pmdsdigit->GetTrackNumber();
837 trpid = pmdsdigit->GetTrackPid();
838 det = pmdsdigit->GetDetector();
839 smn = pmdsdigit->GetSMNumber();
840 irow = pmdsdigit->GetRow();
841 icol = pmdsdigit->GetColumn();
842 edep = pmdsdigit->GetCellEdep();
846 // To decalibrte the adc values
848 Float_t gain1 = Gain(det,smn,irow,icol);
851 Int_t adcDecalib = (Int_t)(adc/gain1);
852 adc = (Float_t) adcDecalib;
858 // Pedestal Decalibration
859 Int_t pedmeanrms = fCalibPed->GetPedMeanRms(det,smn,irow,icol);
860 Int_t pedrms1 = (Int_t) pedmeanrms%100;
861 Float_t pedrms = (Float_t)pedrms1/10.;
862 Float_t pedmean = (Float_t) (pedmeanrms - pedrms1)/1000.0;
865 adc += (pedmean + 3.0*pedrms);
866 AddDigit(trno,trpid,det,smn,irow,icol,adc);
873 fPMDLoader->WriteDigits("OVERWRITE");
876 //____________________________________________________________________________
877 void AliPMDDigitizer::Exec(Option_t *option)
879 // Does the event merging and digitization
880 const char *cdeb = strstr(option,"deb");
883 AliDebug(100," *** PMD Exec is called ***");
886 Int_t ninputs = fManager->GetNinputs();
887 AliDebug(1,Form("Number of files to be processed = %d",ninputs));
890 for (Int_t i = 0; i < ninputs; i++)
892 Int_t troffset = fManager->GetMask(i);
893 MergeSDigits(i, troffset);
896 fRunLoader = AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
897 fPMD = (AliPMD*)gAlice->GetDetector("PMD");
898 fPMDLoader = fRunLoader->GetLoader("PMDLoader");
899 if (fPMDLoader == 0x0)
901 AliError("Can not find PMD or PMDLoader");
903 fPMDLoader->LoadDigits("update");
904 TTree* treeD = fPMDLoader->TreeD();
907 fPMDLoader->MakeTree("D");
908 treeD = fPMDLoader->TreeD();
910 Int_t bufsize = 16000;
911 if (!fDigits) fDigits = new TClonesArray("AliPMDdigit", 1000);
912 treeD->Branch("PMDDigit", &fDigits, bufsize);
920 for (Int_t idet = 0; idet < 2; idet++)
922 for (Int_t ism = 0; ism < fgkTotUM; ism++)
924 for (Int_t jrow = 0; jrow < fgkRow; jrow++)
926 for (Int_t kcol = 0; kcol < fgkCol; kcol++)
930 deltaE = fPRE[ism][jrow][kcol];
931 trno = fPRETrackNo[ism][jrow][kcol];
932 trpid = fPRETrackPid[ism][jrow][kcol];
937 deltaE = fCPV[ism][jrow][kcol];
938 trno = fCPVTrackNo[ism][jrow][kcol];
939 trpid = fCPVTrackPid[ism][jrow][kcol];
947 // Gain decalibration
949 Float_t gain1 = Gain(idet,ism,jrow,kcol);
953 Int_t adcDecalib = (Int_t)(adc/gain1);
954 adc = (Float_t) adcDecalib;
960 // Pedestal Decalibration
962 fCalibPed->GetPedMeanRms(idet,ism,jrow,kcol);
963 Int_t pedrms1 = (Int_t) pedmeanrms%100;
964 Float_t pedrms = (Float_t)pedrms1/10.;
966 (Float_t) (pedmeanrms - pedrms1)/1000.0;
969 adc += (pedmean + 3.0*pedrms);
970 AddDigit(trno,trpid,detno,ism,jrow,kcol,adc);
978 } // supermodule loop
980 fPMDLoader->WriteDigits("OVERWRITE");
981 fPMDLoader->UnloadDigits();
984 //____________________________________________________________________________
985 void AliPMDDigitizer::TrackAssignment2CPVCell()
987 // This block assigns the cell id when there are
988 // multiple tracks in a cell according to the
990 // This method added by Ajay
991 Bool_t jsort = false;
1002 Float_t ****cpvEdep;
1004 cpvTrack = new Int_t ***[fgkTotUM];
1005 cpvEdep = new Float_t ***[fgkTotUM];
1006 for (i=0; i<fgkTotUM; i++)
1008 cpvTrack[i] = new Int_t **[fgkRow];
1009 cpvEdep[i] = new Float_t **[fgkRow];
1012 for (i = 0; i < fgkTotUM; i++)
1014 for (j = 0; j < fgkRow; j++)
1016 cpvTrack[i][j] = new Int_t *[fgkCol];
1017 cpvEdep[i][j] = new Float_t *[fgkCol];
1020 for (i = 0; i < fgkTotUM; i++)
1022 for (j = 0; j < fgkRow; j++)
1024 for (k = 0; k < fgkCol; k++)
1026 Int_t nn = fCPVCounter[i][j][k];
1029 cpvTrack[i][j][k] = new Int_t[nn];
1030 cpvEdep[i][j][k] = new Float_t[nn];
1035 cpvTrack[i][j][k] = new Int_t[nn];
1036 cpvEdep[i][j][k] = new Float_t[nn];
1038 fCPVCounter[i][j][k] = 0;
1044 Int_t nentries = fCPVCell.GetEntries();
1046 Int_t mtrackno, ism, ixp, iyp;
1048 for (i = 0; i < nentries; i++)
1050 AliPMDcell* cpvcell = (AliPMDcell*)fCPVCell.UncheckedAt(i);
1052 mtrackno = cpvcell->GetTrackNumber();
1053 ism = cpvcell->GetSMNumber();
1054 ixp = cpvcell->GetX();
1055 iyp = cpvcell->GetY();
1056 edep = cpvcell->GetEdep();
1057 Int_t nn = fCPVCounter[ism][ixp][iyp];
1058 cpvTrack[ism][ixp][iyp][nn] = (Int_t) mtrackno;
1059 cpvEdep[ism][ixp][iyp][nn] = edep;
1060 fCPVCounter[ism][ixp][iyp]++;
1066 for (im=0; im<fgkTotUM; im++)
1068 for (ix=0; ix<fgkRow; ix++)
1070 for (iy=0; iy<fgkCol; iy++)
1072 nn = fCPVCounter[im][ix][iy];
1075 // This block handles if a cell is fired
1076 // many times by many tracks
1077 status1 = new Int_t[nn];
1078 status2 = new Int_t[nn];
1079 trnarray = new Int_t[nn];
1080 for (iz = 0; iz < nn; iz++)
1082 status1[iz] = cpvTrack[im][ix][iy][iz];
1084 TMath::Sort(nn,status1,status2,jsort);
1085 Int_t trackOld = -99999;
1086 Int_t track, trCount = 0;
1087 for (iz = 0; iz < nn; iz++)
1089 track = status1[status2[iz]];
1090 if (trackOld != track)
1092 trnarray[trCount] = track;
1099 Float_t totEdp = 0.;
1100 trEdp = new Float_t[trCount];
1101 fracEdp = new Float_t[trCount];
1102 for (il = 0; il < trCount; il++)
1105 track = trnarray[il];
1106 for (iz = 0; iz < nn; iz++)
1108 if (track == cpvTrack[im][ix][iy][iz])
1110 trEdp[il] += cpvEdep[im][ix][iy][iz];
1113 totEdp += trEdp[il];
1116 Float_t fracOld = 0.;
1118 for (il = 0; il < trCount; il++)
1120 fracEdp[il] = trEdp[il]/totEdp;
1121 if (fracOld < fracEdp[il])
1123 fracOld = fracEdp[il];
1127 fCPVTrackNo[im][ix][iy] = trnarray[ilOld];
1134 // This only handles if a cell is fired
1135 // by only one track
1137 fCPVTrackNo[im][ix][iy] = cpvTrack[im][ix][iy][0];
1142 // This is if no cell is fired
1143 fCPVTrackNo[im][ix][iy] = -999;
1149 // Delete all the pointers
1151 for (i = 0; i < fgkTotUM; i++)
1153 for (j = 0; j < fgkRow; j++)
1155 for (k = 0; k < fgkCol; k++)
1157 delete []cpvTrack[i][j][k];
1158 delete []cpvEdep[i][j][k];
1163 for (i = 0; i < fgkTotUM; i++)
1165 for (j = 0; j < fgkRow; j++)
1167 delete [] cpvTrack[i][j];
1168 delete [] cpvEdep[i][j];
1172 for (i = 0; i < fgkTotUM; i++)
1174 delete [] cpvTrack[i];
1175 delete [] cpvEdep[i];
1181 // End of the cell id assignment
1184 //____________________________________________________________________________
1186 void AliPMDDigitizer::MergeSDigits(Int_t filenumber, Int_t troffset)
1189 fRunLoader = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(filenumber));
1190 fPMDLoader = fRunLoader->GetLoader("PMDLoader");
1191 fPMDLoader->LoadSDigits("read");
1192 TTree* treeS = fPMDLoader->TreeS();
1193 AliPMDsdigit *pmdsdigit;
1194 TBranch *branch = treeS->GetBranch("PMDSDigit");
1195 if (!fSDigits) fSDigits = new TClonesArray("AliPMDsdigit", 1000);
1196 branch->SetAddress(&fSDigits);
1198 Int_t itrackno, itrackpid, idet, ism;
1201 Int_t nmodules = (Int_t) treeS->GetEntries();
1202 AliDebug(1,Form("Number of Modules in the treeS = %d",nmodules));
1203 AliDebug(1,Form("Track Offset = %d",troffset));
1204 for (Int_t imodule = 0; imodule < nmodules; imodule++)
1206 treeS->GetEntry(imodule);
1207 Int_t nentries = fSDigits->GetLast();
1208 AliDebug(2,Form("Number of Entries per Module = %d",nentries));
1209 for (Int_t ient = 0; ient < nentries+1; ient++)
1211 pmdsdigit = (AliPMDsdigit*)fSDigits->UncheckedAt(ient);
1212 itrackno = pmdsdigit->GetTrackNumber();
1213 itrackpid = pmdsdigit->GetTrackPid();
1214 idet = pmdsdigit->GetDetector();
1215 ism = pmdsdigit->GetSMNumber();
1216 ixp = pmdsdigit->GetRow();
1217 iyp = pmdsdigit->GetColumn();
1218 edep = pmdsdigit->GetCellEdep();
1221 if (fPRE[ism][ixp][iyp] < edep)
1223 fPRETrackNo[ism][ixp][iyp] = troffset + itrackno;
1224 fPRETrackPid[ism][ixp][iyp] = itrackpid;
1226 fPRE[ism][ixp][iyp] += edep;
1230 if (fCPV[ism][ixp][iyp] < edep)
1232 fCPVTrackNo[ism][ixp][iyp] = troffset + itrackno;
1233 fCPVTrackPid[ism][ixp][iyp] = itrackpid;
1235 fCPV[ism][ixp][iyp] += edep;
1241 // ----------------------------------------------------------------------
1242 void AliPMDDigitizer::TrackAssignment2Cell()
1245 // This block assigns the cell id when there are
1246 // multiple tracks in a cell according to the
1247 // energy deposition
1249 Bool_t jsort = false;
1260 Float_t ****pmdEdep;
1262 pmdTrack = new Int_t ***[fgkTotUM];
1263 pmdEdep = new Float_t ***[fgkTotUM];
1264 for (i=0; i<fgkTotUM; i++)
1266 pmdTrack[i] = new Int_t **[fgkRow];
1267 pmdEdep[i] = new Float_t **[fgkRow];
1270 for (i = 0; i < fgkTotUM; i++)
1272 for (j = 0; j < fgkRow; j++)
1274 pmdTrack[i][j] = new Int_t *[fgkCol];
1275 pmdEdep[i][j] = new Float_t *[fgkCol];
1279 for (i = 0; i < fgkTotUM; i++)
1281 for (j = 0; j < fgkRow; j++)
1283 for (k = 0; k < fgkCol; k++)
1285 Int_t nn = fPRECounter[i][j][k];
1288 pmdTrack[i][j][k] = new Int_t[nn];
1289 pmdEdep[i][j][k] = new Float_t[nn];
1294 pmdTrack[i][j][k] = new Int_t[nn];
1295 pmdEdep[i][j][k] = new Float_t[nn];
1297 fPRECounter[i][j][k] = 0;
1303 Int_t nentries = fCell.GetEntries();
1305 Int_t mtrackno, ism, ixp, iyp;
1308 for (i = 0; i < nentries; i++)
1310 AliPMDcell* cell = (AliPMDcell*)fCell.UncheckedAt(i);
1312 mtrackno = cell->GetTrackNumber();
1313 ism = cell->GetSMNumber();
1316 edep = cell->GetEdep();
1317 Int_t nn = fPRECounter[ism][ixp][iyp];
1318 pmdTrack[ism][ixp][iyp][nn] = (Int_t) mtrackno;
1319 pmdEdep[ism][ixp][iyp][nn] = edep;
1320 fPRECounter[ism][ixp][iyp]++;
1327 for (im=0; im<fgkTotUM; im++)
1329 for (ix=0; ix<fgkRow; ix++)
1331 for (iy=0; iy<fgkCol; iy++)
1333 nn = fPRECounter[im][ix][iy];
1336 // This block handles if a cell is fired
1337 // many times by many tracks
1338 status1 = new Int_t[nn];
1339 status2 = new Int_t[nn];
1340 trnarray = new Int_t[nn];
1341 for (iz = 0; iz < nn; iz++)
1343 status1[iz] = pmdTrack[im][ix][iy][iz];
1345 TMath::Sort(nn,status1,status2,jsort);
1346 Int_t trackOld = -99999;
1347 Int_t track, trCount = 0;
1348 for (iz = 0; iz < nn; iz++)
1350 track = status1[status2[iz]];
1351 if (trackOld != track)
1353 trnarray[trCount] = track;
1360 Float_t totEdp = 0.;
1361 trEdp = new Float_t[trCount];
1362 fracEdp = new Float_t[trCount];
1363 for (il = 0; il < trCount; il++)
1366 track = trnarray[il];
1367 for (iz = 0; iz < nn; iz++)
1369 if (track == pmdTrack[im][ix][iy][iz])
1371 trEdp[il] += pmdEdep[im][ix][iy][iz];
1374 totEdp += trEdp[il];
1377 Float_t fracOld = 0.;
1379 for (il = 0; il < trCount; il++)
1381 fracEdp[il] = trEdp[il]/totEdp;
1382 if (fracOld < fracEdp[il])
1384 fracOld = fracEdp[il];
1388 fPRETrackNo[im][ix][iy] = trnarray[ilOld];
1395 // This only handles if a cell is fired
1396 // by only one track
1398 fPRETrackNo[im][ix][iy] = pmdTrack[im][ix][iy][0];
1403 // This is if no cell is fired
1404 fPRETrackNo[im][ix][iy] = -999;
1410 // Delete all the pointers
1412 for (i = 0; i < fgkTotUM; i++)
1414 for (j = 0; j < fgkRow; j++)
1416 for (k = 0; k < fgkCol; k++)
1418 delete [] pmdTrack[i][j][k];
1419 delete [] pmdEdep[i][j][k];
1424 for (i = 0; i < fgkTotUM; i++)
1426 for (j = 0; j < fgkRow; j++)
1428 delete [] pmdTrack[i][j];
1429 delete [] pmdEdep[i][j];
1433 for (i = 0; i < fgkTotUM; i++)
1435 delete [] pmdTrack[i];
1436 delete [] pmdEdep[i];
1441 // End of the cell id assignment
1444 //____________________________________________________________________________
1445 void AliPMDDigitizer::MeV2ADC(Float_t mev, Float_t & adc) const
1447 // This converts the simulated edep to ADC according to the
1449 // PS Test in June 2010, Voltage @ 1300 V
1450 // KeV - ADC conversion for 12bit ADC
1451 // MPV data used for the fit and taken here
1453 const Float_t kConstant = 16.807;
1454 const Float_t kErConstant = 6.62215;
1455 const Float_t kSlope = 128.875;
1456 const Float_t kErSlope = 6.21713;
1458 Float_t cons = gRandom->Gaus(kConstant,kErConstant);
1459 Float_t slop = gRandom->Gaus(kSlope,kErSlope);
1461 Float_t adc12bit = slop*mev*0.001 + cons;
1463 if (adc12bit < 0.) adc12bit = 0.;
1465 if(adc12bit < 1600.0)
1467 adc = (Float_t) adc12bit;
1469 else if (adc12bit >= 1600.0)
1474 //____________________________________________________________________________
1475 void AliPMDDigitizer::AddSDigit(Int_t trnumber, Int_t trpid, Int_t det,
1476 Int_t smnumber, Int_t irow, Int_t icol,
1481 if (!fSDigits) fSDigits = new TClonesArray("AliPMDsdigit", 1000);
1482 TClonesArray &lsdigits = *fSDigits;
1483 new(lsdigits[fNsdigit++]) AliPMDsdigit(trnumber,trpid,det,smnumber,irow,icol,adc);
1485 //____________________________________________________________________________
1487 void AliPMDDigitizer::AddDigit(Int_t trnumber, Int_t trpid, Int_t det,
1488 Int_t smnumber, Int_t irow, Int_t icol,
1493 if (!fDigits) fDigits = new TClonesArray("AliPMDdigit", 1000);
1494 TClonesArray &ldigits = *fDigits;
1495 new(ldigits[fNdigit++]) AliPMDdigit(trnumber,trpid, det,smnumber,irow,icol,adc);
1497 //____________________________________________________________________________
1499 void AliPMDDigitizer::SetZPosition(Float_t zpos)
1503 //____________________________________________________________________________
1504 Float_t AliPMDDigitizer::GetZPosition() const
1508 //____________________________________________________________________________
1510 void AliPMDDigitizer::ResetCell()
1512 // clears the cell array and also the counter
1517 for (Int_t i = 0; i < fgkTotUM; i++)
1519 for (Int_t j = 0; j < fgkRow; j++)
1521 for (Int_t k = 0; k < fgkCol; k++)
1523 fCPVCounter[i][j][k] = 0;
1524 fPRECounter[i][j][k] = 0;
1529 //____________________________________________________________________________
1530 void AliPMDDigitizer::ResetSDigit()
1534 if (fSDigits) fSDigits->Delete();
1536 //____________________________________________________________________________
1537 void AliPMDDigitizer::ResetDigit()
1541 if (fDigits) fDigits->Delete();
1543 //____________________________________________________________________________
1545 void AliPMDDigitizer::ResetCellADC()
1547 // Clears individual cells edep and track number
1548 for (Int_t i = 0; i < fgkTotUM; i++)
1550 for (Int_t j = 0; j < fgkRow; j++)
1552 for (Int_t k = 0; k < fgkCol; k++)
1556 fCPVTrackNo[i][j][k] = 0;
1557 fPRETrackNo[i][j][k] = 0;
1558 fCPVTrackPid[i][j][k] = -1;
1559 fPRETrackPid[i][j][k] = -1;
1564 //____________________________________________________________________________
1566 void AliPMDDigitizer::UnLoad(Option_t *option)
1568 // Unloads all the root files
1570 const char *cS = strstr(option,"S");
1571 const char *cD = strstr(option,"D");
1573 fRunLoader->UnloadgAlice();
1574 fRunLoader->UnloadHeader();
1575 fRunLoader->UnloadKinematics();
1579 fPMDLoader->UnloadHits();
1583 fPMDLoader->UnloadHits();
1584 fPMDLoader->UnloadSDigits();
1588 //----------------------------------------------------------------------
1589 Float_t AliPMDDigitizer::Gain(Int_t det, Int_t smn, Int_t row, Int_t col) const
1591 // returns of the gain of the cell
1592 // Added this method by ZA
1594 //cout<<" I am here in gain "<<fCalibData<< "smn,row, col "<<smn
1595 //<<" "<<row<<" "<<col<<endl;
1598 AliError("No calibration data loaded from CDB!!!");
1603 GainFact = fCalibGain->GetGainFact(det,smn,row,col);
1606 //----------------------------------------------------------------------
1607 AliPMDCalibData* AliPMDDigitizer::GetCalibGain() const
1609 // The run number will be centralized in AliCDBManager,
1610 // you don't need to set it here!
1611 // Added this method by ZA
1612 // Cleaned up by Alberto
1613 AliCDBEntry *entry = AliCDBManager::Instance()->Get("PMD/Calib/Gain");
1615 if(!entry) AliFatal("Calibration object retrieval failed!");
1617 AliPMDCalibData *calibdata=0;
1618 if (entry) calibdata = (AliPMDCalibData*) entry->GetObject();
1620 if (!calibdata) AliFatal("No calibration data from calibration database !");
1624 //----------------------------------------------------------------------
1625 AliPMDPedestal* AliPMDDigitizer::GetCalibPed() const
1627 // The run number will be centralized in AliCDBManager,
1628 // you don't need to set it here!
1630 AliCDBEntry *entry = AliCDBManager::Instance()->Get("PMD/Calib/Ped");
1632 if(!entry) AliFatal("Pedestal object retrieval failed!");
1634 AliPMDPedestal *pedestal=0;
1635 if (entry) pedestal = (AliPMDPedestal*) entry->GetObject();
1637 if (!pedestal) AliFatal("No pedestal data from calibration database !");