]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDCreateDummyCDB.C
Added a commented out version with new digitizers.
[u/mrichter/AliRoot.git] / TRD / AliTRDCreateDummyCDB.C
1 #if !defined( __CINT__) || defined(__MAKECINT__)
2
3 #include <iostream>
4
5 #include <AliCDBManager.h>
6 #include <AliCDBStorage.h>
7 #include <AliCDBEntry.h>
8 #include <AliCDBMetaData.h>
9
10 #include "AliTRDgeometry.h"
11
12 #include "AliTRDCalROC.h"
13 #include "AliTRDCalChamber.h"
14 #include "AliTRDCalStack.h"
15 #include "AliTRDCalPad.h"
16 #include "AliTRDCalDet.h"
17 #include "AliTRDCalGlobals.h"
18
19 #endif
20
21 // run number for the dummy file
22 const Int_t gkDummyRun = 0;
23 AliCDBStorage* gStorLoc = 0;
24
25 TObject* CreatePadObject(const char* shortName, const char* description, Float_t value)
26 {
27   AliTRDCalPad *calPad = new AliTRDCalPad(shortName, description);
28   for (Int_t det=0; det<AliTRDgeometry::kNdet; ++det)
29   {
30     AliTRDCalROC *calROC = calPad->GetCalROC(det);
31     for (Int_t channel=0; channel<calROC->GetNchannels(); ++channel)
32       calROC->SetValue(channel, value);
33   }
34   return calPad;
35 }
36
37 TObject* CreateDetObject(const char* shortName, const char* description, Float_t value)
38 {
39   AliTRDCalDet *object = new AliTRDCalDet(shortName, description);
40   for (Int_t det=0; det<AliTRDgeometry::kNdet; ++det)
41     object->SetValue(det, value);
42   return object;
43 }
44
45 TObject* CreateGlobalsObject()
46 {
47   AliTRDCalGlobals *object = new AliTRDCalGlobals("Globals", "Global TRD calibration parameters");
48   
49   object->SetSamplingFrequency(10.0);
50   object->SetNumberOfTimeBins(22);
51   
52   return object;
53 }
54
55 TObject* CreateChamberObject()
56 {
57   AliTRDCalChamber *object = new AliTRDCalChamber("Chamber", "TRD chamber positions");
58   
59   for (Int_t det=0; det<AliTRDgeometry::kNdet; ++det)
60   {
61     object->SetPos(det, 0, 0, 0);
62     object->SetRot(det, 0, 0, 0);
63   }
64   
65   return object;
66 }
67
68 TObject* CreateStackObject()
69 {
70   AliTRDCalStack *object = new AliTRDCalStack("Stack", "TRD stack positions");
71   
72   for (Int_t sect=0; sect<AliTRDgeometry::kNsect; ++sect)
73   {
74     for (Int_t chamber=0; chamber<AliTRDgeometry::kNcham; ++chamber)
75     {
76       object->SetPos(chamber, sect, 0, 0, 0);
77       object->SetRot(chamber, sect, 0, 0, 0);
78     }
79   }
80     
81   return object;
82 }
83
84 TObject* CreatePRFWidthObject()
85 {
86   AliTRDCalPad *calPad = new AliTRDCalPad("PRFWidth","PRFWidth");
87   for (Int_t plane=0; plane<AliTRDgeometry::kNplan; ++plane)
88   {
89     Float_t value = 0;
90     switch (plane)
91     {
92       case 0: value = 0.515; break;
93       case 1: value = 0.502; break;
94       case 2: value = 0.491; break;
95       case 3: value = 0.481; break;
96       case 4: value = 0.471; break;
97       case 5: value = 0.463; break;
98       default: cout << "CreatePRFWidthObject: UNEXPECTED" << endl; return 0;
99     }
100     for (Int_t chamber=0; chamber<AliTRDgeometry::kNcham; ++chamber)
101     {
102       for (Int_t sector=0; sector<AliTRDgeometry::kNsect; ++sector)
103       {
104         AliTRDCalROC *calROC = calPad->GetCalROC(plane, chamber, sector);
105         for (Int_t channel=0; channel<calROC->GetNchannels(); ++channel)
106           calROC->SetValue(channel, value);
107       }
108     }
109   }
110       
111   return calPad;
112 }
113
114 AliTRDCalPIDLQ* CreatePIDLQObject()
115 {
116   AliTRDCalPIDLQ* pid = new AliTRDCalPIDLQ();
117   pid->ReadData("$ALICE_ROOT/TRD/TRDdEdxHistogramsV1.root");
118   pid->SetMeanChargeRatio(1.0); // The factor is the ratio of Mean of pi charge dist.
119                     // for the New TRD code divided by the Mean of pi charge
120                     // dist. given in AliTRDCalPIDLQ object
121   
122   return pid;
123 }
124
125 AliCDBMetaData* CreateMetaObject(const char* objectClassName)
126 {
127   AliCDBMetaData *md1= new AliCDBMetaData(); 
128   md1->SetObjectClassName(objectClassName);
129   md1->SetResponsible("Jan Fiete Grosse-Oetringhaus");
130   md1->SetBeamPeriod(1);
131   md1->SetAliRootVersion("05-06-00"); //root version
132   md1->SetComment("The dummy values in this calibration file are for testing only");
133   
134   return md1;
135 }
136
137 void StoreObject(const char* cdbPath, TObject* object, AliCDBMetaData* metaData)
138 {
139   AliCDBId id1(cdbPath, gkDummyRun, gkDummyRun); 
140   gStorLoc->Put(object, id1, metaData); 
141 }
142     
143
144 void AliTRDCreateDummyCDB()
145 {
146   cout << endl << "TRD :: Creating dummy CDB with event number " << gkDummyRun << endl;
147   
148   AliCDBManager *man = AliCDBManager::Instance();
149   gStorLoc = man->GetStorage("local://$ALICE_ROOT");
150   if (!gStorLoc)
151     return;
152
153   TObject* obj = 0;
154   AliCDBMetaData* metaData = 0;
155   
156   metaData = CreateMetaObject("AliTRDCalPad");
157   
158   obj = CreatePadObject("LocalVdrift","TRD drift velocities (local variations)", 1);
159   StoreObject("TRD/Calib/LocalVdrift", obj, metaData);
160   
161   obj = CreatePadObject("LocalT0","T0 (local variations)", 1);
162   StoreObject("TRD/Calib/LocalT0", obj, metaData);
163   
164   obj = CreatePadObject("GainFactor","GainFactor", 1);
165   StoreObject("TRD/Calib/GainFactor", obj, metaData);
166
167   obj = CreatePRFWidthObject();
168   StoreObject("TRD/Calib/PRFWidth", obj, metaData);
169
170   metaData = CreateMetaObject("AliTRDCalDet");
171   
172   obj = CreateDetObject("ChamberVdrift","TRD drift velocities (detector value)", 1.5);
173   StoreObject("TRD/Calib/ChamberVdrift", obj, metaData);
174   
175   obj = CreateDetObject("ChamberT0","T0 (detector value)", 0);
176   StoreObject("TRD/Calib/ChamberT0", obj, metaData);
177   
178   metaData = CreateMetaObject("AliTRDCalGlobals");
179   obj = CreateGlobalsObject();
180   StoreObject("TRD/Calib/Globals", obj, metaData);
181   
182   metaData = CreateMetaObject("AliTRDCalChamber");
183   obj = CreateChamberObject();
184   StoreObject("TRD/Calib/Chamber", obj, metaData);
185   
186   metaData = CreateMetaObject("AliTRDCalStack");
187   obj = CreateStackObject();
188   StoreObject("TRD/Calib/Stack", obj, metaData);
189   
190   metaData = CreateMetaObject("AliTRDCalPIDLQ");
191   obj = CreatePIDLQObject();
192   StoreObject("TRD/Calib/PIDLQ", obj, metaData);
193 }