]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDDigitizer.cxx
Geometry V2 updated with the final frames. Now completed.
[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 "AliRunDigitizer.h"
23 #include <AliLoader.h>
24 #include <AliLog.h>
25 #include "AliHMPIDDigitizer.h"
26 #include "AliHMPIDDigit.h"
27 #include "AliHMPID.h"
28 #include "AliHMPIDParam.h"
29 #include <TRandom.h>
30
31 ClassImp(AliHMPIDDigitizer)
32
33 Bool_t AliHMPIDDigitizer::fgDoNoise=kFALSE;
34 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
35 void AliHMPIDDigitizer::Exec(Option_t*)
36 {
37 // This methode is responsible for merging sdigits to a list of digits
38 //Disintegration leeds to the fact that one hit affects several neighbouring pads, which means that the same pad might be affected by few hits.     
39   AliDebug(1,Form("Start with %i input(s) for event %i",fManager->GetNinputs(),fManager->GetOutputEventNr()));
40 //First we read all sdigits from all inputs  
41   AliRunLoader *pInRunLoader=0;//in and out Run loaders
42   AliLoader    *pInRichLoader=0;//in and out HMPID loaders  
43   TClonesArray sdigs("AliHMPIDDigit");//tmp storage for sdigits sum up from all input files
44   Int_t total=0;
45   for(Int_t inFileN=0;inFileN<fManager->GetNinputs();inFileN++){//files loop
46     pInRunLoader  = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(inFileN));          //get run loader from current input 
47     pInRichLoader = pInRunLoader->GetLoader("HMPIDLoader"); if(pInRichLoader==0) continue;       //no HMPID in this input, check the next input
48     if (!pInRunLoader->GetAliRun()) pInRunLoader->LoadgAlice();
49     AliHMPID* pInRich=(AliHMPID*)pInRunLoader->GetAliRun()->GetDetector("HMPID");                  //take HMPID from current input
50     pInRichLoader->LoadSDigits(); pInRichLoader->TreeS()->GetEntry(0);                          //take list of HMPID sdigits from current input 
51     AliDebug(1,Form("input %i has %i sdigits",inFileN,pInRich->SdiLst()->GetEntries()));
52     for(Int_t i=0;i<pInRich->SdiLst()->GetEntries();i++){                                        //collect sdigits from current input
53       AliHMPIDDigit *pSDig=(AliHMPIDDigit*)pInRich->SdiLst()->At(i);
54       pSDig->AddTidOffset(fManager->GetMask(inFileN));                                          //apply TID shift since all inputs count tracks independently starting from 0
55       new(sdigs[total++]) AliHMPIDDigit(*pSDig);       
56     }
57     pInRichLoader->UnloadSDigits();   pInRich->SdiReset(); //close current input and reset 
58   }//files loop
59
60   //PH  if(sdigs.GetEntries()==0) return;                                                              //no sdigits collected, nothing to convert  
61   
62   AliRunLoader *pOutRunLoader  = AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());    //open output stream (only 1 possible)
63   AliLoader    *pOutRichLoader = pOutRunLoader->GetLoader("HMPIDLoader");                         //take output HMPID loader
64   AliRun *pArun = pOutRunLoader->GetAliRun();
65   AliHMPID      *pOutRich       = (AliHMPID*)pArun->GetDetector("HMPID");      //take output HMPID
66   pOutRichLoader->MakeTree("D");   pOutRich->MakeBranch("D");                                    //create TreeD in output stream
67
68   Sdi2Dig(&sdigs,pOutRich->DigLst());
69   
70   pOutRichLoader->TreeD()->Fill();              //fill the output tree with the list of digits
71   pOutRichLoader->WriteDigits("OVERWRITE");     //serialize them to file
72   
73   sdigs.Clear();                      //remove all tmp sdigits
74   pOutRichLoader->UnloadDigits();   pOutRich->DigReset();
75 }//Exec()
76 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
77 void AliHMPIDDigitizer::Sdi2Dig(TClonesArray *pSdiLst,TObjArray *pDigLst)
78 {
79 // Converts list of sdigits to 7 lists of digits, one per each chamber
80 // Arguments: pSDigLst - list of all sdigits
81 //            pDigLst  - list of 7 lists of digits        
82 //   Returns: none  
83   
84   TClonesArray *pLst[7]; Int_t iCnt[7];
85
86   TRandom *rnd = new TRandom();  
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                                                                       
90   }
91
92   // make noise array
93   Float_t arrNoise[7][6][80][48];
94   if(fgDoNoise) {
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] = rnd->Gaus(0,1);
99   }  
100   
101   pSdiLst->Sort();  
102                      
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 
109       continue;
110     }
111     if(i!=0 && AliHMPIDParam::IsOverTh(q))  new((*pLst[iCh])[iCnt[iCh]++]) AliHMPIDDigit(iPad,(Int_t)q,aTids);   //do not create digit for the very first sdigit 
112     iPad=pSdig->Pad(); iCh=AliHMPIDParam::A2C(iPad);                                                            //new sdigit comes, reset collectors
113     iNdigPad=1;
114     aTids[0]=pSdig->GetTrack(0);aTids[1]=aTids[2]=-1; 
115     q=pSdig->Q();
116     if(fgDoNoise) q+=arrNoise[iCh][pSdig->Pc()][pSdig->PadPcX()][pSdig->PadPcY()];
117     arrNoise[iCh][pSdig->Pc()][pSdig->PadPcX()][pSdig->PadPcY()]=0;
118   }//sdigits loop (sorted)
119   
120   if(AliHMPIDParam::IsOverTh(q))  new((*pLst[iCh])[iCnt[iCh]++]) AliHMPIDDigit(iPad,(Int_t)q,aTids);           //add the last one, in case of empty sdigits list q=-1
121 // add noise pad above threshold with no signal merged...if any
122   if(!fgDoNoise) return;
123   aTids[0]=aTids[1]=aTids[2]=-1;
124   for (Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++)
125     for (Int_t iPc=AliHMPIDParam::kMinPc;iPc<=AliHMPIDParam::kMaxPc;iPc++)
126       for(Int_t iPx=AliHMPIDParam::kMinPx;iPx<=AliHMPIDParam::kMaxPx;iPx++)
127         for(Int_t iPy=AliHMPIDParam::kMinPy;iPy<=AliHMPIDParam::kMaxPy;iPy++) {
128           Float_t q = arrNoise[iCh][iPc][iPx][iPy];
129           if(AliHMPIDParam::IsOverTh(q)) new((*pLst[iCh])[iCnt[iCh]++]) AliHMPIDDigit(AliHMPIDParam::Abs(iCh,iPc,iPx,iPy),(Int_t)q,aTids);
130         }
131         
132 }//Sdi2Dig()
133 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++