]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/macros/CalibrateTPC.C
Don't call stage function (Marian)
[u/mrichter/AliRoot.git] / TPC / macros / CalibrateTPC.C
1 /*
2   Simple calibration analysis
3   //
4   //0. Setup memory chcecker if you want 
5   //
6   gSystem->Load("$ROOTSYS/lib/libGui.so");
7   gSystem->Load("$ROOTSYS/lib/libTree.so");
8   gSystem->Load("$MEMSTAT/libMemStat.so");
9   TMemStat *memstat = new TMemStat(100000000,10000000,kTRUE);
10   AliSysInfo::AddCallBack(TMemStatManager::GetInstance()->fStampCallBack);
11
12   AliSysInfo::AddStamp("Start");  
13   //1. Load needed libraries
14   gSystem->Load("libANALYSIS");
15   gSystem->Load("libTPCcalib");
16   //
17   // Setup analysis manager
18   //
19   .L $ALICE_ROOT/TPC/macros/CalibrateTPC.C
20   AliAnalysisManager * mgr = SetupCalibTask();
21   //
22   // Process data - chain
23   //
24   gSystem->AddIncludePath("-I$ALICE_ROOT/TPC/macros");
25   gROOT->LoadMacro("$ALICE_ROOT/TPC/macros/AliXRDPROOFtoolkit.cxx+")
26   AliXRDPROOFtoolkit tool;
27   TChain * chain = tool.MakeChain("chain.txt","esdTree",0,1000);
28   chain->Lookup();
29   // memory
30   mgr->SetNSysInfo(100); 
31   //
32   mgr->SetDebugLevel(1);
33   mgr->StartAnalysis("local",chain);
34   // delete manager
35   //
36   delete mgr;
37   AliSysInfo::AddStamp("End");
38   //
39   // analyze memstat report
40   //
41   delete memstat;
42   TMemStat draw("memstat.root");
43   draw.MakeReport(0,0,"order 0 sortstat 3 sortstamp 0 sortdeep 10 stackdeep 15 maxlength 50")   
44 */
45
46
47 AliAnalysisManager * SetupCalibTask() {
48   //
49   //
50   //
51   TStopwatch stopwatch;
52   stopwatch.Start();
53   //
54   // set magnetic field form the cosmos - it should be provided by framework
55   AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., 2);
56   AliTracker::SetFieldMap(field,0);
57   //
58   AliAnalysisManager *mgr=new AliAnalysisManager("TestManager");
59
60   AliESDInputHandler* esdH=new AliESDInputHandler;
61   esdH->SetActiveBranches("ESDfriend");
62   mgr->SetInputEventHandler(esdH);  
63   //
64   //
65   AliCDBManager::Instance()->SetRun(1) ;
66   AliCDBManager::Instance()->SetDefaultStorage("local://$ALICE_ROOT");
67   AliTPCClusterParam * clusterParam = AliTPCcalibDB::Instance()->GetClusterParam();
68
69   AliTPCAnalysisTaskcalib *task1=new AliTPCAnalysisTaskcalib("TPC calibration task");
70   
71   AliTPCcalibTracksCuts *cuts = new AliTPCcalibTracksCuts(20, 0.4, 0.5, 0.13, 0.018);
72
73   //
74   AliTPCcalibTracks *calibTracks =  new AliTPCcalibTracks("calibTracks", "Resolution calibration object for tracks", clusterParam, cuts); 
75   AliTPCcalibTracksGain *calibTracksGain =  new AliTPCcalibTracksGain("calibTracksGain","Gain calibration using tracks",cuts); 
76   AliTPCcalibAlign *calibAlign = new AliTPCcalibAlign("alignTPC","Alignment of the TPC sectors");
77   calibTracks->SetDebugLevel(0);
78   calibTracks->SetStreamLevel(0);
79   calibTracksGain->SetDebugLevel(0);
80   calibTracksGain->SetStreamLevel(0);
81   calibAlign->SetDebugLevel(20);
82   calibAlign->SetStreamLevel(2);
83   //
84  // ---*---*-----*-*-----*----------*---
85   // ADD CALIB JOBS HERE!!!!!!!!!!!!!!!!
86   task1->AddJob(calibAlign);
87   //task1->AddJob(calibTracksGain);
88   //task1->AddJob(calibTracks);
89   //  task1->AddJob(new AliTPCcalibBase);
90   // task1->AddJob(new AliTPCcalibV0);
91   // -*----*----*---*-*------*-------**--
92   // -------*--*---------*-----*-------*-
93
94   mgr->AddTask(task1);
95
96   AliAnalysisDataContainer *cinput1
97     =mgr->CreateContainer("cchain1",TChain::Class(),
98                           AliAnalysisManager::kInputContainer);
99   AliAnalysisDataContainer *coutput1
100     =mgr->CreateContainer("asdofhaw",TObjArray::Class(),
101                           AliAnalysisManager::kOutputContainer,
102                           "CalibObjects.root");
103
104   mgr->ConnectInput(task1,0,cinput1);
105   mgr->ConnectOutput(task1,0,coutput1);
106
107   if (!mgr->InitAnalysis()) return;
108   mgr->PrintStatus(); 
109   
110   stopwatch.Stop();
111   stopwatch.Print();
112   return mgr;
113 }