ROOT logo
/**************************************************************************
 * Copyright(c) 1998-2000, ALICE Experiment at CERN, All rights reserved. *
 *                                                                        *
 * Author: The ALICE Off-line Project.                                    *
 * Contributors are mentioned in the code where appropriate.              *
 *                                                                        *
 * Permission to use, copy, modify and distribute this software and its   *
 * documentation strictly for non-commercial purposes is hereby granted   *
 * without fee, provided that the above copyright notice appears in all   *
 * copies and that both the copyright notice and this permission notice   *
 * appear in the supporting documentation. The authors make no claims     *
 * about the suitability of this software for any purpose. It is          *
 * provided "as is" without express or implied warranty.                  *
 **************************************************************************/
//.
//.
//.
//.
//.
#include <AliRun.h>
#include <AliRunLoader.h>
#include "AliDigitizationInput.h"
#include <AliLoader.h>
#include "AliLog.h"
#include <AliCDBEntry.h> 
#include <AliCDBManager.h>
#include "AliHMPIDDigitizer.h"
#include "AliHMPIDReconstructor.h"
#include "AliHMPIDDigit.h"
#include "AliHMPID.h"
#include "AliHMPIDParam.h"
#include <TRandom.h>
#include <TMath.h>
#include <TTree.h>
#include <TObjArray.h>

ClassImp(AliHMPIDDigitizer)

Bool_t AliHMPIDDigitizer::fgDoNoise=kTRUE;
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
void AliHMPIDDigitizer::Digitize(Option_t*)
{
// This methode is responsible for merging sdigits to a list of digits
//Disintegration leeds to the fact that one hit affects several neighbouring pads, which means that the same pad might be affected by few hits.     
  AliDebug(1,Form("Start with %i input(s) for event %i",fDigInput->GetNinputs(),fDigInput->GetOutputEventNr()));
//First we read all sdigits from all inputs  
  AliRunLoader *pInRunLoader=0;//in and out Run loaders
  AliLoader    *pInRichLoader=0;//in and out HMPID loaders  
  static TClonesArray sdigs("AliHMPIDDigit");//tmp storage for sdigits sum up from all input files
  Int_t total=0;
  for(Int_t inFileN=0;inFileN<fDigInput->GetNinputs();inFileN++){//files loop
    pInRunLoader  = AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(inFileN));          //get run loader from current input 
    pInRichLoader = pInRunLoader->GetLoader("HMPIDLoader"); if(pInRichLoader==0) continue;       //no HMPID in this input, check the next input
    if (!pInRunLoader->GetAliRun()) pInRunLoader->LoadgAlice();
    AliHMPID* pInRich=(AliHMPID*)pInRunLoader->GetAliRun()->GetDetector("HMPID");                  //take HMPID from current input
    pInRichLoader->LoadSDigits(); pInRichLoader->TreeS()->GetEntry(0);                          //take list of HMPID sdigits from current input 
    AliDebug(1,Form("input %i has %i sdigits",inFileN,pInRich->SdiLst()->GetEntries()));
    for(Int_t i=0;i<pInRich->SdiLst()->GetEntries();i++){                                        //collect sdigits from current input
      AliHMPIDDigit *pSDig=(AliHMPIDDigit*)pInRich->SdiLst()->At(i);
      pSDig->AddTidOffset(fDigInput->GetMask(inFileN));                                          //apply TID shift since all inputs count tracks independently starting from 0
      new(sdigs[total++]) AliHMPIDDigit(*pSDig);       
    }
    pInRichLoader->UnloadSDigits();   pInRich->SdiReset(); //close current input and reset 
  }//files loop

  //PH  if(sdigs.GetEntries()==0) return;                                                              //no sdigits collected, nothing to convert  
  
  AliRunLoader *pOutRunLoader  = AliRunLoader::GetRunLoader(fDigInput->GetOutputFolderName());    //open output stream (only 1 possible)
  AliLoader    *pOutRichLoader = pOutRunLoader->GetLoader("HMPIDLoader");                         //take output HMPID loader
  AliRun *pArun = pOutRunLoader->GetAliRun();
  AliHMPID      *pOutRich       = (AliHMPID*)pArun->GetDetector("HMPID");      //take output HMPID
  pOutRichLoader->MakeTree("D");   pOutRich->MakeBranch("D");                                    //create TreeD in output stream

  Sdi2Dig(&sdigs,pOutRich->DigLst());
  
  pOutRichLoader->TreeD()->Fill();              //fill the output tree with the list of digits
  pOutRichLoader->WriteDigits("OVERWRITE");     //serialize them to file
  
  sdigs.Clear();                      //remove all tmp sdigits
  pOutRichLoader->UnloadDigits();   pOutRich->DigReset();
}//Exec()
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
void AliHMPIDDigitizer::Sdi2Dig(TClonesArray *pSdiLst,TObjArray *pDigLst)
{
// Converts list of sdigits to 7 lists of digits, one per each chamber
// Arguments: pSDigLst - list of all sdigits
//            pDigLst  - list of 7 lists of digits        
//   Returns: none  
  
  TClonesArray *pLst[7]; Int_t iCnt[7];

  for(Int_t i=0;i<7;i++){
    pLst[i]=(TClonesArray*)(*pDigLst)[i];
    iCnt[i]=0; if(pLst[i]->GetEntries()!=0) AliErrorClass("Some of digits lists is not empty");         //in principle those lists should be empty                                                                       
  }
  
  // make noise array
  Float_t arrNoise[7][6][80][48], arrSigmaPed[7][6][80][48];
  if(fgDoNoise) {
    for (Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++)
      for (Int_t iPc=AliHMPIDParam::kMinPc;iPc<=AliHMPIDParam::kMaxPc;iPc++)
        for(Int_t iPx=AliHMPIDParam::kMinPx;iPx<=AliHMPIDParam::kMaxPx;iPx++)
          for(Int_t iPy=AliHMPIDParam::kMinPy;iPy<=AliHMPIDParam::kMaxPy;iPy++){
            arrNoise[iCh][iPc][iPx][iPy] = gRandom->Gaus(0,1.);
            arrSigmaPed[iCh][iPc][iPx][iPy] = 1.;
          }
          
    AliCDBEntry *pDaqSigEnt = AliCDBManager::Instance()->Get("HMPID/Calib/DaqSig");  //contains TObjArray of TObjArray 14 TMatrixF sigmas values for pads 
   
    if(pDaqSigEnt){
      TObjArray *pDaqSig = (TObjArray*)pDaqSigEnt->GetObject();
      for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){                  //chambers loop    
	TMatrixF *pM = (TMatrixF*)pDaqSig->At(iCh);     
	for (Int_t iPc=AliHMPIDParam::kMinPc;iPc<=AliHMPIDParam::kMaxPc;iPc++)
	  for(Int_t iPx=AliHMPIDParam::kMinPx;iPx<=AliHMPIDParam::kMaxPx;iPx++)
	    for(Int_t iPy=AliHMPIDParam::kMinPy;iPy<=AliHMPIDParam::kMaxPy;iPy++){
	      Int_t padX = (iPc%2)*AliHMPIDParam::kPadPcX+iPx; 
	      Int_t padY = (iPc/2)*AliHMPIDParam::kPadPcY+iPy;
	      if((*pM)(padX,padY)>0.){
		arrNoise[iCh][iPc][iPx][iPy] = gRandom->Gaus(0,(*pM)(padX,padY));
		arrSigmaPed[iCh][iPc][iPx][iPy] = (*pM)(padX,padY);}
	      else{
		arrNoise[iCh][iPc][iPx][iPy] = gRandom->Gaus(0,1.);
		arrSigmaPed[iCh][iPc][iPx][iPy] = 1.;}
	    } 
      }
    }
  }  
  
  pSdiLst->Sort();  
                     
  Int_t iPad=-1,iCh=-1,iNdigPad=-1,aTids[3]={-1,-1,-1}; Float_t q=-1;
  for(Int_t i=0;i<pSdiLst->GetEntries();i++){                                                                  //sdigits loop (sorted)
    AliHMPIDDigit *pSdig=(AliHMPIDDigit*)pSdiLst->At(i);                                                       //take current sdigit
    if(pSdig->Pad()==iPad){                                                                                    //if the same pad 
      q+=pSdig->Q();                                                                                           //sum up charge
      iNdigPad++; if(iNdigPad<=3) aTids[iNdigPad-1]=pSdig->GetTrack(0);                                        //collect TID 
      continue;
    }
    if(i!=0 && iCh>=AliHMPIDParam::kMinCh && iCh<=AliHMPIDParam::kMaxCh){
      AliHMPIDParam::Instance()->SetThreshold((TMath::Nint(arrSigmaPed[iCh][pSdig->Pc()][pSdig->PadPcX()][pSdig->PadPcY()])*AliHMPIDParam::Nsig()));
      if(AliHMPIDParam::IsOverTh(q)) new((*pLst[iCh])[iCnt[iCh]++]) AliHMPIDDigit(iPad,(Int_t)q,aTids);}  //do not create digit for the very first sdigit 
    
    iPad=pSdig->Pad(); iCh=AliHMPIDParam::A2C(iPad);                                                            //new sdigit comes, reset collectors
    iNdigPad=1;
    aTids[0]=pSdig->GetTrack(0);aTids[1]=aTids[2]=-1; 
    q=pSdig->Q();
    if(fgDoNoise) q+=arrNoise[iCh][pSdig->Pc()][pSdig->PadPcX()][pSdig->PadPcY()];
    arrNoise[iCh][pSdig->Pc()][pSdig->PadPcX()][pSdig->PadPcY()]=0;
  }//sdigits loop (sorted)
  
  if(iCh>=AliHMPIDParam::kMinCh && iCh<=AliHMPIDParam::kMaxCh){
    Int_t pc = AliHMPIDParam::A2P(iPad);
    Int_t px = AliHMPIDParam::A2X(iPad); 
    Int_t py = AliHMPIDParam::A2Y(iPad);
    AliHMPIDParam::Instance()->SetThreshold((TMath::Nint(arrSigmaPed[iCh][pc][px][py])*AliHMPIDParam::Nsig()));
    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 
  
// add noise pad above threshold with no signal merged...if any
  if(!fgDoNoise) return;
  
  aTids[0]=aTids[1]=aTids[2]=-1;
  for (Int_t iChCurr=AliHMPIDParam::kMinCh;iChCurr<=AliHMPIDParam::kMaxCh;iChCurr++){
    for (Int_t iPc=AliHMPIDParam::kMinPc;iPc<=AliHMPIDParam::kMaxPc;iPc++)
      for(Int_t iPx=AliHMPIDParam::kMinPx;iPx<=AliHMPIDParam::kMaxPx;iPx++)
        for(Int_t iPy=AliHMPIDParam::kMinPy;iPy<=AliHMPIDParam::kMaxPy;iPy++) {
          Float_t qNoise = arrNoise[iChCurr][iPc][iPx][iPy];
          AliHMPIDParam::Instance()->SetThreshold((TMath::Nint(arrSigmaPed[iChCurr][iPc][iPx][iPy])*AliHMPIDParam::Nsig()));
          if(AliHMPIDParam::IsOverTh(qNoise)) new((*pLst[iChCurr])[iCnt[iChCurr]++]) AliHMPIDDigit(AliHMPIDParam::Abs(iChCurr,iPc,iPx,iPy),(Int_t)qNoise,aTids);
        }
  }        
}//Sdi2Dig()
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
 AliHMPIDDigitizer.cxx:1
 AliHMPIDDigitizer.cxx:2
 AliHMPIDDigitizer.cxx:3
 AliHMPIDDigitizer.cxx:4
 AliHMPIDDigitizer.cxx:5
 AliHMPIDDigitizer.cxx:6
 AliHMPIDDigitizer.cxx:7
 AliHMPIDDigitizer.cxx:8
 AliHMPIDDigitizer.cxx:9
 AliHMPIDDigitizer.cxx:10
 AliHMPIDDigitizer.cxx:11
 AliHMPIDDigitizer.cxx:12
 AliHMPIDDigitizer.cxx:13
 AliHMPIDDigitizer.cxx:14
 AliHMPIDDigitizer.cxx:15
 AliHMPIDDigitizer.cxx:16
 AliHMPIDDigitizer.cxx:17
 AliHMPIDDigitizer.cxx:18
 AliHMPIDDigitizer.cxx:19
 AliHMPIDDigitizer.cxx:20
 AliHMPIDDigitizer.cxx:21
 AliHMPIDDigitizer.cxx:22
 AliHMPIDDigitizer.cxx:23
 AliHMPIDDigitizer.cxx:24
 AliHMPIDDigitizer.cxx:25
 AliHMPIDDigitizer.cxx:26
 AliHMPIDDigitizer.cxx:27
 AliHMPIDDigitizer.cxx:28
 AliHMPIDDigitizer.cxx:29
 AliHMPIDDigitizer.cxx:30
 AliHMPIDDigitizer.cxx:31
 AliHMPIDDigitizer.cxx:32
 AliHMPIDDigitizer.cxx:33
 AliHMPIDDigitizer.cxx:34
 AliHMPIDDigitizer.cxx:35
 AliHMPIDDigitizer.cxx:36
 AliHMPIDDigitizer.cxx:37
 AliHMPIDDigitizer.cxx:38
 AliHMPIDDigitizer.cxx:39
 AliHMPIDDigitizer.cxx:40
 AliHMPIDDigitizer.cxx:41
 AliHMPIDDigitizer.cxx:42
 AliHMPIDDigitizer.cxx:43
 AliHMPIDDigitizer.cxx:44
 AliHMPIDDigitizer.cxx:45
 AliHMPIDDigitizer.cxx:46
 AliHMPIDDigitizer.cxx:47
 AliHMPIDDigitizer.cxx:48
 AliHMPIDDigitizer.cxx:49
 AliHMPIDDigitizer.cxx:50
 AliHMPIDDigitizer.cxx:51
 AliHMPIDDigitizer.cxx:52
 AliHMPIDDigitizer.cxx:53
 AliHMPIDDigitizer.cxx:54
 AliHMPIDDigitizer.cxx:55
 AliHMPIDDigitizer.cxx:56
 AliHMPIDDigitizer.cxx:57
 AliHMPIDDigitizer.cxx:58
 AliHMPIDDigitizer.cxx:59
 AliHMPIDDigitizer.cxx:60
 AliHMPIDDigitizer.cxx:61
 AliHMPIDDigitizer.cxx:62
 AliHMPIDDigitizer.cxx:63
 AliHMPIDDigitizer.cxx:64
 AliHMPIDDigitizer.cxx:65
 AliHMPIDDigitizer.cxx:66
 AliHMPIDDigitizer.cxx:67
 AliHMPIDDigitizer.cxx:68
 AliHMPIDDigitizer.cxx:69
 AliHMPIDDigitizer.cxx:70
 AliHMPIDDigitizer.cxx:71
 AliHMPIDDigitizer.cxx:72
 AliHMPIDDigitizer.cxx:73
 AliHMPIDDigitizer.cxx:74
 AliHMPIDDigitizer.cxx:75
 AliHMPIDDigitizer.cxx:76
 AliHMPIDDigitizer.cxx:77
 AliHMPIDDigitizer.cxx:78
 AliHMPIDDigitizer.cxx:79
 AliHMPIDDigitizer.cxx:80
 AliHMPIDDigitizer.cxx:81
 AliHMPIDDigitizer.cxx:82
 AliHMPIDDigitizer.cxx:83
 AliHMPIDDigitizer.cxx:84
 AliHMPIDDigitizer.cxx:85
 AliHMPIDDigitizer.cxx:86
 AliHMPIDDigitizer.cxx:87
 AliHMPIDDigitizer.cxx:88
 AliHMPIDDigitizer.cxx:89
 AliHMPIDDigitizer.cxx:90
 AliHMPIDDigitizer.cxx:91
 AliHMPIDDigitizer.cxx:92
 AliHMPIDDigitizer.cxx:93
 AliHMPIDDigitizer.cxx:94
 AliHMPIDDigitizer.cxx:95
 AliHMPIDDigitizer.cxx:96
 AliHMPIDDigitizer.cxx:97
 AliHMPIDDigitizer.cxx:98
 AliHMPIDDigitizer.cxx:99
 AliHMPIDDigitizer.cxx:100
 AliHMPIDDigitizer.cxx:101
 AliHMPIDDigitizer.cxx:102
 AliHMPIDDigitizer.cxx:103
 AliHMPIDDigitizer.cxx:104
 AliHMPIDDigitizer.cxx:105
 AliHMPIDDigitizer.cxx:106
 AliHMPIDDigitizer.cxx:107
 AliHMPIDDigitizer.cxx:108
 AliHMPIDDigitizer.cxx:109
 AliHMPIDDigitizer.cxx:110
 AliHMPIDDigitizer.cxx:111
 AliHMPIDDigitizer.cxx:112
 AliHMPIDDigitizer.cxx:113
 AliHMPIDDigitizer.cxx:114
 AliHMPIDDigitizer.cxx:115
 AliHMPIDDigitizer.cxx:116
 AliHMPIDDigitizer.cxx:117
 AliHMPIDDigitizer.cxx:118
 AliHMPIDDigitizer.cxx:119
 AliHMPIDDigitizer.cxx:120
 AliHMPIDDigitizer.cxx:121
 AliHMPIDDigitizer.cxx:122
 AliHMPIDDigitizer.cxx:123
 AliHMPIDDigitizer.cxx:124
 AliHMPIDDigitizer.cxx:125
 AliHMPIDDigitizer.cxx:126
 AliHMPIDDigitizer.cxx:127
 AliHMPIDDigitizer.cxx:128
 AliHMPIDDigitizer.cxx:129
 AliHMPIDDigitizer.cxx:130
 AliHMPIDDigitizer.cxx:131
 AliHMPIDDigitizer.cxx:132
 AliHMPIDDigitizer.cxx:133
 AliHMPIDDigitizer.cxx:134
 AliHMPIDDigitizer.cxx:135
 AliHMPIDDigitizer.cxx:136
 AliHMPIDDigitizer.cxx:137
 AliHMPIDDigitizer.cxx:138
 AliHMPIDDigitizer.cxx:139
 AliHMPIDDigitizer.cxx:140
 AliHMPIDDigitizer.cxx:141
 AliHMPIDDigitizer.cxx:142
 AliHMPIDDigitizer.cxx:143
 AliHMPIDDigitizer.cxx:144
 AliHMPIDDigitizer.cxx:145
 AliHMPIDDigitizer.cxx:146
 AliHMPIDDigitizer.cxx:147
 AliHMPIDDigitizer.cxx:148
 AliHMPIDDigitizer.cxx:149
 AliHMPIDDigitizer.cxx:150
 AliHMPIDDigitizer.cxx:151
 AliHMPIDDigitizer.cxx:152
 AliHMPIDDigitizer.cxx:153
 AliHMPIDDigitizer.cxx:154
 AliHMPIDDigitizer.cxx:155
 AliHMPIDDigitizer.cxx:156
 AliHMPIDDigitizer.cxx:157
 AliHMPIDDigitizer.cxx:158
 AliHMPIDDigitizer.cxx:159
 AliHMPIDDigitizer.cxx:160
 AliHMPIDDigitizer.cxx:161
 AliHMPIDDigitizer.cxx:162
 AliHMPIDDigitizer.cxx:163
 AliHMPIDDigitizer.cxx:164
 AliHMPIDDigitizer.cxx:165
 AliHMPIDDigitizer.cxx:166
 AliHMPIDDigitizer.cxx:167
 AliHMPIDDigitizer.cxx:168
 AliHMPIDDigitizer.cxx:169
 AliHMPIDDigitizer.cxx:170
 AliHMPIDDigitizer.cxx:171
 AliHMPIDDigitizer.cxx:172
 AliHMPIDDigitizer.cxx:173
 AliHMPIDDigitizer.cxx:174