]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Bug fixes as suggested by Matus and updates for pass0 making also the postcorrection...
authorcblume <cblume@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 23 Jun 2010 14:32:09 +0000 (14:32 +0000)
committercblume <cblume@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 23 Jun 2010 14:32:09 +0000 (14:32 +0000)
TRD/AliTRDCalibTask.cxx
TRD/AliTRDCalibTask.h
TRD/AliTRDCalibraFillHisto.cxx
TRD/AliTRDCalibraFillHisto.h
TRD/AliTRDCalibraFit.cxx
TRD/AliTRDCalibraFit.h
TRD/AliTRDPreprocessorOffline.cxx
TRD/AliTRDPreprocessorOffline.h

index dafebb0c33c7a09ff0bfe70e9363fd21aa84e8ed..6d238c586160b77e63dd5098bb0749056dbd2809 100644 (file)
@@ -126,6 +126,10 @@ ClassImp(AliTRDCalibTask)
       fNbMaxCluster(2),
       fOfflineTracks(kFALSE),
       fStandaloneTracks(kFALSE),
+      fVersionGainUsed(0),
+      fSubVersionGainUsed(0),
+      fVersionVdriftUsed(0), 
+      fSubVersionVdriftUsed(0),
       fCalDetGain(0x0),
       fMaxEvent(0),
       fCounter(0),
@@ -252,6 +256,10 @@ void AliTRDCalibTask::UserCreateOutputObjects()
   fTRDCalibraFillHisto->SetPRF2dOn(); // choose to look at the PRF
   fTRDCalibraFillHisto->SetLinearFitterOn(fVdriftLinear); // Other possibility vdrift VDRIFT
   fTRDCalibraFillHisto->SetLinearFitterDebugOn(fVdriftLinear); // Other possibility vdrift
+  fTRDCalibraFillHisto->SetVersionGainUsed(fVersionGainUsed); // Gain Used
+  fTRDCalibraFillHisto->SetSubVersionGainUsed(fSubVersionGainUsed); // Gain Used
+  fTRDCalibraFillHisto->SetVersionVdriftUsed(fVersionVdriftUsed); // Vdrift Used
+  fTRDCalibraFillHisto->SetSubVersionVdriftUsed(fSubVersionVdriftUsed); // Vdrift Used
   for(Int_t k = 0; k < 3; k++){
     if(((fNz[k] != 10) && (fNrphi[k] != 10)) && ((fNz[k] != 100) && (fNrphi[k] != 100))) {
       fTRDCalibraFillHisto->SetNz(k,fNz[k]);                                    // Mode calibration
@@ -526,6 +534,7 @@ void AliTRDCalibTask::UserExec(Option_t *)
     return;
   }
   
+  
   fESDfriend = dynamic_cast<AliESDfriend*> (fESD->FindListObject("AliESDfriend"));
   if(!fESDfriend){
     AliError("fESDfriend not available");
@@ -612,13 +621,17 @@ void AliTRDCalibTask::UserExec(Option_t *)
     //ULong_t status = fkEsdTrack->GetStatus();
     
     fFriendTrack = fESDfriend->GetTrack(itrk);
-    if(!fFriendTrack)  continue;
+    if(!fFriendTrack)  {
+      //printf("No friend track %d\n",itrk);
+      continue;
+    }
     //////////////////////////////////////
     // Loop on calibration objects
     //////////////////////////////////////
     Int_t icalib=0;
     Int_t nTRDtrackV1=0;
     while((fCalibObject = (TObject *)(fFriendTrack->GetCalibObject(icalib++)))){
+      //printf("Name %s\n",fCalibObject->IsA()->GetName());
       if(strcmp(fCalibObject->IsA()->GetName(), "AliTRDtrackV1") != 0) continue;
       //printf("Find the calibration object\n");
       ++nTRDtrackV1;
index 39a9f42c5ee0c104a7c4df81c987216c0ba513b7..423403ab448670e6e9d1c167b9e3f6ddf3aef82c 100644 (file)
@@ -64,6 +64,10 @@ class AliTRDCalibTask : public AliAnalysisTaskSE {
   void SetUseSPDVertex()                                               {fVtxTPC=kFALSE; fVtxSPD=kTRUE ;} 
   void SetMinNbOfContributors(Int_t minNbOfContributors)               {fMinNbContributors = minNbOfContributors;};  
   void SetRangePrimaryVertexZ(Double_t rangePrimaryVertexZ)            {fRangePrimaryVertexZ = TMath::Abs(rangePrimaryVertexZ);};  
+  void SetVersionGainUsed(Int_t versionGainUsed)                       { fVersionGainUsed = versionGainUsed;   }
+  void SetSubVersionGainUsed(Int_t subVersionGainUsed)                 { fSubVersionGainUsed = subVersionGainUsed;   }
+  void SetVersionVdriftUsed(Int_t versionVdriftUsed)                   { fVersionVdriftUsed = versionVdriftUsed;   }
+  void SetSubVersionVdriftUsed(Int_t subVersionVdriftUsed)             { fSubVersionVdriftUsed = subVersionVdriftUsed;   }
   
   void SetLow(Int_t low)                                               {fLow=low;};
   void SetHigh(Int_t high)                                             {fHigh=high;};
@@ -73,8 +77,8 @@ class AliTRDCalibTask : public AliAnalysisTaskSE {
   void SetNbMaxCluster(Short_t nbMaxCluster)                           {fNbMaxCluster =  nbMaxCluster; }; 
   void SetOfflineTracks()                                              {fOfflineTracks=kTRUE; fStandaloneTracks=kFALSE; };
   void SetStandaloneTracks()                                           {fStandaloneTracks=kTRUE; fOfflineTracks=kFALSE; };
-
-  void SetCalDetGain(AliTRDCalDet * const calDetGain)                         {fCalDetGain = calDetGain;};  
+  
+  void SetCalDetGain(AliTRDCalDet * const calDetGain)                  {fCalDetGain = calDetGain;};  
 
   void SetMaxEvent(Int_t nbevents)                                     { fMaxEvent = nbevents; };
   void SetDebug(Int_t debug)                                           { fDebug = debug; };
@@ -144,6 +148,11 @@ class AliTRDCalibTask : public AliAnalysisTaskSE {
   Bool_t      fOfflineTracks;                    // Only Offline refitted tracks
   Bool_t      fStandaloneTracks;                 // Take only standalone tracks
 
+  Int_t       fVersionGainUsed;                  // VersionGainUsed 
+  Int_t       fSubVersionGainUsed;               // SubVersionGainUsed
+  Int_t       fVersionVdriftUsed;                // VersionVdriftUsed 
+  Int_t       fSubVersionVdriftUsed;             // SubVersionVdriftUsed
+
   AliTRDCalDet *fCalDetGain;                     // Calib object gain
 
   Int_t       fMaxEvent;                         // max events
index 57a240a81e55394f8a87b931b572497498e7da24..be1b40249c72530ab810d1271859d3329051b5ae 100644 (file)
@@ -142,6 +142,10 @@ AliTRDCalibraFillHisto::AliTRDCalibraFillHisto()
   ,fNormalizeNbOfCluster(kFALSE)
   ,fMaxCluster(0)
   ,fNbMaxCluster(0)
+  ,fVersionGainUsed(0)
+  ,fSubVersionGainUsed(0)
+  ,fVersionVdriftUsed(0) 
+  ,fSubVersionVdriftUsed(0)
   ,fCalibraMode(new AliTRDCalibraMode())
   ,fDebugStreamer(0)
   ,fDebugLevel(0)
@@ -215,6 +219,10 @@ AliTRDCalibraFillHisto::AliTRDCalibraFillHisto(const AliTRDCalibraFillHisto &c)
   ,fNormalizeNbOfCluster(c.fNormalizeNbOfCluster)
   ,fMaxCluster(c.fMaxCluster)
   ,fNbMaxCluster(c.fNbMaxCluster)
+  ,fVersionGainUsed(c.fVersionGainUsed)
+  ,fSubVersionGainUsed(c.fSubVersionGainUsed)
+  ,fVersionVdriftUsed(c.fVersionVdriftUsed) 
+  ,fSubVersionVdriftUsed(c.fSubVersionVdriftUsed)
   ,fCalibraMode(0x0)
   ,fDebugStreamer(0)
   ,fDebugLevel(c.fDebugLevel)
@@ -3201,7 +3209,11 @@ void AliTRDCalibraFillHisto::CreatePH2d(Int_t nn)
   // Create the 2D histos
   //
 
-  TString name("Nz");
+  TString name("Ver");
+  name += fVersionVdriftUsed;
+  name += "Subver";
+  name += fSubVersionVdriftUsed;
+  name += "Nz";
   name += fCalibraMode->GetNz(1);
   name += "Nrphi";
   name += fCalibraMode->GetNrphi(1);
@@ -3222,11 +3234,15 @@ void AliTRDCalibraFillHisto::CreateCH2d(Int_t nn)
   // Create the 2D histos
   //
 
-  TString name("Nz");
+  TString name("Ver");
+  name += fVersionGainUsed;
+  name += "Subver";
+  name += fSubVersionGainUsed;
+  name += "Nz";
   name += fCalibraMode->GetNz(0);
   name += "Nrphi";
   name += fCalibraMode->GetNrphi(0);
-
+  
   fCH2d = new TH2I("CH2d",(const Char_t *) name
                   ,fNumberBinCharge,0,300,nn,0,nn);
   fCH2d->SetYTitle("Det/pad groups");
index 1d1cd27effe9c3e4f72c4bc3157e67cfc9fe2645..e221853a1417e1a4bd6398a1a46a56e018e2bd37 100644 (file)
@@ -104,7 +104,11 @@ class AliTRDCalibraFillHisto : public TObject {
           void     SetVector2d(Bool_t vector2d = kTRUE)                      { fVector2d        = vector2d;          }
          void     SetLinearFitterOn(Bool_t linearfitteron = kTRUE)          { fLinearFitterOn      = linearfitteron;}
          void     SetLinearFitterDebugOn(Bool_t debug = kTRUE)              { fLinearFitterDebugOn = debug;         }
-                 
+         void     SetVersionGainUsed(Int_t versionGainUsed)                 { fVersionGainUsed = versionGainUsed;   }
+         void     SetSubVersionGainUsed(Int_t subVersionGainUsed)           { fSubVersionGainUsed = subVersionGainUsed;   }
+         void     SetVersionVdriftUsed(Int_t versionVdriftUsed)             { fVersionVdriftUsed = versionVdriftUsed;   }
+         void     SetSubVersionVdriftUsed(Int_t subVersionVdriftUsed)       { fSubVersionVdriftUsed = subVersionVdriftUsed;   }
+       
   
          Bool_t   GetPH2dOn() const                                         { return fPH2dOn;                 }
           Bool_t   GetCH2dOn() const                                         { return fCH2dOn;                 }
@@ -112,7 +116,11 @@ class AliTRDCalibraFillHisto : public TObject {
           Bool_t   GetHisto2d() const                                        { return fHisto2d;                }
           Bool_t   GetVector2d() const                                       { return fVector2d;               }
           Bool_t   GetLinearFitterOn() const                                 { return fLinearFitterOn;         }
-         Bool_t   GetLinearFitterDebugOn() const                            { return fLinearFitterDebugOn; }
+         Bool_t   GetLinearFitterDebugOn() const                            { return fLinearFitterDebugOn;    }
+         Int_t    GetVersionGainUsed() const                                { return fVersionGainUsed;        }
+         Int_t    GetSubVersionGainUsed() const                             { return fSubVersionGainUsed;     }
+         Int_t    GetVersionVdriftUsed() const                              { return fVersionVdriftUsed;      }
+         Int_t    GetSubVersionVdriftUsed() const                           { return fSubVersionVdriftUsed;   }
 
 
   // Get stuff that are filled
@@ -197,6 +205,11 @@ AliTRDCalibraVector *GetCalibraVector() const                                { r
          Bool_t   fNormalizeNbOfCluster;   // Normalize with the number of cluster for the gain
          Float_t  fMaxCluster;             // Max amplitude of one cluster
          Short_t  fNbMaxCluster;           // Number of tb at the end
+  // Back correction
+         Int_t    fVersionGainUsed;        // VersionGainUsed 
+         Int_t    fSubVersionGainUsed;     // SubVersionGainUsed
+         Int_t    fVersionVdriftUsed;      // VersionVdriftUsed 
+         Int_t    fSubVersionVdriftUsed;   // SubVersionVdriftUsed
   // Calibration mode
          AliTRDCalibraMode *fCalibraMode;  // Calibration mode
 
index 9debbfd709284a0bfe2e4c40ea5115a93d296d87..a6e67efbbd01242f4f5711c6d2b557da17bee29b 100644 (file)
@@ -63,6 +63,7 @@
 #include <TLinearFitter.h>
 #include <TVectorD.h>
 #include <TROOT.h>
+#include <TString.h>
 
 #include "AliLog.h"
 #include "AliMathBase.h"
@@ -361,7 +362,8 @@ Bool_t AliTRDCalibraFit::AnalyseCH(const TH2I *ch)
   //
 
   // Set the calibration mode
-  const char *name = ch->GetTitle();
+  //const char *name = ch->GetTitle();
+  TString name = ch->GetTitle();
   if(!SetModeCalibration(name,0)) return kFALSE;
 
   // Number of Ybins (detectors or groups of pads)
@@ -459,7 +461,8 @@ Bool_t AliTRDCalibraFit::AnalyseCH(AliTRDCalibraVector *calvect)
   //
 
   // Set the calibraMode
-  const char *name = calvect->GetNameCH();
+  //const char *name = calvect->GetNameCH();
+  TString name = calvect->GetNameCH();
   if(!SetModeCalibration(name,0)) return kFALSE;  
 
   // Number of Xbins (detectors or groups of pads)
@@ -552,7 +555,8 @@ Bool_t AliTRDCalibraFit::AnalysePH(const TProfile2D *ph)
   //
 
   // Set the calibration mode
-  const char *name = ph->GetTitle();
+  //const char *name = ph->GetTitle();
+  TString name = ph->GetTitle();
   if(!SetModeCalibration(name,1)) return kFALSE;
   
   //printf("Mode calibration set\n");
@@ -653,7 +657,8 @@ Bool_t AliTRDCalibraFit::AnalysePH(AliTRDCalibraVector *calvect)
   //
 
   // Set the calibration mode
-  const char *name = calvect->GetNamePH();
+  //const char *name = calvect->GetNamePH();
+  TString name = calvect->GetNamePH();
   if(!SetModeCalibration(name,1)) return kFALSE;
 
   // Number of Xbins (detectors or groups of pads)
@@ -732,7 +737,8 @@ Bool_t AliTRDCalibraFit::AnalysePRF(const TProfile2D *prf)
   //
 
   // Set the calibration mode
-  const char *name = prf->GetTitle();
+  //const char *name = prf->GetTitle();
+  TString name = prf->GetTitle();
   if(!SetModeCalibration(name,2)) return kFALSE;
 
   // Number of Ybins (detectors or groups of pads)
@@ -819,7 +825,8 @@ Bool_t AliTRDCalibraFit::AnalysePRFMarianFit(const TProfile2D *prf)
   //
 
   // Set the calibration mode
-  const char *name = prf->GetTitle();
+  //const char *name = prf->GetTitle();
+  TString name = prf->GetTitle();
   if(!SetModeCalibration(name,2)) return kFALSE;
 
   // Number of Ybins (detectors or groups of pads)
@@ -914,7 +921,8 @@ Bool_t AliTRDCalibraFit::AnalysePRF(AliTRDCalibraVector *calvect)
   //
 
   // Set the calibra mode
-  const char *name = calvect->GetNamePRF();
+  //const char *name = calvect->GetNamePRF();
+  TString name = calvect->GetNamePRF();
   if(!SetModeCalibration(name,2)) return kFALSE;
   //printf("test0 %s\n",name);
 
@@ -990,7 +998,8 @@ Bool_t AliTRDCalibraFit::AnalysePRFMarianFit(AliTRDCalibraVector *calvect)
   //
 
   // Set the calibra mode
-  const char *name = calvect->GetNamePRF();
+  //const char *name = calvect->GetNamePRF();
+  TString name = calvect->GetNamePRF();
   if(!SetModeCalibration(name,2)) return kFALSE;
   //printf("test0 %s\n",name);
   Int_t    nbg     = GetNumberOfGroupsPRF((const char *)name);
@@ -1166,7 +1175,7 @@ Bool_t AliTRDCalibraFit::AnalyseLinearFitters(AliTRDCalibraVdriftLinearFit *cali
 }
 //____________Functions for seeing if the pad is really okey___________________
 //_____________________________________________________________________________
-Int_t AliTRDCalibraFit::GetNumberOfGroupsPRF(const char* nametitle)
+Int_t AliTRDCalibraFit::GetNumberOfGroupsPRF(TString nametitle)
 {
   //
   // Get numberofgroupsprf
@@ -1182,25 +1191,25 @@ Int_t AliTRDCalibraFit::GetNumberOfGroupsPRF(const char* nametitle)
   const Char_t *pattern6 = "Ngp6";
 
   // Nrphi mode
-  if (strstr(nametitle,pattern0)) {
+  if (strstr(nametitle.Data(),pattern0)) {
     return 0;
   }
-  if (strstr(nametitle,pattern1)) {
+  if (strstr(nametitle.Data(),pattern1)) {
     return 1;
   }
-  if (strstr(nametitle,pattern2)) {
+  if (strstr(nametitle.Data(),pattern2)) {
     return 2;
   }
-  if (strstr(nametitle,pattern3)) {
+  if (strstr(nametitle.Data(),pattern3)) {
     return 3;
   }
-  if (strstr(nametitle,pattern4)) {
+  if (strstr(nametitle.Data(),pattern4)) {
     return 4;
   }
-  if (strstr(nametitle,pattern5)) {
+  if (strstr(nametitle.Data(),pattern5)) {
     return 5;
   }
-  if (strstr(nametitle,pattern6)){
+  if (strstr(nametitle.Data(),pattern6)){
     return 6;
   }
   else return -1;
@@ -1208,7 +1217,7 @@ Int_t AliTRDCalibraFit::GetNumberOfGroupsPRF(const char* nametitle)
 
 }
 //_____________________________________________________________________________
-Bool_t AliTRDCalibraFit::SetModeCalibration(const char *name, Int_t i)
+Bool_t AliTRDCalibraFit::SetModeCalibration(TString name, Int_t i)
 {
   //
   // Set fNz[i] and fNrphi[i] of the AliTRDCalibraFit::Instance()
@@ -1222,7 +1231,7 @@ Bool_t AliTRDCalibraFit::SetModeCalibration(const char *name, Int_t i)
 
 }
 //_____________________________________________________________________________
-Bool_t AliTRDCalibraFit::SetNrphiFromTObject(const char *name, Int_t i)
+Bool_t AliTRDCalibraFit::SetNrphiFromTObject(TString name, Int_t i)
 {
   //
   // Set fNrphi[i] of the AliTRDCalibraFit::Instance()
@@ -1245,7 +1254,7 @@ Bool_t AliTRDCalibraFit::SetNrphiFromTObject(const char *name, Int_t i)
   const Char_t *patternz100 = "Nz100";
 
   // Nrphi mode
-  if ((strstr(name,patternrphi100)) && (strstr(name,patternz100))) {
+  if ((strstr(name.Data(),patternrphi100)) && (strstr(name.Data(),patternz100))) {
     fCalibraMode->SetAllTogether(i);
     fNbDet = 540;
     if (fDebugLevel > 1) {
@@ -1253,7 +1262,7 @@ Bool_t AliTRDCalibraFit::SetNrphiFromTObject(const char *name, Int_t i)
     }
     return kTRUE;
   }
-  if ((strstr(name,patternrphi10)) && (strstr(name,patternz10))) {
+  if ((strstr(name.Data(),patternrphi10)) && (strstr(name.Data(),patternz10))) {
     fCalibraMode->SetPerSuperModule(i);
     fNbDet = 30;
     if (fDebugLevel > 1) {
@@ -1262,49 +1271,49 @@ Bool_t AliTRDCalibraFit::SetNrphiFromTObject(const char *name, Int_t i)
     return kTRUE;
   }
   
-  if (strstr(name,patternrphi0)) {
+  if (strstr(name.Data(),patternrphi0)) {
     fCalibraMode->SetNrphi(i ,0);
     if (fDebugLevel > 1) {
       AliInfo(Form("fNbDet %d and 0",fNbDet));
     }
     return kTRUE;
   }
-  if (strstr(name,patternrphi1)) {
+  if (strstr(name.Data(),patternrphi1)) {
     fCalibraMode->SetNrphi(i, 1);
     if (fDebugLevel > 1) {
       AliInfo(Form("fNbDet %d and 1",fNbDet));
     }
     return kTRUE;
   }
-  if (strstr(name,patternrphi2)) {
+  if (strstr(name.Data(),patternrphi2)) {
     fCalibraMode->SetNrphi(i, 2);
     if (fDebugLevel > 1) {
       AliInfo(Form("fNbDet %d and 2",fNbDet));
     }    
     return kTRUE;
   }
-  if (strstr(name,patternrphi3)) {
+  if (strstr(name.Data(),patternrphi3)) {
     fCalibraMode->SetNrphi(i, 3);
     if (fDebugLevel > 1) {
       AliInfo(Form("fNbDet %d and 3",fNbDet));
     }   
     return kTRUE;
   }
-  if (strstr(name,patternrphi4)) {
+  if (strstr(name.Data(),patternrphi4)) {
     fCalibraMode->SetNrphi(i, 4);
     if (fDebugLevel > 1) {
       AliInfo(Form("fNbDet %d and 4",fNbDet));
     }   
     return kTRUE;
   }
-  if (strstr(name,patternrphi5)) {
+  if (strstr(name.Data(),patternrphi5)) {
     fCalibraMode->SetNrphi(i, 5);
     if (fDebugLevel > 1) {
       AliInfo(Form("fNbDet %d and 5",fNbDet));
     }
     return kTRUE;
   }
-  if (strstr(name,patternrphi6)) {
+  if (strstr(name.Data(),patternrphi6)) {
     fCalibraMode->SetNrphi(i, 6);
     if (fDebugLevel > 1) {
       AliInfo(Form("fNbDet %d and 6",fNbDet));
@@ -1320,7 +1329,7 @@ Bool_t AliTRDCalibraFit::SetNrphiFromTObject(const char *name, Int_t i)
   
 }
 //_____________________________________________________________________________
-Bool_t AliTRDCalibraFit::SetNzFromTObject(const char *name, Int_t i)
+Bool_t AliTRDCalibraFit::SetNzFromTObject(TString name, Int_t i)
 {
   //
   // Set fNz[i] of the AliTRDCalibraFit::Instance()
@@ -1339,7 +1348,7 @@ Bool_t AliTRDCalibraFit::SetNzFromTObject(const char *name, Int_t i)
   const Char_t *patternz10 = "Nz10";
   const Char_t *patternz100 = "Nz100";
 
-  if ((strstr(name,patternrphi100)) && (strstr(name,patternz100))) {
+  if ((strstr(name.Data(),patternrphi100)) && (strstr(name.Data(),patternz100))) {
     fCalibraMode->SetAllTogether(i);
     fNbDet = 540;
     if (fDebugLevel > 1) {
@@ -1347,7 +1356,7 @@ Bool_t AliTRDCalibraFit::SetNzFromTObject(const char *name, Int_t i)
     }
     return kTRUE;
   }
-  if ((strstr(name,patternrphi10)) && (strstr(name,patternz10))) {
+  if ((strstr(name.Data(),patternrphi10)) && (strstr(name.Data(),patternz10))) {
     fCalibraMode->SetPerSuperModule(i);
     fNbDet = 30;
     if (fDebugLevel > 1) {
@@ -1355,35 +1364,35 @@ Bool_t AliTRDCalibraFit::SetNzFromTObject(const char *name, Int_t i)
     }
     return kTRUE;
   }
-  if (strstr(name,patternz0)) {
+  if (strstr(name.Data(),patternz0)) {
     fCalibraMode->SetNz(i, 0);
     if (fDebugLevel > 1) {
       AliInfo(Form("fNbDet %d and 0",fNbDet));
     }
     return kTRUE;
   }
-  if (strstr(name,patternz1)) {
+  if (strstr(name.Data(),patternz1)) {
     fCalibraMode->SetNz(i ,1);
     if (fDebugLevel > 1) {
       AliInfo(Form("fNbDet %d and 1",fNbDet));
     }
     return kTRUE;
   }
-  if (strstr(name,patternz2)) {
+  if (strstr(name.Data(),patternz2)) {
     fCalibraMode->SetNz(i ,2);
     if (fDebugLevel > 1) {    
       AliInfo(Form("fNbDet %d and 2",fNbDet));
     }
     return kTRUE;
   }
-  if (strstr(name,patternz3)) {
+  if (strstr(name.Data(),patternz3)) {
     fCalibraMode->SetNz(i ,3);
     if (fDebugLevel > 1) {
       AliInfo(Form("fNbDet %d and 3",fNbDet));
     }
     return kTRUE;  
   }
-  if (strstr(name,patternz4)) {
+  if (strstr(name.Data(),patternz4)) {
     fCalibraMode->SetNz(i ,4);
     if (fDebugLevel > 1) {    
       AliInfo(Form("fNbDet %d and 4",fNbDet));
@@ -6203,3 +6212,4 @@ Double_t AliTRDCalibraFit::GausConstant(const Double_t *x, const Double_t *par)
   return gauss;
 
 }
+
index 407ec7119e829beea3e8254ee94e04a830ee0c49..2c91c850eea1cef66c15b65e3fdc8b2e82357875 100644 (file)
@@ -18,6 +18,7 @@
 # include <TVectorD.h>
 #endif 
 
+class TString;
 class TTree;
 class TProfile2D;
 class TGraphErrors;
@@ -75,7 +76,7 @@ class AliTRDCalibraFit : public TObject {
   Bool_t   AnalyseLinearFitters(AliTRDCalibraVdriftLinearFit *calivdli);
   
   // Pad Calibration
-  Bool_t   SetModeCalibration(const char *name, Int_t i);
+  Bool_t   SetModeCalibration(TString name, Int_t i);
   
   //Reset Function
   void     ResetVectorFit();
@@ -293,9 +294,9 @@ class AliTRDCalibraFit : public TObject {
        void     ReconstructFitRowMinRowMax(Int_t idect, Int_t i);
        Bool_t   CheckFitVoir();
        void     NormierungCharge();
-       Bool_t   SetNrphiFromTObject(const char *name, Int_t i);
-       Bool_t   SetNzFromTObject(const char *name, Int_t i);
-       Int_t    GetNumberOfGroupsPRF(const char* nametitle);
+       Bool_t   SetNrphiFromTObject(TString name, Int_t i);
+       Bool_t   SetNzFromTObject(TString name, Int_t i);
+       Int_t    GetNumberOfGroupsPRF(TString nametitle);
        
        // Calculate the mean coefs from the database
        Bool_t   CalculVdriftCoefMean();
@@ -354,3 +355,4 @@ class AliTRDCalibraFit : public TObject {
 };
   
 #endif
+
index 2cd0c3f0a8c44e972d96c14c01e13652fa54fe01..c0ccdfc892dc4842df99ed30eba40ec535cfee93 100644 (file)
 
 
 /*
-  Responsible: marian.ivanov@cern.ch 
-  Code to analyze the TPC calibration and to produce OCDB entries  
+  Responsible: Raphaelle Bailhache (rbailhache@ikf.uni-frankfurt.de) 
+  Code to analyze the TRD calibration and to produce OCDB entries  
 
 
    .x ~/rootlogon.C
    gSystem->Load("libANALYSIS");
-   gSystem->Load("libTPCcalib");
+   gSystem->Load("libTRDcalib");
 
    AliTRDPreprocessorOffline proces;
    TString ocdbPath="local:////"
@@ -65,6 +65,7 @@ AliTRDPreprocessorOffline::AliTRDPreprocessorOffline():
   fMethodSecond(kTRUE),
   fNameList("TRDCalib"),
   fCalDetGainUsed(0x0),
+  fCalDetVdriftUsed(0x0),
   fCH2d(0x0),
   fPH2d(0x0),
   fPRF2d(0x0),
@@ -72,7 +73,11 @@ AliTRDPreprocessorOffline::AliTRDPreprocessorOffline():
   fNEvents(0x0),
   fAbsoluteGain(0x0),
   fPlots(new TObjArray(8)),
-  fCalibObjects(new TObjArray(8))
+  fCalibObjects(new TObjArray(8)),
+  fVersionGainUsed(0),
+  fSubVersionGainUsed(0),
+  fVersionVdriftUsed(0), 
+  fSubVersionVdriftUsed(0)
 {
   //
   // default constructor
@@ -85,6 +90,7 @@ AliTRDPreprocessorOffline::~AliTRDPreprocessorOffline() {
   //
 
   if(fCalDetGainUsed) delete fCalDetGainUsed;
+  if(fCalDetVdriftUsed) delete fCalDetVdriftUsed;
   if(fCH2d) delete fCH2d;
   if(fPH2d) delete fPH2d;
   if(fPRF2d) delete fPRF2d;
@@ -146,6 +152,7 @@ void AliTRDPreprocessorOffline::CalibGain(const Char_t* file, Int_t startRunNumb
   //
   AnalyzeGain();
   if(fCalDetGainUsed) CorrectFromDetGainUsed();
+  if(fCalDetVdriftUsed) CorrectFromDetVdriftUsed();
   //
   // 3. Append QA plots
   //
@@ -190,11 +197,41 @@ void AliTRDPreprocessorOffline::CalibPRF(const Char_t* file, Int_t startRunNumbe
   
 }
 
+Bool_t AliTRDPreprocessorOffline::Init(const Char_t* fileName){
+  //
+  // read the calibration used during the reconstruction
+  // 
+
+  if(ReadVdriftT0Global(fileName)) {
+    
+    TString nameph = fPH2d->GetTitle();
+    fVersionVdriftUsed = GetVersion(nameph);  
+    fSubVersionVdriftUsed = GetSubVersion(nameph);    
+
+    //printf("Found Version %d, Subversion %d for vdrift\n",fVersionVdriftUsed,fSubVersionVdriftUsed);
+  
+  }
+
+  if(ReadGainGlobal(fileName)) {
+
+    TString namech = fCH2d->GetTitle();
+    fVersionGainUsed = GetVersion(namech);  
+    fSubVersionGainUsed = GetSubVersion(namech);    
+
+    //printf("Found Version %d, Subversion %d for gain\n",fVersionGainUsed,fSubVersionGainUsed);
+
+  }
+   
+  return kTRUE;
+  
+}
+
 
 Bool_t AliTRDPreprocessorOffline::ReadGainGlobal(const Char_t* fileName){
   //
   // read calibration entries from file
   // 
+  if(fCH2d) return kTRUE;
   TFile fcalib(fileName);
   TObjArray * array = (TObjArray*)fcalib.Get(fNameList);
   if (array){
@@ -221,6 +258,7 @@ Bool_t AliTRDPreprocessorOffline::ReadVdriftT0Global(const Char_t* fileName){
   //
   // read calibration entries from file
   // 
+  if(fPH2d) return kTRUE;
   TFile fcalib(fileName);
   TObjArray * array = (TObjArray*)fcalib.Get(fNameList);
   if (array){
@@ -245,6 +283,7 @@ Bool_t AliTRDPreprocessorOffline::ReadVdriftLinearFitGlobal(const Char_t* fileNa
   //
   // read calibration entries from file
   // 
+  if(fAliTRDCalibraVdriftLinearFit) return kTRUE;
   TFile fcalib(fileName);
   TObjArray * array = (TObjArray*)fcalib.Get(fNameList);
   if (array){
@@ -266,6 +305,7 @@ Bool_t AliTRDPreprocessorOffline::ReadPRFGlobal(const Char_t* fileName){
   //
   // read calibration entries from file
   // 
+  if(fPRF2d) return kTRUE;
   TFile fcalib(fileName);
   TObjArray * array = (TObjArray*)fcalib.Get(fNameList);
   if (array){
@@ -304,13 +344,17 @@ Bool_t AliTRDPreprocessorOffline::AnalyzeGain(){
 
 
   Bool_t ok = kFALSE;
+  Bool_t meanother = kFALSE;
   // enough statistics
   if ((nbtg >                  0) && 
       (nbfit        >= 0.5*nbE) && (nbE > 30)) {
     // create the cal objects
-    calibra->PutMeanValueOtherVectorFit(1,kTRUE);
+    if(!fCalDetGainUsed) {
+      calibra->PutMeanValueOtherVectorFit(1,kTRUE);
+      meanother = kTRUE;
+    }
     TObjArray object           = calibra->GetVectorFit();
-    AliTRDCalDet *calDetGain   = calibra->CreateDetObjectGain(&object);
+    AliTRDCalDet *calDetGain   = calibra->CreateDetObjectGain(&object,meanother);
     TH1F *coefGain  = calDetGain->MakeHisto1DAsFunctionOfDet();
     // Put them in the array
     fCalibObjects->AddAt(calDetGain,kGain);
@@ -463,22 +507,72 @@ Bool_t AliTRDPreprocessorOffline::AnalyzePRF(){
 }
 
 void AliTRDPreprocessorOffline::CorrectFromDetGainUsed() {
-
+  //
+  // Correct from the gas gain used afterwards
+  //
   AliTRDCalDet *calDetGain = (AliTRDCalDet *) fCalibObjects->At(kGain);
   if(!calDetGain) return;
 
+  // Calculate mean
+  Double_t mean = 0.0;
+  Int_t nbdet = 0;
+  
   for(Int_t det = 0; det < 540; det++) {
     
     Float_t gaininit = fCalDetGainUsed->GetValue(det);
     Float_t gainout = calDetGain->GetValue(det);
-    
-    calDetGain->SetValue(det,gaininit*gainout);
 
+   
+    if(TMath::Abs(gainout-1.0) > 0.000001) {
+      mean += (gaininit*gainout);
+      nbdet++;
+    }  
+  }
+  if(nbdet > 0) mean = mean/nbdet;
+   
+  for(Int_t det = 0; det < 540; det++) {
+    
+    Float_t gaininit = fCalDetGainUsed->GetValue(det);
+    Float_t gainout = calDetGain->GetValue(det);
+  
+    if(TMath::Abs(gainout-1.0) > 0.000001) calDetGain->SetValue(det,gaininit*gainout);
+    else calDetGain->SetValue(det,mean);
   }
 
 
 }
 
+void AliTRDPreprocessorOffline::CorrectFromDetVdriftUsed() {
+  //
+  // Correct from the drift velocity
+  //
+
+  //printf("Correct for vdrift\n");
+
+  AliTRDCalDet *calDetGain = (AliTRDCalDet *) fCalibObjects->At(kGain);
+  if(!calDetGain) return;
+
+  Int_t detVdrift = kVdriftPHDet;
+  if(fMethodSecond) detVdrift = kVdriftLinear;
+  AliTRDCalDet *calDetVdrift = (AliTRDCalDet *) fCalibObjects->At(detVdrift);
+  if(!calDetVdrift) return;
+
+  // Calculate mean
+  
+  for(Int_t det = 0; det < 540; det++) {
+    
+    Float_t vdriftinit = fCalDetVdriftUsed->GetValue(det);
+    Float_t vdriftout = calDetVdrift->GetValue(det);
+
+    Float_t gain = calDetGain->GetValue(det);
+    if(vdriftout > 0.0) gain = gain*vdriftinit/vdriftout;
+    calDetGain->SetValue(det,gain);
+
+    
+  }
+}
+
 void AliTRDPreprocessorOffline::UpdateOCDBGain(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
   //
   // Update OCDB entry
@@ -493,10 +587,6 @@ void AliTRDPreprocessorOffline::UpdateOCDBGain(Int_t startRunNumber, Int_t endRu
   AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
   AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kGain);
   if(calDet) gStorage->Put(calDet, id1, metaData);
-  //else {
-  //  printf("No calDet object for Gain\n");
-  //}
-
     
 
 }
@@ -519,10 +609,7 @@ void AliTRDPreprocessorOffline::UpdateOCDBVdrift(Int_t startRunNumber, Int_t end
   AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
   AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(detVdrift);
   if(calDet) gStorage->Put(calDet, id1, metaData);
-  //else {
-  //  printf("No calDet object for Vdrift\n");
-  //}
-
+  
   //
 
   if(!fMethodSecond) {
@@ -535,11 +622,7 @@ void AliTRDPreprocessorOffline::UpdateOCDBVdrift(Int_t startRunNumber, Int_t end
     AliCDBId id1Pad("TRD/Calib/LocalVdrift", startRunNumber, endRunNumber);
     AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kVdriftPHPad);
     if(calPad) gStorage->Put(calPad, id1Pad, metaDataPad);
-    //else {
-    //  printf("No calPad object for Vdrift\n");
-    //}
-
-  
+     
   }
 
 }
@@ -558,11 +641,7 @@ void AliTRDPreprocessorOffline::UpdateOCDBT0(Int_t startRunNumber, Int_t endRunN
   AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
   AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kT0PHDet);
   if(calDet) gStorage->Put(calDet, id1, metaData);
-  //else {
-  //  printf("No calDet object for T0\n");
-  //}
  
-
   //
 
   AliCDBMetaData *metaDataPad= new AliCDBMetaData();
@@ -573,9 +652,7 @@ void AliTRDPreprocessorOffline::UpdateOCDBT0(Int_t startRunNumber, Int_t endRunN
   AliCDBId id1Pad("TRD/Calib/LocalT0", startRunNumber, endRunNumber);
   AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kT0PHPad);
   if(calPad) gStorage->Put(calPad, id1Pad, metaDataPad);
-  //else {
-  //  printf("No calPad object for T0\n");
-  //}
 
 
 }
@@ -594,10 +671,60 @@ void AliTRDPreprocessorOffline::UpdateOCDBPRF(Int_t startRunNumber, Int_t endRun
   AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
   AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kPRF);
   if(calPad) gStorage->Put(calPad, id1, metaData);
-  //else {
-  //  printf("No calPad object for PRF\n");
-  //}
  
 
 }
+//_____________________________________________________________________________
+Int_t AliTRDPreprocessorOffline::GetVersion(TString name) const
+{
+  //
+  // Get version from the title
+  //
+  
+  // Some patterns
+  const Char_t *version = "Ver";
+  if(!strstr(name.Data(),version)) return -1;
+  
+  for(Int_t ver = 0; ver < 999999999; ver++) {
+
+    TString vertry(version);
+    vertry += ver;
+    vertry += "Subver";
+
+    //printf("vertry %s and name %s\n",vertry.Data(),name.Data());
+
+    if(strstr(name.Data(),vertry.Data())) return ver;
+    
+  }
+  
+  return -1;
+
+}
+
+//_____________________________________________________________________________
+Int_t AliTRDPreprocessorOffline::GetSubVersion(TString name) const
+{
+  //
+  // Get subversion from the title
+  //
+  
+  // Some patterns
+  const Char_t *subversion = "Subver";
+  if(!strstr(name.Data(),subversion)) return -1;
+  
+  for(Int_t ver = 0; ver < 999999999; ver++) {
+    
+    TString vertry(subversion);
+    vertry += ver;
+    vertry += "Nz";
+
+    //printf("vertry %s and name %s\n",vertry.Data(),name.Data());
+
+    if(strstr(name.Data(),vertry.Data())) return ver;
+    
+  }
+  
+  return -1;
+
+}
 
index e15b66eebc0c7814b505076be0990b2cef775b14..e44f42b53850f4006a176028e719c8ca1d63ca46 100644 (file)
@@ -42,8 +42,12 @@ public:
   void SetNameList(TString nameList) { fNameList = nameList;};
   TString GetNameList() const { return fNameList;}; 
   void SetCalDetGain(AliTRDCalDet *calDetGainUsed) {fCalDetGainUsed = calDetGainUsed;};
+  void SetCalDetVdrift(AliTRDCalDet *calDetVdriftUsed) {fCalDetVdriftUsed = calDetVdriftUsed;};
   AliTRDCalDet *GetCalDetGain() const { return fCalDetGainUsed;};
+  AliTRDCalDet *GetCalDetVdrift() const { return fCalDetVdriftUsed;};
 
+  Bool_t Init(const Char_t* fileName);
+  
   void CalibVdriftT0(const Char_t* file, Int_t startRunNumber, Int_t endRunNumber, TString ocdbStorage="");
   void CalibGain(const Char_t* file, Int_t startRunNumber, Int_t endRunNumber,  TString  ocdbStorage="");
   void CalibPRF(const Char_t* file, Int_t startRunNumber, Int_t endRunNumber,  TString  ocdbStorage="");
@@ -59,17 +63,24 @@ public:
   Bool_t AnalyzePRF(); 
   
   void CorrectFromDetGainUsed();
+  void CorrectFromDetVdriftUsed();
 
   void UpdateOCDBT0(Int_t startRunNumber, Int_t endRunNumber, const char* storagePath);
   void UpdateOCDBVdrift(Int_t startRunNumber, Int_t endRunNumber, const char* storagePath);
   void UpdateOCDBGain(Int_t  startRunNumber, Int_t endRunNumber, const char* storagePath);
   void UpdateOCDBPRF(Int_t  startRunNumber, Int_t endRunNumber, const char* storagePath);
 
+  Int_t    GetVersionGainUsed() const                                { return fVersionGainUsed;        }
+  Int_t    GetSubVersionGainUsed() const                             { return fSubVersionGainUsed;     }
+  Int_t    GetVersionVdriftUsed() const                              { return fVersionVdriftUsed;      }
+  Int_t    GetSubVersionVdriftUsed() const                           { return fSubVersionVdriftUsed;   }
   
-private:
+  
+ private:
   Bool_t fMethodSecond;                   // Second Method for drift velocity   
   TString fNameList;                      // Name of the list
   AliTRDCalDet *fCalDetGainUsed;          // CalDet used and to be corrected for
+  AliTRDCalDet *fCalDetVdriftUsed;        // CalDet used and to be corrected for
   TH2I *fCH2d;                            // Gain
   TProfile2D *fPH2d;                      // Drift velocity first method
   TProfile2D *fPRF2d;                     // PRF
@@ -78,7 +89,13 @@ private:
   TH2F *fAbsoluteGain;                    // Absolute Gain calibration
   TObjArray * fPlots;                     // array with some plots to check
   TObjArray * fCalibObjects;              // array with calibration objects 
+  Int_t    fVersionGainUsed;              // VersionGainUsed 
+  Int_t    fSubVersionGainUsed;           // SubVersionGainUsed
+  Int_t    fVersionVdriftUsed;            // VersionVdriftUsed 
+  Int_t    fSubVersionVdriftUsed;         // SubVersionVdriftUsed
 
+  Int_t GetSubVersion(TString name) const;
+  Int_t GetVersion(TString name) const;