]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
First version of the SDD DA calibration classes. AliITSOnlineSDDBase - for measuremen...
authorcvetan <cvetan@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 23 May 2007 11:36:05 +0000 (11:36 +0000)
committercvetan <cvetan@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 23 May 2007 11:36:05 +0000 (11:36 +0000)
12 files changed:
ITS/AliITSOnlineSDD.cxx [new file with mode: 0644]
ITS/AliITSOnlineSDD.h [new file with mode: 0644]
ITS/AliITSOnlineSDDBTP.cxx [new file with mode: 0644]
ITS/AliITSOnlineSDDBTP.h [new file with mode: 0644]
ITS/AliITSOnlineSDDBase.cxx [new file with mode: 0644]
ITS/AliITSOnlineSDDBase.h [new file with mode: 0644]
ITS/AliITSOnlineSDDInjectors.cxx [new file with mode: 0644]
ITS/AliITSOnlineSDDInjectors.h [new file with mode: 0644]
ITS/AliITSOnlineSDDTP.cxx [new file with mode: 0644]
ITS/AliITSOnlineSDDTP.h [new file with mode: 0644]
ITS/ITSrecLinkDef.h
ITS/libITSrec.pkg

diff --git a/ITS/AliITSOnlineSDD.cxx b/ITS/AliITSOnlineSDD.cxx
new file mode 100644 (file)
index 0000000..893f8da
--- /dev/null
@@ -0,0 +1,37 @@
+/**************************************************************************
+ * Copyright(c) 2007-2009, 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.                  *
+ **************************************************************************/
+
+///////////////////////////////////////////////////////////////////
+//                                                               //
+// Implementation of the base class for SDD detector algorithms  //
+// Origin: F.Prino, Torino, prino@to.infn.it                     //
+//                                                               //
+///////////////////////////////////////////////////////////////////
+
+#include "AliITSOnlineSDD.h"
+
+ClassImp(AliITSOnlineSDD)
+//______________________________________________________________________
+  AliITSOnlineSDD::AliITSOnlineSDD():TObject(),fModuleId(0),fSide(0)
+{
+  // default constructor
+}
+//______________________________________________________________________
+  AliITSOnlineSDD::AliITSOnlineSDD(Int_t mod, Int_t sid):TObject(),fModuleId(0),fSide(0)
+{
+  // standard constructor
+  SetModule(mod);
+  SetDetectorSide(sid);
+}
diff --git a/ITS/AliITSOnlineSDD.h b/ITS/AliITSOnlineSDD.h
new file mode 100644 (file)
index 0000000..41911ac
--- /dev/null
@@ -0,0 +1,35 @@
+#ifndef ALIITSONLINESDD_H
+#define ALIITSONLINESDD_H
+
+
+///////////////////////////////////////////////////////////////////
+//                                                               //
+// Base class for SDD detector algorithms                        //
+// Origin: F.Prino, Torino, prino@to.infn.it                     //
+//                                                               //
+///////////////////////////////////////////////////////////////////
+
+#include<TObject.h>
+#include<TMath.h>
+
+class AliITSOnlineSDD : public TObject {
+
+ public:
+  AliITSOnlineSDD();
+  AliITSOnlineSDD(Int_t mod, Int_t sid);
+  virtual ~AliITSOnlineSDD(){};
+
+  void SetModule(Int_t mod){fModuleId=mod;}
+  void SetDetectorSide(Int_t sid){fSide=sid;}
+
+  Int_t GetModuleId() const {return fModuleId;}
+  Int_t GetDetectorSide() const {return fSide;}
+
+ protected:
+  static const Int_t fgkNAnodes = 256;
+  Int_t fModuleId; // module number from 0 to 255
+  Int_t fSide;     // detector side (0-1)
+
+  ClassDef(AliITSOnlineSDD,1);
+};
+#endif
diff --git a/ITS/AliITSOnlineSDDBTP.cxx b/ITS/AliITSOnlineSDDBTP.cxx
new file mode 100644 (file)
index 0000000..969826e
--- /dev/null
@@ -0,0 +1,209 @@
+/**************************************************************************
+ * Copyright(c) 2007-2009, 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 "AliITSOnlineSDDBTP.h"
+#include <TH2F.h>
+#include <Riostream.h>
+#include <TMath.h>
+
+
+///////////////////////////////////////////////////////////////////
+//                                                               //
+// Class used for SDD baseline, noise and gain analysis          //
+// Origin: F.Prino, Torino, prino@to.infn.it                     //
+//                                                               //
+///////////////////////////////////////////////////////////////////
+
+
+ClassImp(AliITSOnlineSDDBTP)
+//______________________________________________________________________
+  AliITSOnlineSDDBTP::AliITSOnlineSDDBTP():AliITSOnlineSDD(),fNBaseEvents(0),fNTPEvents(0),fMinBaseline(0.),fMaxBaseline(0.),fMinRawNoise(0.),fMaxRawNoise(0.),fNSigmaNoise(0.),fNSigmaGain(0.)
+{
+  // default constructor
+  Reset();
+  SetMinBaseline();
+  SetMaxBaseline();
+  SetMinRawNoise();
+  SetMaxRawNoise();
+  SetNSigmaNoise();
+  SetNSigmaGain();
+}
+//______________________________________________________________________
+  AliITSOnlineSDDBTP::AliITSOnlineSDDBTP(Int_t mod, Int_t sid):AliITSOnlineSDD(mod,sid),fNBaseEvents(0),fNTPEvents(0),fMinBaseline(0.),fMaxBaseline(0.),fMinRawNoise(0.),fMaxRawNoise(0.),fNSigmaNoise(0.),fNSigmaGain(0.)
+{
+  // default constructor
+  Reset();
+  SetMinBaseline();
+  SetMaxBaseline();
+  SetMinRawNoise();
+  SetMaxRawNoise();
+  SetNSigmaNoise();
+  SetNSigmaGain();
+}
+//______________________________________________________________________
+AliITSOnlineSDDBTP::~AliITSOnlineSDDBTP(){
+  // Destructor
+}
+//______________________________________________________________________
+void AliITSOnlineSDDBTP::Reset(){
+  fNBaseEvents=0;
+  fNTPEvents=0;
+  for(Int_t i=0;i<fgkNAnodes;i++){
+    fGoodAnode[i]=1;
+    fSumBaseline[i]=0.;
+    fSumRawNoise[i]=0.;
+    fSumCMN[i]=0.;
+    fSumTPPeak[i]=0.;
+    fTPPos[i]=0.;
+  }
+}
+//______________________________________________________________________
+void  AliITSOnlineSDDBTP::ValidateAnodes(){
+  for(Int_t ian=0;ian<fgkNAnodes;ian++){
+    fGoodAnode[ian]=1;
+    if(GetAnodeBaseline(ian)>fMaxBaseline || GetAnodeBaseline(ian)<fMinBaseline) fGoodAnode[ian]=0;
+    if(GetAnodeRawNoise(ian)>fMaxRawNoise || GetAnodeRawNoise(ian)<fMinRawNoise) fGoodAnode[ian]=0;
+    if(GetAnodeRawNoise(ian)>fNSigmaNoise*CalcMeanRawNoise()) fGoodAnode[ian]=0;
+  }
+  if(fNTPEvents>0){
+    Float_t meang,rmsg;
+    StatGain(meang,rmsg);
+    Float_t lowlim=meang-fNSigmaGain*rmsg;
+    Float_t hilim=meang+fNSigmaGain*rmsg;
+    for(Int_t ian=0;ian<fgkNAnodes;ian++){
+      if(!fGoodAnode[ian]) continue;
+      if(GetChannelGain(ian)<lowlim||GetChannelGain(ian)>hilim) fGoodAnode[ian]=0;
+    }
+  }
+}
+
+//______________________________________________________________________
+void AliITSOnlineSDDBTP::AddTPEvent(TH2F* hrawd, Float_t xDAC){
+  // 
+  if(fNBaseEvents==0) return;
+  fNTPEvents++;
+  Float_t tbmax=(Float_t)hrawd->GetNbinsX();
+  for(Int_t ian=0;ian<fgkNAnodes;ian++){
+    Float_t auxmax=0.;
+    Int_t auxtb=0;
+    if(!fGoodAnode[ian]) continue;
+    for(Int_t itb=0;itb<tbmax;itb++){
+      Float_t cnt=hrawd->GetBinContent(itb+1,ian+1);
+      if(cnt>auxmax){ 
+       auxmax=cnt;
+       auxtb=itb;
+      }
+    }
+    fSumTPPeak[ian]+=(auxmax-GetAnodeBaseline(ian))/xDAC;
+    fTPPos[ian]+=auxtb;
+  }
+}
+//______________________________________________________________________
+void AliITSOnlineSDDBTP::AddBaseEvent(TH2F* hrawd){
+  // 
+  fNBaseEvents++;
+  Float_t tbmax=(Float_t)hrawd->GetNbinsX();
+  Float_t sum[fgkNAnodes];
+  for(Int_t ian=0;ian<fgkNAnodes;ian++){
+    Float_t sumQ=0.;
+    sum[ian]=0.;
+    for(Int_t itb=0;itb<tbmax;itb++){
+      sum[ian]+=hrawd->GetBinContent(itb+1,ian+1);
+      sumQ+=TMath::Power(hrawd->GetBinContent(itb+1,ian+1),2);      
+    }
+    sum[ian]/=tbmax;
+    sumQ/=tbmax;
+    fSumBaseline[ian]+=sum[ian];
+    fSumRawNoise[ian]+=sumQ;
+    if(fNBaseEvents==1) ValidateAnodes();
+  }
+
+
+  const Int_t itbmax=int(tbmax);
+  Float_t *cmnEven = new Float_t[itbmax];
+  Float_t *cmnOdd  = new Float_t[itbmax];
+  for(Int_t itb=0;itb<tbmax;itb++){
+    Float_t sumEven=0., sumOdd=0.;
+    Int_t countEven=0,countOdd=0;
+    for(Int_t ian=0;ian<fgkNAnodes;ian+=2){
+      if(!fGoodAnode[ian]) continue;
+      sumEven+=hrawd->GetBinContent(itb+1,ian+1)-sum[ian];
+      countEven++;
+    }
+    for(Int_t ian=1;ian<fgkNAnodes;ian+=2){
+      if(!fGoodAnode[ian]) continue;
+      sumOdd+=hrawd->GetBinContent(itb+1,ian+1)-sum[ian];
+      countOdd++;
+    }
+    cmnEven[itb]=sumEven/countEven;
+    cmnOdd[itb]=sumOdd/countOdd;
+  }
+  for(Int_t ian=0;ian<fgkNAnodes;ian++){
+    Float_t num=0.,den=0.;
+    if(!fGoodAnode[ian]) continue;
+    for(Int_t itb=0;itb<tbmax;itb++){
+      Float_t cmnCoef=cmnOdd[itb];
+      if(ian%2==0) cmnCoef=cmnEven[itb];
+      num+=(hrawd->GetBinContent(itb+1,ian+1)-sum[ian])*cmnCoef;
+      den+=TMath::Power(cmnCoef,2);
+    }
+    if(den!=0) fSumCMN[ian]+=num/den;
+  }
+  delete [] cmnEven;
+  delete [] cmnOdd;
+}
+//______________________________________________________________________
+Float_t AliITSOnlineSDDBTP::CalcMeanRawNoise(){
+  //
+  Float_t meanns=0.;
+  Int_t cnt=0;
+  for(Int_t ian=0;ian<fgkNAnodes;ian++){
+    if(!fGoodAnode[ian]) continue;  
+    meanns+=GetAnodeRawNoise(ian);
+    cnt++;
+  }
+  if(cnt>0) meanns/=(Float_t)cnt;
+  return meanns;
+}
+//______________________________________________________________________
+void AliITSOnlineSDDBTP::StatGain(Float_t &mean, Float_t  &rms){
+  Float_t sum=0.,sumq=0.;
+  Int_t cnt=0;
+  for(Int_t ian=0;ian<fgkNAnodes;ian++){
+    if(!fGoodAnode[ian]) continue;
+    sum+=GetChannelGain(ian);
+    sumq+=TMath::Power(GetChannelGain(ian),2);
+    cnt++;
+  }
+  if(cnt>0){ 
+    mean=sum/(Float_t)cnt;
+    rms=TMath::Sqrt(sumq/(Float_t)cnt-mean*mean);
+  }else{ 
+    mean=0.;
+    rms=0.;
+  }
+  return;
+}
+//______________________________________________________________________
+void AliITSOnlineSDDBTP::WriteToFXS(){
+  //
+  Char_t outfilnam[100];
+  sprintf(outfilnam,"SDDbase_mod%03d_sid%d.data",fModuleId,fSide);
+  FILE* outf=fopen(outfilnam,"w");
+  Float_t corrnoise=2.;
+  for(Int_t ian=0;ian<fgkNAnodes;ian++){
+    fprintf(outf,"%d %d %8.3f %8.3f %8.3f %8.3f %8.3f\n",ian,IsAnodeGood(ian),GetAnodeBaseline(ian),GetAnodeRawNoise(ian),GetAnodeCommonMode(ian),corrnoise,GetChannelGain(ian));
+  }
+  fclose(outf);  
+}
diff --git a/ITS/AliITSOnlineSDDBTP.h b/ITS/AliITSOnlineSDDBTP.h
new file mode 100644 (file)
index 0000000..a551fe8
--- /dev/null
@@ -0,0 +1,75 @@
+#ifndef ALIITSONLINESDDBTP_H
+#define ALIITSONLINESDDBTP_H
+
+///////////////////////////////////////////////////////////////////
+//                                                               //
+// Class used for SDD baseline, noise and gain analysis          //
+// Origin: F.Prino, Torino, prino@to.infn.it                     //
+//                                                               //
+///////////////////////////////////////////////////////////////////
+#include"AliITSOnlineSDD.h"
+
+class TH2F;
+class TGraph;
+class AliITSOnlineSDDBTP : public AliITSOnlineSDD {
+
+ public:
+  AliITSOnlineSDDBTP();
+  AliITSOnlineSDDBTP(Int_t mod, Int_t sid);
+  virtual ~AliITSOnlineSDDBTP();
+  void Reset();
+  void AddBaseEvent(TH2F* hrawd);
+  void AddTPEvent(TH2F* hrawd, Float_t xDAC);
+  void ValidateAnodes();
+  void SetMinBaseline(Float_t bas=10.){fMinBaseline=bas;}
+  void SetMaxBaseline(Float_t bas=150.){fMaxBaseline=bas;}
+  void SetMinRawNoise(Float_t ns=0.001){fMinRawNoise=ns;}
+  void SetMaxRawNoise(Float_t ns=9.){fMaxRawNoise=ns;}
+  void SetNSigmaNoise(Float_t ns=4.){fNSigmaNoise=ns;}
+  void SetNSigmaGain(Float_t sig=3.){fNSigmaGain=sig;}
+  Bool_t IsAnodeGood(Int_t iAnode)const{ return fGoodAnode[iAnode];}
+  Float_t GetAnodeBaseline(Int_t iAnode) const{
+    if(fNBaseEvents>0) return fSumBaseline[iAnode]/fNBaseEvents;
+    else return 0;
+  }
+  Float_t GetAnodeRawNoise(Int_t iAnode) const{
+    if(fNBaseEvents>0) return TMath::Sqrt(fSumRawNoise[iAnode]/fNBaseEvents-TMath::Power(GetAnodeBaseline(iAnode),2));
+    
+    else return 0;
+  }
+  Float_t CalcMeanRawNoise();
+  void StatGain(Float_t &mean, Float_t  &rms);
+  Float_t GetAnodeCommonMode(Int_t iAnode) const{
+    if(fNBaseEvents>0) return fSumCMN[iAnode]/fNBaseEvents;
+    else return 0;
+  }
+  Float_t GetChannelGain(Int_t iAnode)const{
+    if(fNTPEvents>0) return fSumTPPeak[iAnode]/fNTPEvents;
+    else return 0;
+  }
+  Int_t GetNBaseEvents() const {return fNBaseEvents;}
+  Int_t GetNTPEvents() const {return fNTPEvents;}
+  void WriteToFXS();
+
+ protected:
+
+ private:
+
+  Int_t fNBaseEvents;
+  Int_t fNTPEvents;
+  Bool_t fGoodAnode[fgkNAnodes];
+  Float_t fSumBaseline[fgkNAnodes];
+  Float_t fSumRawNoise[fgkNAnodes];
+  Float_t fSumCMN[fgkNAnodes];
+  Float_t fSumTPPeak[fgkNAnodes];
+  Float_t fTPPos[fgkNAnodes];
+  Float_t fMinBaseline;
+  Float_t fMaxBaseline;
+  Float_t fMinRawNoise;
+  Float_t fMaxRawNoise;
+  Float_t fNSigmaNoise;
+  Float_t fNSigmaGain;
+
+  ClassDef(AliITSOnlineSDDBTP,1);
+};
+#endif
diff --git a/ITS/AliITSOnlineSDDBase.cxx b/ITS/AliITSOnlineSDDBase.cxx
new file mode 100644 (file)
index 0000000..c317890
--- /dev/null
@@ -0,0 +1,155 @@
+/**************************************************************************
+ * Copyright(c) 2007-2009, 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 "AliITSOnlineSDDBase.h"
+#include <TH2F.h>
+#include <TMath.h>
+
+
+///////////////////////////////////////////////////////////////////
+//                                                               //
+// Implementation of the class used for SDD baselines            //
+// and noise analysis                                            //
+// Origin: F.Prino, Torino, prino@to.infn.it                     //
+//                                                               //
+///////////////////////////////////////////////////////////////////
+
+
+ClassImp(AliITSOnlineSDDBase)
+//______________________________________________________________________
+  AliITSOnlineSDDBase::AliITSOnlineSDDBase():AliITSOnlineSDD(),fNEvents(0),fMinBaseline(0.),fMaxBaseline(0.),fMinRawNoise(0.),fMaxRawNoise(0.),fNSigmaNoise(0.)
+{
+  // default constructor
+  Reset();
+  SetMinBaseline();
+  SetMaxBaseline();
+  SetMinRawNoise();
+  SetMaxRawNoise();
+  SetNSigmaNoise();
+}
+//______________________________________________________________________
+  AliITSOnlineSDDBase::AliITSOnlineSDDBase(Int_t mod, Int_t sid):AliITSOnlineSDD(mod,sid),fNEvents(0),fMinBaseline(0.),fMaxBaseline(0.),fMinRawNoise(0.),fMaxRawNoise(0.),fNSigmaNoise(0.)
+{
+  // default constructor
+  Reset();
+  SetMinBaseline();
+  SetMaxBaseline();
+  SetMinRawNoise();
+  SetMaxRawNoise();
+  SetNSigmaNoise();
+}
+//______________________________________________________________________
+AliITSOnlineSDDBase::~AliITSOnlineSDDBase(){
+  // Destructor
+}
+//______________________________________________________________________
+void AliITSOnlineSDDBase::Reset(){
+  fNEvents=0;
+  for(Int_t i=0;i<fgkNAnodes;i++){
+    fGoodAnode[i]=1;
+    fSumBaseline[i]=0.;
+    fSumRawNoise[i]=0.;
+    fSumCMN[i]=0.;
+  }
+}
+//______________________________________________________________________
+void  AliITSOnlineSDDBase::ValidateAnodes(){
+  for(Int_t ian=0;ian<fgkNAnodes;ian++){
+    fGoodAnode[ian]=1;
+    if(GetAnodeBaseline(ian)>fMaxBaseline || GetAnodeBaseline(ian)<fMinBaseline) fGoodAnode[ian]=0;
+    if(GetAnodeRawNoise(ian)>fMaxRawNoise || GetAnodeRawNoise(ian)<fMinRawNoise) fGoodAnode[ian]=0;
+    if(GetAnodeRawNoise(ian)>fNSigmaNoise*CalcMeanRawNoise()) fGoodAnode[ian]=0;
+  }
+}
+
+//______________________________________________________________________
+void AliITSOnlineSDDBase::AddEvent(TH2F* hrawd){
+  // 
+  fNEvents++;
+  Float_t tbmax=(Float_t)hrawd->GetNbinsX();
+  Float_t sum[fgkNAnodes];
+  for(Int_t ian=0;ian<fgkNAnodes;ian++){
+    Float_t sumQ=0.;
+    sum[ian]=0.;
+    for(Int_t itb=0;itb<tbmax;itb++){
+      sum[ian]+=hrawd->GetBinContent(itb+1,ian+1);
+      sumQ+=TMath::Power(hrawd->GetBinContent(itb+1,ian+1),2);      
+    }
+    sum[ian]/=tbmax;
+    sumQ/=tbmax;
+    fSumBaseline[ian]+=sum[ian];
+    fSumRawNoise[ian]+=sumQ;
+    if(fNEvents==1) ValidateAnodes();
+  }
+
+
+  const Int_t itbmax=int(tbmax);
+  Float_t *cmnEven = new Float_t[itbmax];
+  Float_t *cmnOdd  = new Float_t[itbmax];
+  for(Int_t itb=0;itb<tbmax;itb++){
+    Float_t sumEven=0., sumOdd=0.;
+    Int_t countEven=0,countOdd=0;
+    for(Int_t ian=0;ian<fgkNAnodes;ian+=2){
+      if(!fGoodAnode[ian]) continue;
+      sumEven+=hrawd->GetBinContent(itb+1,ian+1)-sum[ian];
+      countEven++;
+    }
+    for(Int_t ian=1;ian<fgkNAnodes;ian+=2){
+      if(!fGoodAnode[ian]) continue;
+      sumOdd+=hrawd->GetBinContent(itb+1,ian+1)-sum[ian];
+      countOdd++;
+    }
+    cmnEven[itb]=sumEven/countEven;
+    cmnOdd[itb]=sumOdd/countOdd;
+  }
+  for(Int_t ian=0;ian<fgkNAnodes;ian++){
+    Float_t num=0.,den=0.;
+    if(!fGoodAnode[ian]) continue;
+    for(Int_t itb=0;itb<tbmax;itb++){
+      Float_t cmnCoef=cmnOdd[itb];
+      if(ian%2==0) cmnCoef=cmnEven[itb];
+      num+=(hrawd->GetBinContent(itb+1,ian+1)-sum[ian])*cmnCoef;
+      den+=TMath::Power(cmnCoef,2);
+    }
+    if(den!=0) fSumCMN[ian]+=num/den;
+  }
+
+  delete [] cmnEven;
+  delete [] cmnOdd;
+}
+//______________________________________________________________________
+Float_t AliITSOnlineSDDBase::CalcMeanRawNoise(){
+  //
+  Float_t meanns=0.;
+  Int_t cnt=0;
+  for(Int_t ian=0;ian<fgkNAnodes;ian++){
+    if(!fGoodAnode[ian]) continue;  
+    meanns+=GetAnodeRawNoise(ian);
+    cnt++;
+  }
+  if(cnt>0) meanns/=(Float_t)cnt;
+  return meanns;
+}
+//______________________________________________________________________
+void AliITSOnlineSDDBase::WriteToFXS(){
+  //
+  Char_t outfilnam[100];
+  sprintf(outfilnam,"SDDbase_step1_mod%03d_sid%d.data",fModuleId,fSide);
+  FILE* outf=fopen(outfilnam,"w");
+  Float_t corrnoise=2.;
+  for(Int_t ian=0;ian<fgkNAnodes;ian++){
+    fprintf(outf,"%d %d %8.3f %8.3f %8.3f %8.3f\n",ian,IsAnodeGood(ian),GetAnodeBaseline(ian),GetAnodeRawNoise(ian),GetAnodeCommonMode(ian),corrnoise);
+  }
+  fclose(outf);  
+}
diff --git a/ITS/AliITSOnlineSDDBase.h b/ITS/AliITSOnlineSDDBase.h
new file mode 100644 (file)
index 0000000..c1cafb5
--- /dev/null
@@ -0,0 +1,66 @@
+#ifndef ALIITSONLINESDDBASE_H
+#define ALIITSONLINESDDBASE_H
+
+///////////////////////////////////////////////////////////////////
+//                                                               //
+// Class used for SDD baseline and noise analysis                //
+// Origin: F.Prino, Torino, prino@to.infn.it                     //
+//                                                               //
+///////////////////////////////////////////////////////////////////
+#include"AliITSOnlineSDD.h"
+#include<TMath.h>
+
+class TH2F;
+class TGraph;
+class AliITSOnlineSDDBase : public AliITSOnlineSDD {
+
+ public:
+  AliITSOnlineSDDBase();
+  AliITSOnlineSDDBase(Int_t mod, Int_t sid);
+  virtual ~AliITSOnlineSDDBase();
+  void Reset();
+  void AddEvent(TH2F* hrawd);
+  void ValidateAnodes();
+
+  void SetMinBaseline(Float_t bas=10.){fMinBaseline=bas;}
+  void SetMaxBaseline(Float_t bas=150.){fMaxBaseline=bas;}
+  void SetMinRawNoise(Float_t ns=0.001){fMinRawNoise=ns;}
+  void SetMaxRawNoise(Float_t ns=9.){fMaxRawNoise=ns;}
+  void SetNSigmaNoise(Float_t ns=4.){fNSigmaNoise=ns;}
+
+  Bool_t IsAnodeGood(Int_t iAnode)const{ return fGoodAnode[iAnode];}
+  Float_t GetAnodeBaseline(Int_t iAnode) const{
+    if(fNEvents>0) return fSumBaseline[iAnode]/fNEvents;
+    else return 0;
+  }
+  Float_t GetAnodeRawNoise(Int_t iAnode) const{
+    if(fNEvents>0) return TMath::Sqrt(fSumRawNoise[iAnode]/fNEvents-TMath::Power(GetAnodeBaseline(iAnode),2));
+    
+    else return 0;
+  }
+
+  Float_t CalcMeanRawNoise();
+  Float_t GetAnodeCommonMode(Int_t iAnode) const{
+    if(fNEvents>0) return fSumCMN[iAnode]/fNEvents;
+    else return 0;
+  }
+  Int_t GetNEvents() const {return fNEvents;}
+  void WriteToFXS();
+
+ protected:
+
+ private:
+  Int_t fNEvents;
+  Bool_t fGoodAnode[fgkNAnodes];
+  Float_t fSumBaseline[fgkNAnodes];
+  Float_t fSumRawNoise[fgkNAnodes];
+  Float_t fSumCMN[fgkNAnodes];
+  Float_t fMinBaseline;
+  Float_t fMaxBaseline;
+  Float_t fMinRawNoise;
+  Float_t fMaxRawNoise;
+  Float_t fNSigmaNoise;
+
+  ClassDef(AliITSOnlineSDDBase,1);
+};
+#endif
diff --git a/ITS/AliITSOnlineSDDInjectors.cxx b/ITS/AliITSOnlineSDDInjectors.cxx
new file mode 100644 (file)
index 0000000..47387b4
--- /dev/null
@@ -0,0 +1,391 @@
+/**************************************************************************
+ * Copyright(c) 2007-2009, 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 "AliITSOnlineSDDInjectors.h"
+#include <TH2F.h>
+#include <TGraphErrors.h>
+#include <TMath.h>
+
+///////////////////////////////////////////////////////////////////
+//                                                               //
+// Implementation of the class used for SDD injector analysis    //
+// Origin: F.Prino, Torino, prino@to.infn.it                     //
+//                                                               //
+///////////////////////////////////////////////////////////////////
+
+ClassImp(AliITSOnlineSDDInjectors)
+
+const Float_t AliITSOnlineSDDInjectors::fgkSaturation = 1008.;
+const Float_t AliITSOnlineSDDInjectors::fgkJitterTB = 8.;
+
+//______________________________________________________________________
+  AliITSOnlineSDDInjectors::AliITSOnlineSDDInjectors():AliITSOnlineSDD(),fHisto(),fTbZero(0.),fParam(),fPolOrder(0),fMinDriftVel(0.),fMaxDriftVel(0.),fThreshold(0.)
+{
+  // default constructor
+  SetMinDriftVel();
+  SetMaxDriftVel();
+  SetRangeLine1();
+  SetRangeLine2();
+  SetRangeLine3();
+  SetPositions();
+  SetPolOrder();
+  SetThreshold();
+}
+//______________________________________________________________________
+AliITSOnlineSDDInjectors::AliITSOnlineSDDInjectors(Int_t mod, Int_t sid):AliITSOnlineSDD(mod,sid),fHisto(),fTbZero(0.),fParam(),fPolOrder(0),fMinDriftVel(0.),fMaxDriftVel(0.),fThreshold(0.)
+{ // standard constructor
+  SetMinDriftVel();
+  SetMaxDriftVel();
+  SetRangeLine1();
+  SetRangeLine2();
+  SetRangeLine3();
+  SetPositions();
+  SetPolOrder();
+  SetThreshold();
+}
+//______________________________________________________________________
+AliITSOnlineSDDInjectors::~AliITSOnlineSDDInjectors(){
+  // Destructor
+  if(fHisto) delete fHisto;  
+  if(fParam) delete [] fParam;
+}
+//______________________________________________________________________
+void AliITSOnlineSDDInjectors::SetPositions(){
+  // 
+  Float_t kLinFromCenterUm[3]={31860.,17460.,660.};
+  Float_t kAnodeFromCenterUm=35085;
+  for(Int_t i=0;i<3;i++){
+    fPosition[i]=kAnodeFromCenterUm-kLinFromCenterUm[i];
+    fPosition[i]/=10000.; // from microns to cm
+  }
+}
+//______________________________________________________________________
+void AliITSOnlineSDDInjectors::Reset(){
+  for(Int_t i=0;i<kNInjectors;i++){ 
+    fDriftVel[i]=0.;
+    fSigmaDriftVel[i]=0.;
+  }
+  for(Int_t i=0;i<kNInjectors;i++){
+    for(Int_t j=0;j<3;j++){
+      fGoodInj[i][j]=0;
+      fCentroid[i][j]=0.;
+      fRMSCentroid[i][j]=0.;
+    }
+  }
+}
+//______________________________________________________________________
+void AliITSOnlineSDDInjectors::AnalyzeEvent(TH2F* his){
+  //
+  Reset();
+  fHisto=his;
+  FindGoodInjectors();
+  FindCentroids();
+  CalcTimeBinZero();
+  for(Int_t j=0;j<kNInjectors;j++) CalcDriftVelocity(j);
+  FitDriftVelocityVsAnode();
+}
+//______________________________________________________________________
+TGraphErrors* AliITSOnlineSDDInjectors::GetLineGraph(Int_t jlin){
+  // 
+  Float_t x[4],y[4],ex[4],ey[4];
+  x[0]=0.;
+  ex[0]=0.;
+  y[0]=fTbZero;
+  ey[0]=0.;
+  for(Int_t i=0;i<3;i++){
+    x[i+1]=fPosition[i];
+    ex[i+1]=0.;
+    y[i+1]=fCentroid[jlin][i];
+    ey[i+1]=fRMSCentroid[jlin][i];
+  }
+  TGraphErrors *g=new TGraphErrors(4,x,y,ex,ey);
+  return g;
+}
+//______________________________________________________________________
+Float_t AliITSOnlineSDDInjectors::GetDriftCoordinate(Float_t cAnode, Float_t cTimeBin){
+  Float_t vel=0;
+  for(Int_t i=0;i<=fPolOrder;i++) vel+=fParam[i]*TMath::Power(cAnode,(Float_t)i);
+  return vel*(cTimeBin-(fTbZero-fgkJitterTB))*25/1000.; 
+}
+//______________________________________________________________________
+TGraphErrors* AliITSOnlineSDDInjectors::GetDriftVelocityGraph(){
+  // 
+
+  Int_t ipt=0;
+  TGraphErrors *g=new TGraphErrors(0);
+  for(Int_t i=0;i<kNInjectors;i++){
+    if(fDriftVel[i]>0){ 
+      g->SetPoint(ipt,GetAnodeNumber(i),fDriftVel[i]);
+      g->SetPointError(ipt,0,fSigmaDriftVel[i]);
+      ipt++;
+    }
+  }
+  return g;
+}
+//______________________________________________________________________
+void AliITSOnlineSDDInjectors::CalcTimeBinZero(){
+  Float_t tzero=0.,intCont=0.;
+  for(Int_t ian=0;ian<fgkNAnodes;ian++){
+    for(Int_t itb=1;itb<fTbMin[0];itb++){
+      Float_t cont=fHisto->GetBinContent(itb,ian+1);
+      if(cont>fThreshold){
+       tzero+=cont*float(itb);
+       intCont+=cont;
+      }
+    }
+  }
+  if(intCont>0) fTbZero=tzero/intCont;
+}
+//______________________________________________________________________
+void AliITSOnlineSDDInjectors::FitDriftVelocityVsAnode(){
+  const Int_t nn=fPolOrder+1;
+  Float_t **mat = new Float_t*[nn];
+  for(Int_t i=0; i < nn; i++) mat[i] = new Float_t[nn];
+  Float_t *vect = new Float_t[nn];
+  for(Int_t k1=0;k1<nn;k1++){
+    vect[k1]=0;
+    for(Int_t k2=0;k2<nn;k2++){
+      mat[k1][k2]=0;
+      for(Int_t n=0; n<kNInjectors;n++){
+       Float_t x=(Float_t)GetAnodeNumber(n);
+       if(fDriftVel[n]>0) mat[k1][k2]+=TMath::Power(x,k1+k2)/TMath::Power(fSigmaDriftVel[n],2);
+      }
+    }
+  }
+  for(Int_t k1=0;k1<nn;k1++){
+    for(Int_t n=0; n<kNInjectors;n++){
+      Float_t x=(Float_t)GetAnodeNumber(n);
+      if(fDriftVel[n]>0) vect[k1]+=fDriftVel[n]*TMath::Power(x,k1)/TMath::Power(fSigmaDriftVel[n],2);
+    }
+  }
+  Int_t *iPivot = new Int_t[nn];
+  Int_t *indxR = new Int_t[nn];
+  Int_t *indxC = new Int_t[nn];
+  for(Int_t i=0;i<nn;i++) iPivot[i]=0;
+  Int_t iCol=-1,iRow=-1;
+  for(Int_t i=0;i<nn;i++){
+    Float_t big=0.;
+    for(Int_t j=0;j<nn;j++){
+      if(iPivot[j]!=1){
+       for(Int_t k=0;k<nn;k++){
+          if(iPivot[k]==0){
+            if(TMath::Abs(mat[j][k])>=big){
+              big=TMath::Abs(mat[j][k]);
+              iRow=j;
+              iCol=k;
+            }
+         }
+       }
+      }
+    }
+    iPivot[iCol]++;
+    Float_t aux;
+    if(iRow!=iCol){
+      for(Int_t l=0;l<nn;l++){
+       aux=mat[iRow][l];
+       mat[iRow][l]=mat[iCol][l];
+       mat[iCol][l]=aux;
+      }
+      aux=vect[iRow];
+      vect[iRow]=vect[iCol];
+      vect[iCol]=aux;
+    }
+    indxR[i]=iRow;
+    indxC[i]=iCol;
+    if(mat[iCol][iCol]==0) break;
+    Float_t pivinv=1./mat[iCol][iCol];
+    mat[iCol][iCol]=1;
+    for(Int_t l=0;l<nn;l++) mat[iCol][l]*=pivinv;
+    vect[iCol]*=pivinv;
+    for(Int_t m=0;m<nn;m++){
+      if(m!=iCol){
+       aux=mat[m][iCol];
+       mat[m][iCol]=0;
+       for(Int_t n=0;n<nn;n++) mat[m][n]-=mat[iCol][n]*aux;
+       vect[m]-=vect[iCol]*aux;
+      }
+    }    
+  }
+  delete [] iPivot;
+  delete [] indxR;
+  delete [] indxC;
+
+  if(fParam) delete [] fParam;
+  fParam=new Float_t[nn];
+  for(Int_t i=0; i<nn;i++)fParam[i]=vect[i];
+
+  for(Int_t i=0; i < nn; i++) delete [] mat[i];
+  delete [] mat;
+  delete [] vect;
+}
+//______________________________________________________________________
+void AliITSOnlineSDDInjectors::CalcDriftVelocity(Int_t jlin){
+  // 
+  Float_t sumY=0,sumX=0,sumXX=0,sumYY=0.,sumXY=0,sumWEI=0.;
+  Int_t npt=0;
+  Float_t y[3],ey[3];
+  Float_t tzero=0,erry=0;
+  for(Int_t i=0;i<3;i++){ 
+    y[i]=fCentroid[jlin][i];
+    ey[i]=fRMSCentroid[jlin][i];
+  }
+  for(Int_t i=0;i<3;i++){
+    if(fGoodInj[jlin][i]){
+      sumY+=y[i]/ey[i]/ey[i];
+      sumX+=fPosition[i]/ey[i]/ey[i];
+      sumXX+=fPosition[i]*fPosition[i]/ey[i]/ey[i];
+      sumYY+=y[i]*y[i]/ey[i]/ey[i];
+      sumXY+=fPosition[i]*y[i]/ey[i]/ey[i];
+      sumWEI+=1./ey[i]/ey[i];
+      tzero=fTbZero/ey[i]/ey[i];
+      erry=ey[i];
+      npt++;
+    }
+  }
+  Float_t vel=0,evel=0;
+  if(npt>1){ 
+    Float_t slope=(sumWEI*sumXY-sumY*sumX)/(sumWEI*sumXX-sumX*sumX);
+    Float_t eslope=TMath::Sqrt(sumWEI/(sumWEI*sumXX-sumX*sumX));
+    vel=1./slope*10000./25.;// micron/ns
+    evel=eslope/slope/slope*10000./25.;// micron/ns
+  }
+  if(npt==1){
+    Float_t slope=(sumY-tzero)/sumX;
+    Float_t eslope=erry/sumX;
+    vel=1./slope*10000./25.;// micron/ns    
+    evel=eslope/slope/slope*10000./25.;// micron/ns
+  }
+  if(vel>fMaxDriftVel||vel<fMinDriftVel){ 
+    vel=0.;
+    evel=0.;
+  }
+  fDriftVel[jlin]=vel;
+  fSigmaDriftVel[jlin]=evel;
+}
+//______________________________________________________________________
+Int_t AliITSOnlineSDDInjectors::GetAnodeNumber(Int_t iInjLine){
+  Int_t ian=-1;
+  if(iInjLine>32) return ian;
+  if(!fSide){
+    ian=iInjLine*8;
+    if(iInjLine==32) ian--;
+  }else{
+    ian=iInjLine*8-1;
+    if(iInjLine==0) ian=0;
+  }
+  return ian;
+}
+
+//______________________________________________________________________
+void AliITSOnlineSDDInjectors::FindGoodInjectors(){
+  // 
+  for(Int_t iii=0;iii<kNInjectors;iii++){
+    Int_t ian=GetAnodeNumber(iii);
+    for(Int_t ninj=0;ninj<3;ninj++){
+      for(Int_t jjj=fTbMin[ninj];jjj<fTbMax[ninj];jjj++){
+       Float_t c1=fHisto->GetBinContent(jjj,ian+1);
+       Float_t c2=fHisto->GetBinContent(jjj+1,ian+1);
+       Float_t c3=fHisto->GetBinContent(jjj+2,ian+1);
+       if(c1>fThreshold && c2>fThreshold && c3>fThreshold){ 
+         fGoodInj[iii][ninj]=1;
+         break;
+       }
+      }
+      //      for(Int_t jjj=fTbMin[ninj];jjj<fTbMax[ninj];jjj++){
+      //       Float_t c1=fHisto->GetBinContent(jjj,ian+1);
+      //       if(c1>=fgkSaturation){
+      //         fGoodInj[iii][ninj]=0;
+      //         break;
+      //       }
+      //      }
+    }
+  }
+}
+//______________________________________________________________________
+void AliITSOnlineSDDInjectors::FindCentroids(){
+  // 
+  for(Int_t iii=0;iii<kNInjectors;iii++){
+    Int_t ian=GetAnodeNumber(iii);
+    for(Int_t ninj=0;ninj<3;ninj++){
+      if(!fGoodInj[iii][ninj]) continue;
+      Float_t maxcont=0;
+      Int_t ilmax=-1;
+      for(Int_t jjj=fTbMin[ninj];jjj<fTbMax[ninj];jjj++){
+       Float_t cont=fHisto->GetBinContent(jjj,ian+1);
+       if(cont>maxcont){
+         maxcont=cont;
+         ilmax=jjj;
+       }
+      }
+      Float_t intCont=0;
+      Int_t jjj=ilmax;
+      while(1){
+       Float_t cont=fHisto->GetBinContent(jjj,ian+1);
+       if(cont<fThreshold) break;
+       if(cont<fgkSaturation){
+         fCentroid[iii][ninj]+=cont*(Float_t)jjj;
+         fRMSCentroid[iii][ninj]+=cont*TMath::Power((Float_t)jjj,2);
+         intCont+=cont;
+       }
+       jjj--;
+      }
+      jjj=ilmax+1;
+      while(1){
+       Float_t cont=fHisto->GetBinContent(jjj,ian+1);
+       if(cont<fThreshold) break;
+       if(cont<fgkSaturation){
+         fCentroid[iii][ninj]+=cont*float(jjj);
+         fRMSCentroid[iii][ninj]+=cont*TMath::Power((Float_t)jjj,2);
+         intCont+=cont;
+       }
+       jjj++;
+      }
+      if(intCont>0){ 
+       fCentroid[iii][ninj]/=intCont;
+       fRMSCentroid[iii][ninj]=TMath::Sqrt(fRMSCentroid[iii][ninj]/intCont-fCentroid[iii][ninj]*fCentroid[iii][ninj]);
+      }
+      else{ 
+       fCentroid[iii][ninj]=0.;
+       fRMSCentroid[iii][ninj]=0.;
+       fGoodInj[iii][ninj]=0;
+      }
+    }
+  }
+}
+//______________________________________________________________________
+void AliITSOnlineSDDInjectors::PrintInjMap(){
+  //
+  for(Int_t iii=0;iii<kNInjectors;iii++){
+    printf("Line%d-Anode%d: %d %d %d\n",iii,GetAnodeNumber(iii),fGoodInj[iii][0],fGoodInj[iii][1],fGoodInj[iii][2]);
+  }
+}
+//______________________________________________________________________
+void AliITSOnlineSDDInjectors::PrintCentroids(){
+  //
+  for(Int_t iii=0;iii<kNInjectors;iii++){
+    printf("Line%d-Anode%d: %f+-%f %f+-%f %f+-%f\n",iii,GetAnodeNumber(iii),fCentroid[iii][0],fRMSCentroid[iii][0],fCentroid[iii][1],fRMSCentroid[iii][1],fCentroid[iii][2],fRMSCentroid[iii][2]);
+  }
+}
+//______________________________________________________________________
+void AliITSOnlineSDDInjectors::WriteToFXS(){
+  //
+  Char_t outfilnam[100];
+  sprintf(outfilnam,"SDDinj_mod%03d_sid%d.data",fModuleId,fSide);
+  FILE* outf=fopen(outfilnam,"w");
+  for(Int_t ic=0;ic<fPolOrder+1;ic++){
+    fprintf(outf,"%G ",fParam[ic]);
+  }
+  fprintf(outf,"\n");
+  fclose(outf);  
+}
diff --git a/ITS/AliITSOnlineSDDInjectors.h b/ITS/AliITSOnlineSDDInjectors.h
new file mode 100644 (file)
index 0000000..2604286
--- /dev/null
@@ -0,0 +1,89 @@
+#ifndef ALIITSONLINESDDINJECTORS_H
+#define ALIITSONLINESDDINJECTORS_H
+
+
+///////////////////////////////////////////////////////////////////
+//                                                               //
+// Class used for SDD injector analysis                           //
+// Origin: F.Prino, Torino, prino@to.infn.it                     //
+//                                                               //
+///////////////////////////////////////////////////////////////////
+#include"AliITSOnlineSDD.h"
+
+class TH2F;
+class TGraphErrors;
+class AliITSOnlineSDDInjectors : public AliITSOnlineSDD {
+
+ public:
+  AliITSOnlineSDDInjectors();      
+  AliITSOnlineSDDInjectors(Int_t mod, Int_t sid);
+  virtual ~AliITSOnlineSDDInjectors();
+
+  void SetSide(Int_t sid){fSide=sid;}
+  void SetThreshold(Float_t thr=75.){fThreshold=thr;}
+  void SetRangeLine1(Int_t tbmin=40, Int_t tbmax=90){
+    fTbMin[0]=tbmin; fTbMax[0]=tbmax; 
+  }
+  void SetRangeLine2(Int_t tbmin=90, Int_t tbmax=140){
+    fTbMin[1]=tbmin; fTbMax[1]=tbmax; 
+  }
+  void SetRangeLine3(Int_t tbmin=170, Int_t tbmax=220){
+    fTbMin[2]=tbmin; fTbMax[2]=tbmax; 
+  }
+  void SetPolOrder(Int_t n=3){fPolOrder=n;}
+  void SetMinDriftVel(Float_t vmin=4.){fMinDriftVel=vmin;}
+  void SetMaxDriftVel(Float_t vmax=9.){fMaxDriftVel=vmax;}
+
+  TGraphErrors* GetLineGraph(Int_t jlin);
+  TGraphErrors* GetDriftVelocityGraph();
+  Float_t* GetDriftVelFitParam()const{ return fParam;}
+  Float_t GetDriftVelocity(Int_t jlin) const{return fDriftVel[jlin];}
+  Float_t GetSigmaDriftVelocity(Int_t jlin) const{return fSigmaDriftVel[jlin];}
+  Float_t GetTimeBinZero() const{return fTbZero;}
+  Float_t GetDriftCoordinate(Float_t cAnode, Float_t cTimeBin);
+  Int_t GetAnodeNumber(Int_t iInjLine);
+
+  void PrintInjMap();
+  void PrintCentroids();
+  void WriteToFXS();
+
+  void Reset();      
+  void AnalyzeEvent(TH2F* his);      
+  void FindGoodInjectors();
+  void FindCentroids();
+  void CalcDriftVelocity(Int_t jlin);
+  void CalcTimeBinZero();
+  void FitDriftVelocityVsAnode();
+
+ protected:
+  void SetPositions();
+ private:
+
+  enum {
+    kNInjectors = 33
+  };
+
+  AliITSOnlineSDDInjectors(const AliITSOnlineSDDInjectors& source);
+  AliITSOnlineSDDInjectors& operator = (const AliITSOnlineSDDInjectors& source);
+  static const Float_t fgkSaturation;
+  static const Float_t fgkJitterTB;
+
+  TH2F* fHisto;
+  Float_t fTbZero;
+  Float_t fPosition[3];
+  UShort_t fTbMin[3];
+  UShort_t fTbMax[3];
+  Bool_t fGoodInj[kNInjectors][3];
+  Float_t fCentroid[kNInjectors][3];
+  Float_t fRMSCentroid[kNInjectors][3];
+  Float_t fDriftVel[kNInjectors];
+  Float_t fSigmaDriftVel[kNInjectors];
+  Float_t *fParam;
+  Int_t fPolOrder;
+  Float_t fMinDriftVel;
+  Float_t fMaxDriftVel;
+  Float_t fThreshold;
+
+  ClassDef(AliITSOnlineSDDInjectors,1)
+};
+#endif
diff --git a/ITS/AliITSOnlineSDDTP.cxx b/ITS/AliITSOnlineSDDTP.cxx
new file mode 100644 (file)
index 0000000..4b66f52
--- /dev/null
@@ -0,0 +1,164 @@
+/**************************************************************************
+ * Copyright(c) 2007-2009, 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 "AliITSOnlineSDDTP.h"
+#include <TH2F.h>
+#include <TMath.h>
+
+
+///////////////////////////////////////////////////////////////////
+//                                                               //
+// Implemetation of the class SDD Test Pulse analysis            //
+// Origin: F.Prino, Torino, prino@to.infn.it                     //
+//                                                               //
+///////////////////////////////////////////////////////////////////
+
+
+ClassImp(AliITSOnlineSDDTP)
+
+//______________________________________________________________________
+AliITSOnlineSDDTP::AliITSOnlineSDDTP():AliITSOnlineSDD(),fNEvents(0),fDAQ(0.),fNSigmaGain(0.)
+{
+  // default constructor
+  Reset();
+  SetNSigmaGain();
+}
+//______________________________________________________________________
+AliITSOnlineSDDTP::AliITSOnlineSDDTP(Int_t mod, Int_t sid, Float_t xDAQ):AliITSOnlineSDD(mod,sid),fNEvents(0),fDAQ(xDAQ),fNSigmaGain(0.)
+{
+  // standard constructor
+  Reset();
+  SetNSigmaGain();
+}
+//______________________________________________________________________
+AliITSOnlineSDDTP::~AliITSOnlineSDDTP(){
+  // Destructor
+}
+//______________________________________________________________________
+void AliITSOnlineSDDTP::Reset(){
+  fNEvents=0;
+  for(Int_t i=0;i<fgkNAnodes;i++){
+    fGoodAnode[i]=1;
+    fBaseline[i]=0.;
+    fSumTPPeak[i]=0.;
+    fTPPos[i]=0.;
+  }
+  ReadBaselines();
+}
+
+//______________________________________________________________________
+void AliITSOnlineSDDTP::AddEvent(TH2F* hrawd){
+  // 
+  fNEvents++;
+  Double_t tbmax=(Double_t)hrawd->GetNbinsX();
+  for(Int_t ian=0;ian<fgkNAnodes;ian++){
+    Float_t auxmax=0.;
+    Int_t auxtb=0;
+    if(!fGoodAnode[ian]) continue;
+    for(Int_t itb=0;itb<tbmax;itb++){
+      Float_t cnt=hrawd->GetBinContent(itb+1,ian+1);
+      if(cnt>auxmax){ 
+       auxmax=cnt;
+       auxtb=itb;
+      }
+    }
+    fSumTPPeak[ian]+=auxmax-fBaseline[ian];
+    fTPPos[ian]+=auxtb;
+  }
+}
+//______________________________________________________________________
+void AliITSOnlineSDDTP::ReadBaselines(){
+  // assume baselines and good anodes are taken from previous run
+  Char_t basfilnam[100];
+  sprintf(basfilnam,"SDDbase_step1_mod%03d_sid%d.data",fModuleId,fSide);
+  FILE* basf=fopen(basfilnam,"r");
+  if(basf==0){
+    printf("Baselinefile not present, Set all baselines to 50\n");
+    for(Int_t ian=0;ian<fgkNAnodes;ian++){ 
+      fBaseline[ian]=50.;
+      fGoodAnode[ian]=1;
+    }
+    return;
+  }
+  Int_t n,ok;
+  Float_t base,rms,cmn,corrnoi;
+  for(Int_t ian=0;ian<fgkNAnodes;ian++){
+    fscanf(basf,"%d %d %f %f %f %f\n",&n,&ok,&base,&rms,&cmn,&corrnoi);
+    fBaseline[ian]=base;
+    fGoodAnode[ian]=ok;
+  }
+  fclose(basf);
+}
+
+//______________________________________________________________________
+void AliITSOnlineSDDTP::ValidateAnodes(){
+  Float_t meang,rmsg;
+  StatGain(meang,rmsg);
+  printf("<gain>=%f,rms=%f\n",meang,rmsg);
+  Float_t lowlim=meang-fNSigmaGain*rmsg;
+  Float_t hilim=meang+fNSigmaGain*rmsg;
+
+  for(Int_t ian=0;ian<fgkNAnodes;ian++){
+    if(!fGoodAnode[ian]) continue;
+    if(GetChannelGain(ian)<lowlim||GetChannelGain(ian)>hilim) fGoodAnode[ian]=0;
+  }
+}
+
+
+//______________________________________________________________________
+void AliITSOnlineSDDTP::StatGain(Float_t &mean, Float_t  &rms){
+  Float_t sum=0.,sumq=0.;
+  Int_t cnt=0;
+  for(Int_t ian=0;ian<fgkNAnodes;ian++){
+    if(!fGoodAnode[ian]) continue;
+    sum+=GetChannelGain(ian);
+    sumq+=TMath::Power(GetChannelGain(ian),2);
+    cnt++;
+  }
+  if(cnt>0){ 
+    mean=sum/(Float_t)cnt;
+    rms=TMath::Sqrt(sumq/(Float_t)cnt-mean*mean);
+  }else{ 
+    mean=0.;
+    rms=0.;
+  }
+  return;
+}
+
+//______________________________________________________________________
+void AliITSOnlineSDDTP::WriteToFXS(){
+  //
+  Char_t basfilnam[100];
+  sprintf(basfilnam,"SDDbase_step1_mod%03d_sid%d.data",fModuleId,fSide);
+  FILE* basf=fopen(basfilnam,"r");
+  Int_t n,ok;
+  Float_t base,rms,cmn,corrnoi;
+  Float_t noise[fgkNAnodes],cmncoef[fgkNAnodes],corrnoise[fgkNAnodes];
+  for(Int_t ian=0;ian<fgkNAnodes;ian++){
+    fscanf(basf,"%d %d %f %f %f %f\n",&n,&ok,&base,&rms,&cmn,&corrnoi);
+    noise[ian]=rms;
+    cmncoef[ian]=cmn;
+    corrnoise[ian]=corrnoi;
+  }
+  fclose(basf);
+  printf("Read All******************\n");
+  Char_t outfilnam[100];
+  sprintf(outfilnam,"SDDbase_mod%03d_sid%d.data",fModuleId,fSide);
+  FILE* outf=fopen(outfilnam,"w");
+  for(Int_t ian=0;ian<fgkNAnodes;ian++){
+    fprintf(outf,"%d %d %8.3f %8.3f %8.3f %8.3f %8.3f\n",ian,IsAnodeGood(ian),fBaseline[ian], noise[ian],cmncoef[ian],corrnoise[ian],GetChannelGain(ian));
+  }
+  fclose(outf);  
+}
+
diff --git a/ITS/AliITSOnlineSDDTP.h b/ITS/AliITSOnlineSDDTP.h
new file mode 100644 (file)
index 0000000..b1052d7
--- /dev/null
@@ -0,0 +1,50 @@
+#ifndef ALIITSONLINESDDTP_H
+#define ALIITSONLINESDDTP_H
+
+
+///////////////////////////////////////////////////////////////////
+//                                                               //
+// Class used for SDD Test Pulse analysis                        //
+// Origin: F.Prino, Torino, prino@to.infn.it                     //
+//                                                               //
+///////////////////////////////////////////////////////////////////
+
+#include"AliITSOnlineSDD.h"
+
+class TH2F;
+class TGraph;
+class AliITSOnlineSDDTP : public AliITSOnlineSDD {
+
+ public:
+  AliITSOnlineSDDTP();
+  AliITSOnlineSDDTP(Int_t mod, Int_t sid,Float_t xDAQ);
+  virtual ~AliITSOnlineSDDTP();
+  void Reset();
+  void AddEvent(TH2F* hrawd);
+  void ValidateAnodes();
+  void ReadBaselines();
+
+  void SetNSigmaGain(Float_t sig=3.){fNSigmaGain=sig;}
+  Bool_t IsAnodeGood(Int_t iAnode)const{ return fGoodAnode[iAnode];}
+  Int_t GetNEvents() const {return fNEvents;}
+  Float_t GetChannelGain(Int_t iAnode)const{
+    if(fNEvents>0) return fSumTPPeak[iAnode]/fNEvents/fDAQ;
+    else return 0;
+  }
+  void StatGain(Float_t &mean, Float_t  &rms);
+  void WriteToFXS();
+
+ protected:
+
+ private:
+  Int_t fNEvents;
+  Float_t fDAQ;
+  Bool_t fGoodAnode[fgkNAnodes];
+  Float_t fBaseline[fgkNAnodes];
+  Float_t fSumTPPeak[fgkNAnodes];
+  Float_t fTPPos[fgkNAnodes];
+  Float_t fNSigmaGain;
+
+  ClassDef(AliITSOnlineSDDTP,1);
+};
+#endif
index 092e4406c268b42988db662ada5604e776092c97..30b8f1996406ae3c96edad933d8dc221e83c5aef 100644 (file)
 #pragma link C++ class AliITSOnlineSPDscanMultiple+;
 #pragma link C++ class AliITSOnlineSPDscanSingle+;
 #pragma link C++ class AliITSPreprocessorSDD+;
+#pragma link C++ class AliITSOnlineSDD+;
+#pragma link C++ class AliITSOnlineSDDBase+;
+#pragma link C++ class AliITSOnlineSDDTP+;
+#pragma link C++ class AliITSOnlineSDDBTP+;
+#pragma link C++ class AliITSOnlineSDDInjectors+;
 
 #endif
index 168b7f2e9fc500ed0f4b36f05219310fcebdd8b6..1f949d67fdf203aa7054018d37ad1d9c664f394c 100644 (file)
@@ -68,7 +68,12 @@ SRCS =       AliITSDetTypeRec.cxx \
                AliITSBadChannelsAuxSPD.cxx \
                AliITSBadChannelsSPD.cxx  \
                AliITSChannelSPD.cxx  \
-               AliITSPreprocessorSDD.cxx
+               AliITSPreprocessorSDD.cxx \
+               AliITSOnlineSDD.cxx \
+               AliITSOnlineSDDBase.cxx \
+               AliITSOnlineSDDTP.cxx \
+               AliITSOnlineSDDBTP.cxx \
+               AliITSOnlineSDDInjectors.cxx
 
 HDRS:=  $(SRCS:.cxx=.h)