]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Updates necessary because of the changes i the framework
authormarian <marian@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 23 Sep 2011 15:41:44 +0000 (15:41 +0000)
committermarian <marian@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 23 Sep 2011 15:41:44 +0000 (15:41 +0000)
(Mikolaj)

PWG1/AliAnalysisTaskITSTPCalignment.cxx
PWG1/AliRelAlignerKalmanArray.cxx
PWG1/AliRelAlignerKalmanArray.h

index 2c1094b178702191795adec7650155f0653ca6e4..fc2ea39d6dbda33efb916ccecdb252debc6837e8 100644 (file)
@@ -116,9 +116,9 @@ Bool_t AliAnalysisTaskITSTPCalignment::UserNotify()
 void AliAnalysisTaskITSTPCalignment::UserCreateOutputObjects()
 {
   // Create output objects
-  // Called once
+
   fArrayITSglobal = new AliRelAlignerKalmanArray();
-  fArrayITSglobal->SetName("outputArrayITSglobal");
+  fArrayITSglobal->SetName("array ITS global");
   fArrayITSglobal->SetupArray(fT0,fTend,fSlotWidth);
   fArrayITSglobal->SetOutRejSigmaOnMerge(fOutRejSigmaOnMerge);
   //set up the template
@@ -128,7 +128,7 @@ void AliAnalysisTaskITSTPCalignment::UserCreateOutputObjects()
   templ->SetRejectOutliersSigma2Median(fRejectOutliersSigma2Median);
   templ->SetOutRejSigma2Median(fOutRejSigma2Median);
   fArrayITSsa = new AliRelAlignerKalmanArray();
-  fArrayITSsa->SetName("outputArrayITSsa");
+  fArrayITSsa->SetName("array ITS SA");
   fArrayITSsa->SetupArray(fT0,fTend,fSlotWidth);
   fArrayITSsa->SetOutRejSigmaOnMerge(fOutRejSigmaOnMerge);
   //set up the template
@@ -139,6 +139,8 @@ void AliAnalysisTaskITSTPCalignment::UserCreateOutputObjects()
   templ->SetOutRejSigma2Median(fOutRejSigma2Median);
 
   fList = new TList();
+  fList->SetName("QA");
+  fList->SetOwner(kTRUE);
 
   TH2F* pZYAResidualsHistBpos = new TH2F("fZYAResidualsHistBpos","z-r\\phi residuals side A (z>0), Bpos", 100, -4., 4., 100, -1.5, 1.5 );
   pZYAResidualsHistBpos->SetXTitle("\\deltaz [cm]");
@@ -189,6 +191,8 @@ void AliAnalysisTaskITSTPCalignment::UserCreateOutputObjects()
   pLowPtZCResidualsHistBpos->SetXTitle("Pt");
   pLowPtZCResidualsHistBpos->SetYTitle("\\deltaz [cm]");
   TList* listBpos = new TList();
+  listBpos->SetName("B+");
+  listBpos->SetOwner(kTRUE);
   fList->Add(listBpos); //0
   listBpos->Add(pZYAResidualsHistBpos);
   listBpos->Add(pZYCResidualsHistBpos);
@@ -256,6 +260,8 @@ void AliAnalysisTaskITSTPCalignment::UserCreateOutputObjects()
   pLowPtZCResidualsHistBneg->SetXTitle("Pt");
   pLowPtZCResidualsHistBneg->SetYTitle("\\deltaz [cm]");
   TList* listBneg = new TList();
+  listBneg->SetName("B-");
+  listBneg->SetOwner(kTRUE);
   fList->Add(listBneg); //1
   listBneg->Add(pZYAResidualsHistBneg);
   listBneg->Add(pZYCResidualsHistBneg);
@@ -323,6 +329,8 @@ void AliAnalysisTaskITSTPCalignment::UserCreateOutputObjects()
   pLowPtZCResidualsHistBnil->SetXTitle("Pt");
   pLowPtZCResidualsHistBnil->SetYTitle("\\deltaz [cm]");
   TList* listBnil = new TList();
+  listBnil->SetName("B0");
+  listBnil->SetOwner(kTRUE);
   fList->Add(listBnil); //2
   listBnil->Add(pZYAResidualsHistBnil);
   listBnil->Add(pZYCResidualsHistBnil);
@@ -340,18 +348,19 @@ void AliAnalysisTaskITSTPCalignment::UserCreateOutputObjects()
   listBnil->Add(pLowPtYCResidualsHistBnil);
   listBnil->Add(pLowPtZAResidualsHistBnil);
   listBnil->Add(pLowPtZCResidualsHistBnil);
+
   TH1F* pNmatchingEff=new TH1F("pNmatchingEff","matching efficiency",50,0.,1.);
   fList->Add(pNmatchingEff); //3
 
-  TH1I* pChecks = new TH1I("pChecks","various checks",10,0,10);
+  TH1I* pChecks = new TH1I("some counters","some counters",10,0,10);
   pChecks->GetXaxis()->SetBinLabel(1+kNoESD,"no ESD");
   pChecks->GetXaxis()->SetBinLabel(1+kNoESDfriend,"no ESDfriend");
   pChecks->GetXaxis()->SetBinLabel(1+kNoFriendTrack,"no friendtrack");
-  pChecks->GetXaxis()->SetBinLabel(1+kNoITSoutParams,"no params");
-  pChecks->GetXaxis()->SetBinLabel(1+kESDfriend,"ESD friends");
+  pChecks->GetXaxis()->SetBinLabel(1+kNoITSoutParams,"no ITS out params");
+  pChecks->GetXaxis()->SetBinLabel(1+kESDfriend,"ESD friend");
   pChecks->GetXaxis()->SetBinLabel(1+kFriendsSkipBit,"friends+skipbit");
   pChecks->GetXaxis()->SetBinLabel(1+kFriendTrack,"friendtrack");
-  pChecks->GetXaxis()->SetBinLabel(1+kITSoutParams,"params");
+  pChecks->GetXaxis()->SetBinLabel(1+kITSoutParams,"ITS out params");
   fList->Add(pChecks); //4
 
   fAligner = new AliRelAlignerKalman();
@@ -363,6 +372,12 @@ void AliAnalysisTaskITSTPCalignment::UserCreateOutputObjects()
 
   fDebugTree = new TTree("debugTree","tree with debug info");
   fDebugTree->Branch("aligner","AliRelAlignerKalman",&fAligner);
+  
+  // Post output data.
+  PostData(0, fDebugTree);
+  PostData(1, fList);
+  PostData(2, fArrayITSglobal);
+  PostData(3, fArrayITSsa);
 }
 
 //________________________________________________________________________
@@ -393,16 +408,16 @@ void AliAnalysisTaskITSTPCalignment::UserExec(Option_t *)
     pChecks->Fill(kESDfriend);
   }
 
-  if (eventFriend->TestSkipBit()) pChecks->Fill(kFriendsSkipBit);
+  //event->SetESDfriend(eventFriend);
+  
+  if (eventFriend->TestSkipBit()) 
+  {
+    pChecks->Fill(kFriendsSkipBit);
+    return;
+  }
 
   //Update the parmeters
   AnalyzeESDevent(event);
-
-  // Post output data.
-  PostData(0, fDebugTree);
-  PostData(1, fList);
-  PostData(2, fArrayITSglobal);
-  PostData(3, fArrayITSsa);
 }
 
 //________________________________________________________________________
@@ -427,7 +442,7 @@ void AliAnalysisTaskITSTPCalignment::AnalyzeESDevent(AliESDEvent* event)
     TObjArray arrayParamTPC(ntracks);
 
     Int_t n = FindMatchingTracks(arrayParamITS, arrayParamTPC, event);
-    AliInfo(Form("n matched: %i\n",n));
+    AliInfo(Form("matched tracklet pairs: %i\n",n));
 
     TH1F* nMatchingEff=static_cast<TH1F*>(fList->At(3));
     Float_t ratio = (float)n/(float)ntracks;
@@ -454,7 +469,7 @@ void AliAnalysisTaskITSTPCalignment::AnalyzeESDevent(AliESDEvent* event)
         fAligner->SetTimeStamp(event->GetTimeStamp());
       }
 
-      //alignment check
+      //alignment
       if (alignerSA) alignerSA->AddTrackParams(paramsITS, paramsTPC);
     }
     arrayParamITS.Delete();
@@ -466,6 +481,8 @@ void AliAnalysisTaskITSTPCalignment::AnalyzeESDevent(AliESDEvent* event)
 void AliAnalysisTaskITSTPCalignment::Terminate(Option_t *)
 {
   // Called once at the end of the query
+  fArrayITSsa->Print("a");
+  fAligner->Print();
 }
 
 //________________________________________________________________________
index 7f933748028945ac6fda95542da926d0a57ab7de..096a6f5f569977ebb220380ff862f320efa6d6a4 100644 (file)
@@ -33,6 +33,8 @@
 #include <AliESDEvent.h>
 #include <AliRelAlignerKalman.h>
 #include "AliRelAlignerKalmanArray.h"
+#include "TList.h"
+#include "TBrowser.h"
 
 ClassImp(AliRelAlignerKalmanArray)
 
@@ -45,7 +47,8 @@ AliRelAlignerKalmanArray::AliRelAlignerKalmanArray():
     fOutRejSigmaOnMerge(10.),
     fOutRejSigmaOnSmooth(1.),
     fAlignerTemplate(),
-    fPArray(NULL)
+    fPArray(NULL),
+    fListOfGraphs(NULL)
 {
   //ctor
 }
@@ -59,7 +62,8 @@ AliRelAlignerKalmanArray::AliRelAlignerKalmanArray(Int_t t0, Int_t tend, Int_t s
     fOutRejSigmaOnMerge(10.),
     fOutRejSigmaOnSmooth(1.),
     fAlignerTemplate(),
-    fPArray(NULL)
+    fPArray(NULL),
+    fListOfGraphs(new TList)
 {
   //ctor
   if (slotwidth==0) fSize = 1;
@@ -67,6 +71,8 @@ AliRelAlignerKalmanArray::AliRelAlignerKalmanArray(Int_t t0, Int_t tend, Int_t s
   fPArray = new AliRelAlignerKalman* [fSize];
   if (fPArray) for (Int_t i=0;i<fSize;i++){fPArray[i]=NULL;}//fill with zeros
   else fSize=0;
+  fListOfGraphs->SetName("graphs");
+  fListOfGraphs->SetOwner(kTRUE);
 }
 
 //______________________________________________________________________________
@@ -78,7 +84,8 @@ AliRelAlignerKalmanArray::AliRelAlignerKalmanArray( const AliRelAlignerKalmanArr
     fOutRejSigmaOnMerge(in.fOutRejSigmaOnMerge),
     fOutRejSigmaOnSmooth(in.fOutRejSigmaOnSmooth),
     fAlignerTemplate(in.fAlignerTemplate),
-    fPArray(NULL)
+    fPArray(NULL),
+    fListOfGraphs(new TList)
 {
   //copy ctor
   fPArray = new AliRelAlignerKalman* [fSize];
@@ -88,6 +95,8 @@ AliRelAlignerKalmanArray::AliRelAlignerKalmanArray( const AliRelAlignerKalmanArr
     if (in.fPArray[i]) fPArray[i]=new AliRelAlignerKalman(*(in.fPArray[i]));
     else fPArray[i]=NULL;
   }
+  fListOfGraphs->SetName("graphs");
+  fListOfGraphs->SetOwner(kTRUE);
 }
 
 //______________________________________________________________________________
@@ -96,6 +105,7 @@ AliRelAlignerKalmanArray::~AliRelAlignerKalmanArray()
   //dtor
   ClearContents();
   delete [] fPArray;
+  delete fListOfGraphs;
 }
 
 //______________________________________________________________________________
@@ -153,6 +163,7 @@ Long64_t AliRelAlignerKalmanArray::Merge( TCollection* list )
         if (!a1 && a2) fPArray[i] = new AliRelAlignerKalman(*a2);
     }
   }
+  fListOfGraphs->Delete();
   return 0;
 }
 
@@ -316,10 +327,11 @@ TGraphErrors* AliRelAlignerKalmanArray::MakeGraph(Int_t iparam) const
     return NULL;
   }
 
-  TVectorD vx(GetEntries());
-  TVectorD vy(GetEntries());
-  TVectorD vex(GetEntries());
-  TVectorD vey(GetEntries());
+  Int_t n=GetEntries();
+  TVectorD vx(n);
+  TVectorD vy(n);
+  TVectorD vex(n);
+  TVectorD vey(n);
   Int_t entry=0;
   for (Int_t i=0; i<fSize; i++)
   {
@@ -374,6 +386,7 @@ TGraphErrors* AliRelAlignerKalmanArray::MakeGraph(Int_t iparam) const
     graphtitley="dv/dy [cm/\\micros/m]";
     break;
   }
+  graph->SetName(graphtitle);
   graph->SetTitle(graphtitle);
   TAxis* xas = graph->GetXaxis();
   TAxis* yas = graph->GetYaxis();
@@ -429,3 +442,22 @@ void AliRelAlignerKalmanArray::PropagateToTime(AliRelAlignerKalman* al, UInt_t t
   (*cov) += corr;
 }
 
+//______________________________________________________________________________
+void AliRelAlignerKalmanArray::Browse(TBrowser* b )
+{
+   if (!b) return;
+   if (!fListOfGraphs) 
+   {
+     fListOfGraphs=new TList();
+     fListOfGraphs->SetName("graphs");
+     fListOfGraphs->SetOwner(kTRUE);
+   }
+   for (Int_t i=0;i<9;i++)
+   {
+     TGraphErrors* gr = dynamic_cast<TGraphErrors*>(fListOfGraphs->At(i));
+     if (gr) continue;
+     fListOfGraphs->AddAt(MakeGraph(i),i);
+   }
+   b->Add(fListOfGraphs);
+}
+
index 5429f6c308e691c3335614dc60c2548b8ccab335..fd896593ab14819c516e8b2f9211e0d60a61ac81 100644 (file)
@@ -16,6 +16,8 @@ class TNamed;
 class TTree;
 class TCollection;
 class AliESDEvent;
+class TBrowser;
+class TList;
 #include <AliRelAlignerKalman.h>
 
 class AliRelAlignerKalmanArray:public TNamed
@@ -47,6 +49,7 @@ public:
   AliRelAlignerKalmanArray* MakeSmoothArray() const;
   void SetOutRejSigmaOnMerge(Double_t s) {fOutRejSigmaOnMerge=s;}
   void SetOutRejSigmaOnSmooth(Double_t s) {fOutRejSigmaOnSmooth=s;}
+  void Browse(TBrowser *b);
 
 private:
   void ClearContents();
@@ -59,8 +62,9 @@ private:
   Double_t fOutRejSigmaOnSmooth;          //how much outlier rejection on Smooth
   AliRelAlignerKalman fAlignerTemplate;  //template
   AliRelAlignerKalman** fPArray;         //[fSize] an array of aligners
+  TList* fListOfGraphs;                  //!hold the graphs
   
-  ClassDef(AliRelAlignerKalmanArray,4);   //AliRelAlignerKalman class
+  ClassDef(AliRelAlignerKalmanArray,5);   //AliRelAlignerKalman class
 };