]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/macrosSDD/ITSQArecoparam.C
defects from coverity fixed
[u/mrichter/AliRoot.git] / ITS / macrosSDD / ITSQArecoparam.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include <TFile.h>
3 #include <TAlienFile.h>
4 #include <TSystem.h>
5 #include <TGrid.h>
6 #include <TGridResult.h>
7 #include <Riostream.h>
8 #include <TObjArray.h>
9 #include <TKey.h>
10 #include <TStyle.h>
11 #include <TTree.h>
12 #include <TROOT.h>
13 #include <TPluginManager.h>
14 #include <TClass.h>
15 #include <iostream>
16 #include <stdio.h>
17 #include <TGeoGlobalMagField.h>
18 #include <TMap.h>
19 #include "AliRawReader.h"
20 #include "AliRawReaderRoot.h"
21 #include "AliRawReaderDate.h"
22 #include "AliITSQADataMakerRec.h"
23 #include "AliITSQASDDDataMakerRec.h"
24 #include "AliITSQAChecker.h"
25 #include "AliQAChecker.h"
26 #include "AliITSQASDDChecker.h"
27 #include "AliReconstructor.h"
28 #include "AliCDBManager.h"
29 #include "AliQAv1.h"
30 #include "AliMagF.h"
31 #include "AliGeomManager.h"
32 #include "AliITSInitGeometry.h"
33 #include "AliITSgeom.h"
34 #include "AliRecoParam.h"
35 #include "AliCDBPath.h"
36 #include "AliCDBEntry.h"
37 #include "AliRecoParam.h"
38 #include "AliDetectorRecoParam.h"
39 #include "AliITSReconstructor.h"
40 #include "AliITSRecPointContainer.h"
41 #include "AliLog.h"
42 #include "AliGRPObject.h"
43 #include "AliRunInfo.h"
44 #endif
45 void ITSQArecoparam(char *iFile,Int_t idet=2,Int_t FirstEvt=0, Int_t MaxEvts=1000000)
46 {
47   TString namefile(iFile);
48   if(namefile.Contains("alien"))
49     {
50       TGrid::Connect("alien://"); 
51       if(!gGrid) {
52         printf("gGrid not found! exit macro\n");
53         return;
54       }
55     }
56   gStyle->SetPalette(1);
57   Int_t ic=0;
58   AliRawReader *rd=NULL; 
59   if(strstr(iFile,".root")!=0){rd = new AliRawReaderRoot(iFile,FirstEvt);}
60   else{rd=new AliRawReaderDate(iFile,FirstEvt);}
61   Int_t runNumber = rd->GetRunNumber();
62   cout << "ITS Quality Assurance Prototype" << endl; 
63   //TStopwatch mytimer;
64   //TString namefile(iFile);
65   // Set OCDB if needed
66   AliCDBManager* man = AliCDBManager::Instance();
67   if (!man->IsDefaultStorageSet()) {
68     printf("Setting a local default storage and run number \n");
69     if(namefile.Contains("alien")){
70       man->SetDefaultStorage("raw://");
71     }
72     else{man->SetDefaultStorage("local://$ALICE_ROOT/OCDB");}
73     man->SetRun(runNumber);
74     AliQAv1::SetQARefStorage("local://$ALICE_ROOT/QARef") ;
75   }
76
77   AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
78
79   AliGRPObject *fGRPData=NULL;
80   
81   TMap* m = dynamic_cast<TMap*>(entry->GetObject());  // old GRP entry
82
83     if (m) {
84        printf("Found a TMap in GRP/GRP/Data, converting it into an AliGRPObject");
85        m->Print();
86        fGRPData = new AliGRPObject();
87        fGRPData->ReadValuesFromMap(m);
88     }
89
90     else {
91        printf("Found an AliGRPObject in GRP/GRP/Data, reading it");
92        fGRPData = dynamic_cast<AliGRPObject*>(entry->GetObject());  // new GRP entry
93        entry->SetOwner(0);
94     }
95
96
97  if (!fGRPData) {
98      printf("Error   No GRP entry found in OCDB!");
99      return;
100   }
101
102
103   TString lhcState = fGRPData->GetLHCState();
104   if (lhcState==AliGRPObject::GetInvalidString()) {
105     printf("Error  GRP/GRP/Data entry:  missing value for the LHC state ! Using UNKNOWN");
106     lhcState = "UNKNOWN";
107   }
108
109   TString beamType = fGRPData->GetBeamType();
110   if (beamType==AliGRPObject::GetInvalidString()) {
111     printf("Error   GRP/GRP/Data entry:  missing value for the beam type ! Using UNKNOWN");
112     beamType = "UNKNOWN";
113   }
114
115   Float_t beamEnergy = fGRPData->GetBeamEnergy();
116   if (beamEnergy==AliGRPObject::GetInvalidFloat()) {
117     printf("Error   GRP/GRP/Data entry:  missing value for the beam energy ! Using 0");
118     beamEnergy = 0;
119   }
120
121   TString runType = fGRPData->GetRunType();
122   if (runType==AliGRPObject::GetInvalidString()) {
123     printf("Error  GRP/GRP/Data entry:  missing value for the run type ! Using UNKNOWN");
124     runType = "UNKNOWN";
125   }
126
127   Int_t activeDetectors = fGRPData->GetDetectorMask();
128   if (activeDetectors==AliGRPObject::GetInvalidUInt()) {
129     printf("Error  GRP/GRP/Data entry:  missing value for the detector mask ! Using 1074790399");
130     activeDetectors = 1074790399;
131   }
132
133   AliRunInfo *fRunInfo = new AliRunInfo(lhcState, beamType, beamEnergy, runType, activeDetectors);
134   fRunInfo->Dump();
135
136  if ( TGeoGlobalMagField::Instance()->IsLocked() ) {
137     if (TGeoGlobalMagField::Instance()->GetField()->TestBit(AliMagF::kOverrideGRP)) {
138       printf("ExpertMode!!! GRP information will be ignored !");
139       printf("ExpertMode!!! Running with the externally locked B field !");
140     }
141     else {
142       printf("Destroying existing B field instance!");
143       delete TGeoGlobalMagField::Instance();
144     }    
145   }
146   if ( !TGeoGlobalMagField::Instance()->IsLocked() ) {
147     // Construct the field map out of the information retrieved from GRP.
148     Bool_t ok = kTRUE;
149     // L3
150     Float_t l3Current = fGRPData->GetL3Current((AliGRPObject::Stats)0);
151     if (l3Current == AliGRPObject::GetInvalidFloat()) {
152       printf("Error :  GRP/GRP/Data entry:  missing value for the L3 current !");
153       ok = kFALSE;
154     }
155     
156     Char_t l3Polarity = fGRPData->GetL3Polarity();
157     if (l3Polarity == AliGRPObject::GetInvalidChar()) {
158       printf("Error:   GRP/GRP/Data entry:  missing value for the L3 polarity !");
159       ok = kFALSE;
160     }
161
162    // Dipole
163     Float_t diCurrent = fGRPData->GetDipoleCurrent((AliGRPObject::Stats)0);
164     if (diCurrent == AliGRPObject::GetInvalidFloat()) {
165       printf("Error  GRP/GRP/Data entry:  missing value for the dipole current !");
166       ok = kFALSE;
167     }
168
169     Char_t diPolarity = fGRPData->GetDipolePolarity();
170     if (diPolarity == AliGRPObject::GetInvalidChar()) {
171       printf("Error GRP/GRP/Data entry:  missing value for the dipole polarity !");
172       ok = kFALSE;
173     }
174
175
176    Int_t  polConvention = fGRPData->IsPolarityConventionLHC() ? AliMagF::kConvLHC : AliMagF::kConvDCS2008;
177    Bool_t uniformB = fGRPData->IsUniformBMap();
178    
179    if (ok) { 
180      AliMagF* fld = AliMagF::CreateFieldMap(TMath::Abs(l3Current) * (l3Polarity ? -1:1), 
181                                             TMath::Abs(diCurrent) * (diPolarity ? -1:1), 
182                                             polConvention,uniformB,beamEnergy, beamType.Data());
183      if (fld) {
184        TGeoGlobalMagField::Instance()->SetField( fld );
185        TGeoGlobalMagField::Instance()->Lock();
186        printf("Running with the B field constructed out of GRP !");
187      }
188      else printf("Fatal Failed to create a B field map !");
189    }
190    else printf("Fatal  B field is neither set nor constructed from GRP ! Exitig...");
191   }
192   
193
194
195
196   AliITSQADataMakerRec *itsQAdm = new AliITSQADataMakerRec(kTRUE,idet,0); //online kTRUE
197   itsQAdm->SetWriteExpert() ;
198   itsQAdm->SetRunNumber(runNumber);  
199   //________________________For the RecPoints____________________________________
200   /************************************************/
201   TPluginManager* pluginManager=NULL;
202   TPluginHandler* pluginHandler=NULL; 
203   AliReconstructor* reconstructor = NULL;
204   AliITSRecPointContainer* rpcont=NULL;
205   AliGeomManager::LoadGeometry("geometry.root");
206   AliGeomManager::ApplyAlignObjsFromCDB("ITS");
207   //   ITS initializations
208   
209   AliITSInitGeometry initgeom;
210   AliITSgeom *geom = initgeom.CreateAliITSgeom();
211   printf("Geometry name: %s\n",(initgeom.GetGeometryName()).Data());
212   
213   printf("Loading reconstruction parameter objects for detector ITS\n");
214   AliRecoParam fRecoParam; 
215   AliCDBPath path("ITS","Calib","RecoParam");
216   AliCDBEntry *entry2=AliCDBManager::Instance()->Get(path.GetPath());
217   Bool_t cacheStatus = AliCDBManager::Instance()->GetCacheFlag();
218   if(!entry2){printf("Couldn't find RecoParam entry in OCDB for detector ITS");entry2=NULL;}
219   else {
220     TObject *recoParamObj = entry2->GetObject();
221     if (dynamic_cast<TObjArray*>(recoParamObj)) {
222       printf("RecoParam TObjArray\n");
223       fRecoParam.AddDetRecoParamArray(0,dynamic_cast<TObjArray*>(recoParamObj));
224     }
225     else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
226       printf("RecoParam AliDetectorRecoParam\n");
227       printf("Single set of reconstruction parameters found for detector ITS");
228       (dynamic_cast<AliDetectorRecoParam*>(recoParamObj))->SetAsDefault();
229       fRecoParam.AddDetRecoParam(0,(dynamic_cast<AliDetectorRecoParam*>(recoParamObj)));
230     }
231     else {printf("Error: No valid RecoParam object found in the OCDB for detector ITS");}
232     entry2->SetOwner(0);
233   }
234   if(!cacheStatus)entry2->SetObject(NULL);
235   if(!cacheStatus){ delete entry2;}
236   
237   // load the reconstructor object
238  pluginManager = gROOT->GetPluginManager();
239   TString detName = "ITS";
240   TString recName = "Ali" + detName + "Reconstructor";  
241   
242 pluginHandler = pluginManager->FindHandler("AliReconstructor", "ITS");
243   // if not, add a plugin for it
244   if (!pluginHandler) {
245     printf("defining plugin for ITS\n");
246     TString libs = gSystem->GetLibraries();
247     if (libs.Contains("lib" + detName + "base.so") ||
248         (gSystem->Load("lib" + detName + "base.so") >= 0)) {pluginManager->AddHandler("AliReconstructor", detName,recName, detName + "rec", recName + "()");}
249     else {pluginManager->AddHandler("AliReconstructor", detName,recName, detName, recName + "()");}
250     pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
251   }
252   if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);}
253   if (fRecoParam.GetDetRecoParamArray(0) && !AliReconstructor::GetRecoParam(0)) {
254     const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(0);
255     reconstructor->Init();
256     reconstructor->SetRecoParam(par);
257   }
258   
259   /*AliITSRecPointContainer**/ rpcont=AliITSRecPointContainer::Instance();
260   rpcont->PrepareToRead();
261
262   Int_t cycleLength = 5;  
263   //cout << "Processing Run " << runNumber << endl;
264   cout << "Init: " << AliQAv1::kRAWS << endl;
265   
266   TObjArray **objArray= itsQAdm->Init(AliQAv1::kRAWS, cycleLength);
267   cout<<"raw tobjarray :"<<objArray<<"\n"<<endl;
268   for(Int_t spec = 0 ; spec < 5 ; spec++){
269     if(spec==1){
270       itsQAdm->SetEventSpecie(AliRecoParam::kLowMult);
271       AliQAv1::Instance()->SetEventSpecie(AliRecoParam::kLowMult);
272       itsQAdm->InitRaws();
273     }
274     else{continue;}
275   }
276   itsQAdm->StartOfCycle(AliQAv1::kRAWS,runNumber,kFALSE);
277   /*********************************************************************/
278
279     cout << "Init: " << AliQAv1::kRECPOINTS << endl;
280     TObjArray **objArray1=  itsQAdm->Init(AliQAv1::kRECPOINTS, cycleLength);
281     cout<<"recpoint tobjarray :"<<objArray1<<"\n"<<endl;
282     for(Int_t spec = 0 ; spec < 5 ; spec++){
283       if(spec==1){
284         itsQAdm->SetEventSpecie(AliRecoParam::kLowMult);
285         AliQAv1::Instance()->SetEventSpecie(AliRecoParam::kLowMult);
286         itsQAdm->InitRecPoints();
287       }
288       else{continue;}
289     }
290     itsQAdm->StartOfCycle(AliQAv1::kRECPOINTS,runNumber,kTRUE);
291
292   /*********************************************************************/
293   Int_t iev = 0;
294   while(rd->NextEvent() && iev < MaxEvts ) {
295     cout<<">>>>>>>   Processing event number: "<<++iev<<endl;
296     /*******************************************************/
297     if(itsQAdm->IsCycleDone()) {
298       cout << "end of cycle" << endl;
299       AliQAChecker::Instance()->SetRunNumber(AliCDBManager::Instance()->GetRun());      
300       itsQAdm->EndOfCycle(AliQAv1::kRAWS);
301       itsQAdm->StartOfCycle(AliQAv1::kRAWS,ic++,kFALSE);
302     }
303     /******************************************************/
304     /*************************************************/
305
306       if(itsQAdm->IsCycleDone()) {
307         cout << "end of cycle" << endl;
308         AliQAChecker::Instance()->SetRunNumber(AliCDBManager::Instance()->GetRun());
309         itsQAdm->EndOfCycle(AliQAv1::kRECPOINTS);
310         itsQAdm->StartOfCycle(AliQAv1::kRECPOINTS,ic,kTRUE);
311       } 
312  
313     /*************************************************/
314     cout<<"Beginning Exec"<<endl;
315     cout<<"AliQAv1::kRAWS   "<<AliQAv1::kRAWS<<endl;
316     itsQAdm->SetEventSpecie(AliRecoParam::kLowMult);
317     AliQAv1::Instance()->SetEventSpecie(AliRecoParam::kLowMult);
318     itsQAdm->Exec(AliQAv1::kRAWS,rd);
319     /***************************************************/
320     //  if(kRecPoints){
321       cout<<"AliQAv1::kRECPOINTS   "<<AliQAv1::kRECPOINTS<<endl;
322       cout << "DigitsToRecPoints" << endl;
323       TTree* fTreeR = new TTree("TreeR", "Reconstructed Points Container"); //make a tree
324       Char_t option[5];
325
326       if(idet==0)sprintf(option,"ALL");
327       else if(idet==1)sprintf(option,"SPD");
328       else if(idet==2)sprintf(option,"SDD");
329       else if(idet==3)sprintf(option,"SSD");
330       printf("\t\t===========>option is %s\n",option);
331
332       rpcont->PrepareToRead();
333       reconstructor->Reconstruct(rd,fTreeR);
334
335       itsQAdm->SetEventSpecie(AliRecoParam::kLowMult);
336       AliQAv1::Instance()->SetEventSpecie(AliRecoParam::kLowMult)  ; 
337       itsQAdm->Exec(AliQAv1::kRECPOINTS,fTreeR);
338       cout<<"Finishing Exec"<<endl;
339     
340       ((AliITSReconstructor*)reconstructor)->ResetRecPoints();
341       delete fTreeR; 
342       fTreeR=NULL; 
343       /*****************************************************/
344     }
345
346   cout << "end RAWS cycle: " << AliQAv1::kRAWS << endl;
347   cout << "refStorage: " << AliQAv1::GetQARefStorage() << endl;
348   cout << "end of cycle 2" << endl;
349   AliQAChecker::Instance()->SetRunNumber(AliCDBManager::Instance()->GetRun());
350   itsQAdm->EndOfCycle(AliQAv1::kRAWS);   
351   cout << "Raws QA completed for " << iev << " events" << endl;
352   /*******************************************************************/
353
354     AliQAChecker::Instance()->SetRunNumber(AliCDBManager::Instance()->GetRun());
355     itsQAdm->EndOfCycle(AliQAv1::kRECPOINTS);
356     cout << "RecPoints QA completed for " << iev << " events" << endl;
357
358   /*******************************************************************/
359   itsQAdm->Finish(); // write to the output File
360
361   cout << "Call AliITSQASDDDataMakerRec destructor" << endl;
362   delete itsQAdm;
363   itsQAdm=NULL;
364
365 }
366