1 /**************************************************************************
2 * Copyright(c) 1998-2000, 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 **************************************************************************/
21 #include <AliRunLoader.h>
22 #include "AliRunDigitizer.h"
23 #include <AliLoader.h>
25 #include "AliHMPIDDigitizer.h"
26 #include "AliHMPIDDigit.h"
28 #include "AliHMPIDParam.h"
32 ClassImp(AliHMPIDDigitizer)
34 Bool_t AliHMPIDDigitizer::fgDoNoise=kFALSE;
35 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
36 void AliHMPIDDigitizer::Exec(Option_t*)
38 // This methode is responsible for merging sdigits to a list of digits
39 //Disintegration leeds to the fact that one hit affects several neighbouring pads, which means that the same pad might be affected by few hits.
40 AliDebug(1,Form("Start with %i input(s) for event %i",fManager->GetNinputs(),fManager->GetOutputEventNr()));
41 //First we read all sdigits from all inputs
42 AliRunLoader *pInRunLoader=0;//in and out Run loaders
43 AliLoader *pInRichLoader=0;//in and out HMPID loaders
44 static TClonesArray sdigs("AliHMPIDDigit");//tmp storage for sdigits sum up from all input files
46 for(Int_t inFileN=0;inFileN<fManager->GetNinputs();inFileN++){//files loop
47 pInRunLoader = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(inFileN)); //get run loader from current input
48 pInRichLoader = pInRunLoader->GetLoader("HMPIDLoader"); if(pInRichLoader==0) continue; //no HMPID in this input, check the next input
49 if (!pInRunLoader->GetAliRun()) pInRunLoader->LoadgAlice();
50 AliHMPID* pInRich=(AliHMPID*)pInRunLoader->GetAliRun()->GetDetector("HMPID"); //take HMPID from current input
51 pInRichLoader->LoadSDigits(); pInRichLoader->TreeS()->GetEntry(0); //take list of HMPID sdigits from current input
52 AliDebug(1,Form("input %i has %i sdigits",inFileN,pInRich->SdiLst()->GetEntries()));
53 for(Int_t i=0;i<pInRich->SdiLst()->GetEntries();i++){ //collect sdigits from current input
54 AliHMPIDDigit *pSDig=(AliHMPIDDigit*)pInRich->SdiLst()->At(i);
55 pSDig->AddTidOffset(fManager->GetMask(inFileN)); //apply TID shift since all inputs count tracks independently starting from 0
56 new(sdigs[total++]) AliHMPIDDigit(*pSDig);
58 pInRichLoader->UnloadSDigits(); pInRich->SdiReset(); //close current input and reset
61 //PH if(sdigs.GetEntries()==0) return; //no sdigits collected, nothing to convert
63 AliRunLoader *pOutRunLoader = AliRunLoader::GetRunLoader(fManager->GetOutputFolderName()); //open output stream (only 1 possible)
64 AliLoader *pOutRichLoader = pOutRunLoader->GetLoader("HMPIDLoader"); //take output HMPID loader
65 AliRun *pArun = pOutRunLoader->GetAliRun();
66 AliHMPID *pOutRich = (AliHMPID*)pArun->GetDetector("HMPID"); //take output HMPID
67 pOutRichLoader->MakeTree("D"); pOutRich->MakeBranch("D"); //create TreeD in output stream
69 Sdi2Dig(&sdigs,pOutRich->DigLst());
71 pOutRichLoader->TreeD()->Fill(); //fill the output tree with the list of digits
72 pOutRichLoader->WriteDigits("OVERWRITE"); //serialize them to file
74 sdigs.Clear(); //remove all tmp sdigits
75 pOutRichLoader->UnloadDigits(); pOutRich->DigReset();
77 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
78 void AliHMPIDDigitizer::Sdi2Dig(TClonesArray *pSdiLst,TObjArray *pDigLst)
80 // Converts list of sdigits to 7 lists of digits, one per each chamber
81 // Arguments: pSDigLst - list of all sdigits
82 // pDigLst - list of 7 lists of digits
85 TClonesArray *pLst[7]; Int_t iCnt[7];
87 for(Int_t i=0;i<7;i++){
88 pLst[i]=(TClonesArray*)(*pDigLst)[i];
89 iCnt[i]=0; if(pLst[i]->GetEntries()!=0) AliErrorClass("Some of digits lists is not empty"); //in principle those lists should be empty
93 Float_t arrNoise[7][6][80][48];
95 for (Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++)
96 for (Int_t iPc=AliHMPIDParam::kMinPc;iPc<=AliHMPIDParam::kMaxPc;iPc++)
97 for(Int_t iPx=AliHMPIDParam::kMinPx;iPx<=AliHMPIDParam::kMaxPx;iPx++)
98 for(Int_t iPy=AliHMPIDParam::kMinPy;iPy<=AliHMPIDParam::kMaxPy;iPy++) arrNoise[iCh][iPc][iPx][iPy] = gRandom->Gaus(0,1.1);
103 Int_t iPad=-1,iCh=-1,iNdigPad=-1,aTids[3]={-1,-1,-1}; Float_t q=-1;
104 for(Int_t i=0;i<pSdiLst->GetEntries();i++){ //sdigits loop (sorted)
105 AliHMPIDDigit *pSdig=(AliHMPIDDigit*)pSdiLst->At(i); //take current sdigit
106 if(pSdig->Pad()==iPad){ //if the same pad
107 q+=pSdig->Q(); //sum up charge
108 iNdigPad++; if(iNdigPad<=3) aTids[iNdigPad-1]=pSdig->GetTrack(0); //collect TID
111 if(i!=0 && AliHMPIDParam::IsOverTh(q)
112 && iCh>=AliHMPIDParam::kMinCh
113 && iCh<=AliHMPIDParam::kMaxCh) new((*pLst[iCh])[iCnt[iCh]++]) AliHMPIDDigit(iPad,(Int_t)q,aTids); //do not create digit for the very first sdigit
115 iPad=pSdig->Pad(); iCh=AliHMPIDParam::A2C(iPad); //new sdigit comes, reset collectors
117 aTids[0]=pSdig->GetTrack(0);aTids[1]=aTids[2]=-1;
119 if(fgDoNoise) q+=arrNoise[iCh][pSdig->Pc()][pSdig->PadPcX()][pSdig->PadPcY()];
120 arrNoise[iCh][pSdig->Pc()][pSdig->PadPcX()][pSdig->PadPcY()]=0;
121 }//sdigits loop (sorted)
123 if(AliHMPIDParam::IsOverTh(q)
124 && iCh>=AliHMPIDParam::kMinCh
125 && iCh<=AliHMPIDParam::kMaxCh) new((*pLst[iCh])[iCnt[iCh]++]) AliHMPIDDigit(iPad,(Int_t)q,aTids); //add the last one, in case of empty sdigits list q=-1
127 // add noise pad above threshold with no signal merged...if any
128 if(!fgDoNoise) return;
129 aTids[0]=aTids[1]=aTids[2]=-1;
130 for (Int_t iChCurr=AliHMPIDParam::kMinCh;iChCurr<=AliHMPIDParam::kMaxCh;iChCurr++)
131 for (Int_t iPc=AliHMPIDParam::kMinPc;iPc<=AliHMPIDParam::kMaxPc;iPc++)
132 for(Int_t iPx=AliHMPIDParam::kMinPx;iPx<=AliHMPIDParam::kMaxPx;iPx++)
133 for(Int_t iPy=AliHMPIDParam::kMinPy;iPy<=AliHMPIDParam::kMaxPy;iPy++) {
134 Float_t qNoise = arrNoise[iChCurr][iPc][iPx][iPy];
135 if(AliHMPIDParam::IsOverTh(qNoise)) new((*pLst[iChCurr])[iCnt[iChCurr]++]) AliHMPIDDigit(AliHMPIDParam::Abs(iChCurr,iPc,iPx,iPy),(Int_t)qNoise,aTids);
139 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++