]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/TOFPreprocessor.C
New default values for baselines (F.Prino)
[u/mrichter/AliRoot.git] / TOF / TOFPreprocessor.C
1 /* 
2 $Id$ 
3 */
4
5 // This class runs the test TOF preprocessor
6 // It uses AliTestShuttle to simulate a full Shuttle process
7
8 // The input data is created in the functions
9 //   CreateDCSAliasMap() creates input that would in the same way come from DCS
10 //   ReadDCSAliasMap() reads from a file
11 //   CreateInputFilesMap() creates a list of local files, that can be accessed by the shuttle
12
13 extern TBenchmark *gBenchmark;
14 void TOFPreprocessor(Char_t * RunType="PHYSICS")
15 {
16   gSystem->Load("$ALICE_ROOT/SHUTTLE/TestShuttle/libTestShuttle.so");
17
18   AliLog::SetClassDebugLevel("AliTOFPreprocessor",2);
19   // initialize location of CDB
20   AliTestShuttle::SetMainCDB("local://$ALICE_ROOT/SHUTTLE/TestShuttle/TestCDB");
21   AliTestShuttle::SetMainRefStorage("local://$ALICE_ROOT/SHUTTLE/TestShuttle/TestReference");
22
23   // create AliTestShuttle instance
24   Int_t nrun = 6;
25   AliTestShuttle* shuttle = new AliTestShuttle(nrun, 0, 1000);
26   //setting run type to physiscs
27   shuttle->SetInputRunType(RunType);
28
29   // Generation of "fake" input DCS data
30   TMap* dcsAliasMap = CreateDCSAliasMap();  
31
32   // now give the alias map to the shuttle
33   shuttle->SetDCSInput(dcsAliasMap);   
34
35   // processing files. for the time being, the files are local.
36   shuttle->AddInputFile(AliTestShuttle::kDAQ, "TOF", "DELAYS", "MON", "$ALICE_ROOT/TOF/ShuttleInput/Total.root");
37   shuttle->AddInputFile(AliTestShuttle::kDAQ, "TOF", "RUNLevel", "MON", "$ALICE_ROOT/TOF/ShuttleInput/Partial.root");
38   shuttle->AddInputFile(AliTestShuttle::kDCS, "TOF", "TofFeeMap", "", "$ALICE_ROOT/TOF/ShuttleInput/TOFFEE.20080310.163032.4000");
39   char filename[100];
40   char LDCname[5];
41
42   for (Int_t iLDC=0;iLDC<2;iLDC++){
43     sprintf(filename,"$ALICE_ROOT/TOF/ShuttleInput/TOFoutPulserLDC_%02i.root",iLDC*2);
44     sprintf(LDCname,"LDC%i",iLDC*2);
45     shuttle->AddInputFile(AliTestShuttle::kDAQ, "TOF", "PULSER", LDCname, filename);
46     sprintf(filename,"$ALICE_ROOT/TOF/ShuttleInput/TOFoutNoiseLDC_%02i.root",iLDC*2);
47     sprintf(LDCname,"LDC%i",iLDC*2);
48     shuttle->AddInputFile(AliTestShuttle::kDAQ, "TOF", "NOISE", LDCname, filename);
49   }
50
51   // instantiation of the preprocessor
52   AliPreprocessor* pp = new AliTOFPreprocessor(shuttle);
53
54   // preprocessing
55   gBenchmark->Start("process");
56   shuttle->Process();
57   gBenchmark->Stop("process");
58   gBenchmark->Print("process");
59
60   // checking the file which should have been created  
61   AliCDBEntry* chkEntry = AliCDBManager::Instance()->GetStorage(AliShuttleInterface::GetMainCDB())->Get("TOF/Calib/ParOnlineDelay", nrun);
62   if (!chkEntry)
63   {
64     printf("The file is not there. Something went wrong.\n");
65     return;
66   }
67
68   AliTOFDataDCS* output = dynamic_cast<AliTOFDataDCS*> (chkEntry->GetObject());
69   // If everything went fine, draw the result
70   if (output)
71     printf("Output found.\n");
72   //    output->Draw();
73   
74 }
75
76 TMap* CreateDCSAliasMap()
77 {
78   // Creates a DCS structure
79   // The structure is the following:
80   //   TMap (key --> value)
81   //     <DCSAlias> --> <valueList>
82   //     <DCSAlias> is a string
83   //     <valueList> is a TObjArray of AliDCSValue
84   //     An AliDCSValue consists of timestamp and a value in form of a AliSimpleValue
85
86   // In this example 6 aliases exists: DCSAlias1 ... DCSAlias6
87   // Each contains 1000 values randomly generated by TRandom::Gaus + 5*nAlias
88
89   TMap* aliasMap = new TMap;
90   aliasMap->SetOwner(1);
91
92   TRandom random;
93   TDatime *datime = new TDatime();
94   Int_t time = datime->GetTime();
95   Int_t date = datime->GetDate();
96   Int_t pid  = gSystem->GetPid();
97   delete datime;
98   Int_t iseed = TMath::Abs(10000 * pid + time - date); 
99
100   //Float_t thrHVv=7.75, thrHVc=3, thrLVv=2.75, thrLVc=2.5,
101   //thrFEEthr=1.5, thrFEEt=10, thrTemp=35, thrPress=1000;
102   //Float_t tentHVv=6.75, tentHVc=2, tentLVv=1.75, tentLVc=1.5,
103   //  tentFEEthr=0.5, te    result=0;
104   //ntFEEt=20, tentTemp=25, tentPress=900;
105   //Float_t sigmaHVv=1, sigmaHVc=0.25, sigmaLVv=0.25, sigmaLVc=0.25,
106   //  sigmaFEEthr=0.05, sigmaFEEt=5, sigmaTemp=1, sigmaPress=10;
107
108   Float_t tentHVv=6500, tentHVi=80;
109   Float_t sigmaHVv=10, sigmaHVi=10;
110
111   Float_t tent=0, sigma=0, thr=0;
112   // to have all the aliases, decomment the following line:
113   Int_t NAliases=360, NHV=90;
114
115   // if not all the aliases are there, use this:
116   //Int_t NAliases=120, NHV=90;
117
118   for(int nAlias=0;nAlias<NAliases;nAlias++) {
119
120     TObjArray* valueSet = new TObjArray;
121     valueSet->SetOwner(1);
122
123     TString sindex;
124     TString aliasName;
125     if (nAlias<NHV){
126       aliasName = "tof_hv_vp_";
127       sindex.Form("%02i",nAlias);
128       aliasName += sindex;
129       //aliasName += nAlias;
130       tent=tentHVv;
131       sigma=sigmaHVv;
132       //      thr=thrHVv;
133     }
134     else if (nAlias<NHV*2){
135       //      aliasName = "HVvneg";
136       //aliasName += nAlias-NHV;
137       aliasName = "tof_hv_vn_";
138       sindex.Form("%02i",nAlias-NHV);
139       aliasName += sindex;
140       tent=-tentHVv;
141       sigma=-sigmaHVv;
142       //thr=-thrHVv;
143     }
144     else if (nAlias<NHV*3){
145       //      aliasName = "HVcpos";
146       //aliasName += nAlias-2*NHV;
147       aliasName = "tof_hv_ip_";
148       sindex.Form("%02i",nAlias-2*NHV);
149       aliasName += sindex;
150       tent=tentHVi;
151       sigma=sigmaHVi;
152       //thr=thrHVc;
153     }
154     else if (nAlias<NHV*4){
155       //      aliasName = "HVcneg";
156       //aliasName += nAlias-3*NHV;
157       aliasName = "tof_hv_in_";
158       sindex.Form("%02i",nAlias-3*NHV);
159       aliasName += sindex;
160       tent=-tentHVi;
161       sigma=-sigmaHVi;
162       //thr=-thrHVc;
163     }
164     // gauss generation of values 
165     for (int timeStamp=0;timeStamp<1000;timeStamp+=10){
166     //for (int timeStamp=0;timeStamp<1;timeStamp++){
167       Float_t gaussvalue = (Float_t) (random.Gaus(tent,sigma));
168       if (TMath::Abs(gaussvalue-tent)>sigma){
169         AliDCSValue* dcsVal = new AliDCSValue(gaussvalue, timeStamp);
170         valueSet->Add(dcsVal);
171       }
172     }
173
174     aliasMap->Add(new TObjString(aliasName), valueSet);
175   }
176
177   return aliasMap;
178 }
179
180 TMap* ReadDCSAliasMap()
181 {
182   // Open a file that contains DCS input data
183   // The CDB framework is used to open the file, this means the file is located
184   // in $ALICE_ROOT/SHUTTLE/TestShuttle/TestCDB/<detector>/DCS/Data
185   // The file contains an AliCDBEntry that contains a TMap with the DCS structure.
186   // An explanation of the structure can be found in CreateDCSAliasMap()
187
188   AliCDBEntry *entry = AliCDBManager::Instance()->Get("TOF/DCS/Data", 0);
189   return dynamic_cast<TMap*> (entry->GetObject());
190 }
191
192 void WriteDCSAliasMap()
193 {
194   // This writes the output from CreateDCSAliasMap to a CDB file
195
196   TMap* dcsAliasMap = CreateDCSAliasMap();
197
198   AliCDBMetaData metaData;
199         metaData.SetBeamPeriod(0);
200         metaData.SetResponsible("Chiara");
201         metaData.SetComment("Test object for TOFPreprocessor.C");
202
203   AliCDBId id("TOF/DCS/Data", 0, 0);
204
205   // initialize location of CDB
206   AliCDBManager::Instance()->SetDefaultStorage("local://$ALICE_ROOT/SHUTTLE/TestShuttle/TestCDB");
207
208   AliCDBManager::Instance()->Put(dcsAliasMap, id, &metaData);
209 }