]> git.uio.no Git - u/mrichter/AliRoot.git/blob - VZERO/UpdateVZEROTimeDelaysEntries.C
Fixing memory leaks
[u/mrichter/AliRoot.git] / VZERO / UpdateVZEROTimeDelaysEntries.C
1 void UpdateVZEROTimeDelaysEntries(Int_t year = 2010, Int_t period = 0, Bool_t local = kFALSE)
2 {
3   Int_t runRange[2];
4   if (year == 2009) {
5     runRange[0] = 0;
6     runRange[1] = 105268;
7   }
8   else if (year == 2010) {
9     if (period == 0) {
10       runRange[0] = 106031;
11       runRange[1] = 116353;
12     }
13     else if (period == 1) {
14       runRange[0] = 116354;
15       runRange[1] = 118285;
16     }
17     else if (period == 2) {
18       runRange[0] = 118486;
19       runRange[1] = 118556;
20     }
21     else {
22       printf("Invalid run period...\n");
23       return;
24     }
25   }
26   else {
27     printf("Invalid year...\n");
28     return;
29   }
30
31   AliCDBManager *man = AliCDBManager::Instance();
32   man->SetDefaultStorage(Form("alien://?folder=/alice/data/%d/OCDB",year));
33
34   man->SetRun(runRange[1]);
35
36   AliCDBEntry *entry = man->Get("VZERO/Calib/Data");
37   AliVZEROCalibData *calibda = (AliVZEROCalibData*)entry->GetObject();
38   printf("Year %d Period %d:\n",year,period);
39   printf("Calib:\n");
40   for(Int_t i = 0; i < 64; ++i) printf("%.2f ",calibda->GetTimeOffset(i));
41   printf("\n");
42   entry = man->Get("VZERO/Calib/TimeDelays");
43   TH1F *delays = (TH1F*)entry->GetObject();
44   TH1F *delaysNew = new TH1F(*delays);
45   printf("Delay:\n");
46   for(Int_t i = 0; i < 64; ++i) printf("%.2f ",delaysNew->GetBinContent(i+1));
47   printf("\n");
48   for(Int_t i = 0; i < 64; ++i) {
49     if (year == 2009) delaysNew->SetBinContent(i+1,delays->GetBinContent(i+1)+5.0);
50     if (year == 2010 && period == 0) delaysNew->SetBinContent(i+1,5.0);
51     if (year == 2010 && (period == 1 || period == 2)) {
52       Int_t board = i / 8;
53       Int_t channel = i % 8;
54       Int_t j = AliVZEROCalibData::GetOfflineChannelNumber(board,channel);
55       delaysNew->SetBinContent(j+1,
56                                delays->GetBinContent(j+1)+
57                                calibda->GetTimeOffset(i)-
58                                calibda->GetTimeOffset(j));
59     }
60   }
61   printf("CorrDelay:\n");
62   for(Int_t i = 0; i < 64; ++i) printf("%.2f ",delaysNew->GetBinContent(i+1));
63   printf("\n");
64
65   if (local) man->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
66
67   {
68     AliCDBMetaData *md= new AliCDBMetaData(); // metaData describing the object
69     md->SetResponsible("Brigitte Cheynis");
70     md->SetBeamPeriod(0);
71     md->SetAliRootVersion(gSystem->Getenv("ARVERSION"));
72     md->SetComment(Form("Time delays channel by channel (corrected values for year %d, period %d)",year,period));
73     md->PrintMetaData();
74
75     AliCDBStorage *storLoc = man->GetDefaultStorage();
76     AliCDBId id("VZERO/Calib/TimeDelays",runRange[0],runRange[1]);
77
78     storLoc->Put(delaysNew, id, md);
79
80     storLoc->Delete();
81     delete md;
82   }
83
84 }