]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDDigitizer.cxx
Coding rule violation corrected.
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDDigitizer.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2000, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15 //.
16 //.
17 //.
18 //.
19 //.
20 #include <AliRun.h>
21 #include <AliRunLoader.h>
22 #include "AliDigitizationInput.h"
23 #include <AliLoader.h>
24 #include "AliLog.h"
25 #include <AliCDBEntry.h> 
26 #include <AliCDBManager.h>
27 #include "AliHMPIDDigitizer.h"
28 #include "AliHMPIDReconstructor.h"
29 #include "AliHMPIDDigit.h"
30 #include "AliHMPID.h"
31 #include "AliHMPIDParam.h"
32 #include <TRandom.h>
33 #include <TMath.h>
34 #include <TTree.h>
35 #include <TObjArray.h>
36
37 ClassImp(AliHMPIDDigitizer)
38
39 Bool_t AliHMPIDDigitizer::fgDoNoise=kTRUE;
40 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
41 void AliHMPIDDigitizer::Digitize(Option_t*)
42 {
43 // This methode is responsible for merging sdigits to a list of digits
44 //Disintegration leeds to the fact that one hit affects several neighbouring pads, which means that the same pad might be affected by few hits.     
45   AliDebug(1,Form("Start with %i input(s) for event %i",fDigInput->GetNinputs(),fDigInput->GetOutputEventNr()));
46 //First we read all sdigits from all inputs  
47   AliRunLoader *pInRunLoader=0;//in and out Run loaders
48   AliLoader    *pInRichLoader=0;//in and out HMPID loaders  
49   static TClonesArray sdigs("AliHMPIDDigit");//tmp storage for sdigits sum up from all input files
50   Int_t total=0;
51   for(Int_t inFileN=0;inFileN<fDigInput->GetNinputs();inFileN++){//files loop
52     pInRunLoader  = AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(inFileN));          //get run loader from current input 
53     pInRichLoader = pInRunLoader->GetLoader("HMPIDLoader"); if(pInRichLoader==0) continue;       //no HMPID in this input, check the next input
54     if (!pInRunLoader->GetAliRun()) pInRunLoader->LoadgAlice();
55     AliHMPID* pInRich=(AliHMPID*)pInRunLoader->GetAliRun()->GetDetector("HMPID");                  //take HMPID from current input
56     pInRichLoader->LoadSDigits(); pInRichLoader->TreeS()->GetEntry(0);                          //take list of HMPID sdigits from current input 
57     AliDebug(1,Form("input %i has %i sdigits",inFileN,pInRich->SdiLst()->GetEntries()));
58     for(Int_t i=0;i<pInRich->SdiLst()->GetEntries();i++){                                        //collect sdigits from current input
59       AliHMPIDDigit *pSDig=(AliHMPIDDigit*)pInRich->SdiLst()->At(i);
60       pSDig->AddTidOffset(fDigInput->GetMask(inFileN));                                          //apply TID shift since all inputs count tracks independently starting from 0
61       new(sdigs[total++]) AliHMPIDDigit(*pSDig);       
62     }
63     pInRichLoader->UnloadSDigits();   pInRich->SdiReset(); //close current input and reset 
64   }//files loop
65
66   //PH  if(sdigs.GetEntries()==0) return;                                                              //no sdigits collected, nothing to convert  
67   
68   AliRunLoader *pOutRunLoader  = AliRunLoader::GetRunLoader(fDigInput->GetOutputFolderName());    //open output stream (only 1 possible)
69   AliLoader    *pOutRichLoader = pOutRunLoader->GetLoader("HMPIDLoader");                         //take output HMPID loader
70   AliRun *pArun = pOutRunLoader->GetAliRun();
71   AliHMPID      *pOutRich       = (AliHMPID*)pArun->GetDetector("HMPID");      //take output HMPID
72   pOutRichLoader->MakeTree("D");   pOutRich->MakeBranch("D");                                    //create TreeD in output stream
73
74   Sdi2Dig(&sdigs,pOutRich->DigLst());
75   
76   pOutRichLoader->TreeD()->Fill();              //fill the output tree with the list of digits
77   pOutRichLoader->WriteDigits("OVERWRITE");     //serialize them to file
78   
79   sdigs.Clear();                      //remove all tmp sdigits
80   pOutRichLoader->UnloadDigits();   pOutRich->DigReset();
81 }//Exec()
82 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
83 void AliHMPIDDigitizer::Sdi2Dig(TClonesArray *pSdiLst,TObjArray *pDigLst)
84 {
85 // Converts list of sdigits to 7 lists of digits, one per each chamber
86 // Arguments: pSDigLst - list of all sdigits
87 //            pDigLst  - list of 7 lists of digits        
88 //   Returns: none  
89   
90   TClonesArray *pLst[7]; Int_t iCnt[7];
91
92   for(Int_t i=0;i<7;i++){
93     pLst[i]=(TClonesArray*)(*pDigLst)[i];
94     iCnt[i]=0; if(pLst[i]->GetEntries()!=0) AliErrorClass("Some of digits lists is not empty");         //in principle those lists should be empty                                                                       
95   }
96   
97   // make noise array
98   Float_t arrNoise[7][6][80][48], arrSigmaPed[7][6][80][48];
99   if(fgDoNoise) {
100     for (Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++)
101       for (Int_t iPc=AliHMPIDParam::kMinPc;iPc<=AliHMPIDParam::kMaxPc;iPc++)
102         for(Int_t iPx=AliHMPIDParam::kMinPx;iPx<=AliHMPIDParam::kMaxPx;iPx++)
103           for(Int_t iPy=AliHMPIDParam::kMinPy;iPy<=AliHMPIDParam::kMaxPy;iPy++){
104             arrNoise[iCh][iPc][iPx][iPy] = gRandom->Gaus(0,1.);
105             arrSigmaPed[iCh][iPc][iPx][iPy] = 1.;
106           }
107           
108     AliCDBEntry *pDaqSigEnt = AliCDBManager::Instance()->Get("HMPID/Calib/DaqSig");  //contains TObjArray of TObjArray 14 TMatrixF sigmas values for pads 
109    
110     if(pDaqSigEnt){
111       TObjArray *pDaqSig = (TObjArray*)pDaqSigEnt->GetObject();
112       for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){                  //chambers loop    
113         TMatrixF *pM = (TMatrixF*)pDaqSig->At(iCh);     
114         for (Int_t iPc=AliHMPIDParam::kMinPc;iPc<=AliHMPIDParam::kMaxPc;iPc++)
115           for(Int_t iPx=AliHMPIDParam::kMinPx;iPx<=AliHMPIDParam::kMaxPx;iPx++)
116             for(Int_t iPy=AliHMPIDParam::kMinPy;iPy<=AliHMPIDParam::kMaxPy;iPy++){
117               Int_t padX = (iPc%2)*AliHMPIDParam::kPadPcX+iPx; 
118               Int_t padY = (iPc/2)*AliHMPIDParam::kPadPcY+iPy;
119               if((*pM)(padX,padY)>0.){
120                 arrNoise[iCh][iPc][iPx][iPy] = gRandom->Gaus(0,(*pM)(padX,padY));
121                 arrSigmaPed[iCh][iPc][iPx][iPy] = (*pM)(padX,padY);}
122               else{
123                 arrNoise[iCh][iPc][iPx][iPy] = gRandom->Gaus(0,1.);
124                 arrSigmaPed[iCh][iPc][iPx][iPy] = 1.;}
125             } 
126       }
127     }
128   }  
129   
130   pSdiLst->Sort();  
131                      
132   Int_t iPad=-1,iCh=-1,iNdigPad=-1,aTids[3]={-1,-1,-1}; Float_t q=-1;
133   for(Int_t i=0;i<pSdiLst->GetEntries();i++){                                                                  //sdigits loop (sorted)
134     AliHMPIDDigit *pSdig=(AliHMPIDDigit*)pSdiLst->At(i);                                                       //take current sdigit
135     if(pSdig->Pad()==iPad){                                                                                    //if the same pad 
136       q+=pSdig->Q();                                                                                           //sum up charge
137       iNdigPad++; if(iNdigPad<=3) aTids[iNdigPad-1]=pSdig->GetTrack(0);                                        //collect TID 
138       continue;
139     }
140     if(i!=0 && iCh>=AliHMPIDParam::kMinCh && iCh<=AliHMPIDParam::kMaxCh){
141       AliHMPIDParam::Instance()->SetThreshold((TMath::Nint(arrSigmaPed[iCh][pSdig->Pc()][pSdig->PadPcX()][pSdig->PadPcY()])*AliHMPIDParam::Nsig()));
142       if(AliHMPIDParam::IsOverTh(q)) new((*pLst[iCh])[iCnt[iCh]++]) AliHMPIDDigit(iPad,(Int_t)q,aTids);}  //do not create digit for the very first sdigit 
143     
144     iPad=pSdig->Pad(); iCh=AliHMPIDParam::A2C(iPad);                                                            //new sdigit comes, reset collectors
145     iNdigPad=1;
146     aTids[0]=pSdig->GetTrack(0);aTids[1]=aTids[2]=-1; 
147     q=pSdig->Q();
148     if(fgDoNoise) q+=arrNoise[iCh][pSdig->Pc()][pSdig->PadPcX()][pSdig->PadPcY()];
149     arrNoise[iCh][pSdig->Pc()][pSdig->PadPcX()][pSdig->PadPcY()]=0;
150   }//sdigits loop (sorted)
151   
152   if(iCh>=AliHMPIDParam::kMinCh && iCh<=AliHMPIDParam::kMaxCh){
153     Int_t pc = AliHMPIDParam::A2P(iPad);
154     Int_t px = AliHMPIDParam::A2X(iPad); 
155     Int_t py = AliHMPIDParam::A2Y(iPad);
156     AliHMPIDParam::Instance()->SetThreshold((TMath::Nint(arrSigmaPed[iCh][pc][px][py])*AliHMPIDParam::Nsig()));
157     if(AliHMPIDParam::IsOverTh(q)) new((*pLst[iCh])[iCnt[iCh]++]) AliHMPIDDigit(iPad,(Int_t)q,aTids);
158   }  //add the last one, in case of empty sdigits list q=-1 
159   
160 // add noise pad above threshold with no signal merged...if any
161   if(!fgDoNoise) return;
162   
163   aTids[0]=aTids[1]=aTids[2]=-1;
164   for (Int_t iChCurr=AliHMPIDParam::kMinCh;iChCurr<=AliHMPIDParam::kMaxCh;iChCurr++){
165     for (Int_t iPc=AliHMPIDParam::kMinPc;iPc<=AliHMPIDParam::kMaxPc;iPc++)
166       for(Int_t iPx=AliHMPIDParam::kMinPx;iPx<=AliHMPIDParam::kMaxPx;iPx++)
167         for(Int_t iPy=AliHMPIDParam::kMinPy;iPy<=AliHMPIDParam::kMaxPy;iPy++) {
168           Float_t qNoise = arrNoise[iChCurr][iPc][iPx][iPy];
169           AliHMPIDParam::Instance()->SetThreshold((TMath::Nint(arrSigmaPed[iChCurr][iPc][iPx][iPy])*AliHMPIDParam::Nsig()));
170           if(AliHMPIDParam::IsOverTh(qNoise)) new((*pLst[iChCurr])[iCnt[iChCurr]++]) AliHMPIDDigit(AliHMPIDParam::Abs(iChCurr,iPc,iPx,iPy),(Int_t)qNoise,aTids);
171         }
172   }        
173 }//Sdi2Dig()
174 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++