]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/CalibMacros/dumpObject.C
added info on evaluatation on weighted merging of vzero event planes
[u/mrichter/AliRoot.git] / PWGPP / CalibMacros / dumpObject.C
1 /*
2   //
3   // Function to dump recursivally in human readable format, content of the objects + example use case
4   // 
5   //
6   .L  $mcProd/dumpObject.C+
7   ExampleUse(); > AliTPCClusterParam.dump
8   AliTPCRecoParam param; 
9   DumpObjectRecursive(&param); >> AliTPCRecoParam.dump
10
11   Use case examples 
12   1.) compare oontent of alignent OCDB files for differnt yers
13   2.) compare ClusterParam for different periods
14    
15 =================================================================================================================
16   // 1.)
17   // Compare alignment example:
18   // Compare TPC alignemnt 2013 and 2010
19   //
20   DumpOCDBFile("/cvmfs/alice.gsi.de/alice/data/2013/OCDB/TPC/Align/Data/Run0_999999999_v1_s0.root","TPCalign2013.dump",1);
21   DumpOCDBFile("/cvmfs/alice.gsi.de/alice/data/2010/OCDB/TPC/Align/Data/Run0_999999999_v1_s0.root","TPCalign2010.dump",1);
22   diff  TPCalign2013.dump TPCalign2010.dump > TPCalign2013_TPCalign2010.diff
23   //
24   //    
25 =================================================================================================================
26 //  2.) 
27   // Compare CluterParam OCDB etry
28   //
29   DumpOCDBFile("/cvmfs/alice.gsi.de/alice/data/2010/OCDB/TPC/Calib/ClusterParam/Run131541_999999999_v2_s0.root","2010_TPC_Calib_ClusterParam_Run131541_999999999_v2_s0.dump",1);
30   DumpOCDBFile("/cvmfs/alice.gsi.de/alice/data/2010/OCDB/TPC/Calib/ClusterParam/Run0_999999999_v1_s0.root","2010_TPC_Calib_ClusterParam_Run0_999999999_v1_s0.dump",1);
31   DumpOCDBFile("/cvmfs/alice.gsi.de/alice/data/2013/OCDB/TPC/Calib/ClusterParam/Run0_999999999_v1_s0.root","2013_TPC_Calib_ClusterParam_Run0_999999999_v1_s0.dump",1);
32   
33   diff 2010_TPC_Calib_ClusterParam_Run131541_999999999_v2_s0.dump 2010_TPC_Calib_ClusterParam_Run0_999999999_v1_s0.dump
34
35
36 */
37 #if !defined(__CINT__) || defined(__MAKECINT__)
38 #include "TMath.h"
39 #include "TFile.h"
40 #include "TTree.h"
41 #include <TVectorF.h>
42 #include <TLinearFitter.h>
43 #include <TH1F.h>
44 #include <TH3F.h>
45 #include <TProfile2D.h>
46 #include <TVectorD.h>
47 #include <TObjArray.h>
48 #include "THnBase.h"
49 #include "TRealData.h"
50 #include "TDataMember.h"
51 #include "TClass.h"
52 #include "AliCDBEntry.h"
53 #include "TObjArray.h"
54 #include "TNamed.h"
55 #include "TROOT.h"
56 #endif
57
58 void DumpObjectRecursive(TObject *obj);
59 void DumpObjectRecursive(TObject *obj, TString prefix, Int_t &counterRec);
60
61 void ExampleUse(){
62   //
63   //
64   //
65   TFile *f = TFile::Open("/cvmfs/alice.gsi.de/alice/simulation/2008/v4-15-Release/Residual/TPC/Calib/ClusterParam/Run127712_130850_v3_s0.root");
66   //TFile *f = TFile::Open("./OCDB_NoClustersBelowThreshold/TPC/Calib/RecoParam/Run127712_130850_v0_s0.root");
67   AliCDBEntry * entry = (AliCDBEntry*)f->Get("AliCDBEntry");
68   TObject  *obj = ( TObject*)entry->GetObject();
69   DumpObjectRecursive(obj);
70   //
71   //
72 }
73
74 void DumpOCDBFile(const char *finput , const char *foutput, Bool_t dumpMetaData, Bool_t xml){
75   //
76   //  
77   //  DumpOCDBFile("$ALICE_ROOT/OCDB/ITS/Align/Data/Run0_999999999_v0_s0.root", "ITS_Align_Data_Run0_999999999_v0_s0.dump")
78   //
79   if (finput==0) return ;
80   TFile *falignITS  = TFile::Open(finput);
81   AliCDBEntry *entry  = (AliCDBEntry*)falignITS->Get("AliCDBEntry");
82   if (!entry) return; 
83   TObject *obj = ((AliCDBEntry*)falignITS->Get("AliCDBEntry"))->GetObject();  
84
85   //
86   if (!xml){
87     if (dumpMetaData) gROOT->ProcessLine(TString::Format("((TObject*)%p)->Dump(); >%s",entry, foutput).Data());
88     if (!obj) return;
89     gROOT->ProcessLine(TString::Format("DumpObjectRecursive((TObject*)%p); >>%s",obj, foutput).Data());
90   }
91   if (xml){
92     TFile * f = TFile::Open(TString::Format("%s.xml",foutput).Data(),"recreate");
93     if (dumpMetaData) entry->Write("AliCDBEntry");
94     else obj->Write("AliCDBEntry");
95     f->Close();
96   }
97 }
98
99
100
101 void DumpObjectRecursive(TObject *obj){
102   //
103   //
104   //
105   Int_t counterRec=0;
106   printf("==> Dumping object at: %p, name=%s, class=%s)\n", obj, obj->GetName(), (obj->IsA()->GetName()));
107   DumpObjectRecursive(obj, TString(obj->IsA()->GetName())+".",counterRec);
108 }
109  
110 //
111 //
112 //
113 void DumpObjectRecursive(TObject *obj, TString prefix, Int_t &counterRec){
114   //
115   // Recursive dump of the TObject
116   // Dump all basic types and follow pointers to the objects
117   // current limitation:
118   //    a.) clases and structures not derived from TObject not followed (to fix)
119   //    b.) dynamic arrays not followed
120   //    c.) std maps,array ....  not followed
121   //    
122   //
123   if (!obj) return;
124   //
125   // Special case of Collection classes
126   //
127   if (obj->IsA()->InheritsFrom(TCollection::Class())) {
128     TIter myiter((TCollection*)obj);
129     TObject  *arObject=0;
130     Int_t counter=0;
131     while ((arObject = (TObject*)myiter.Next())) {
132       TString prefixArr = TString::Format("%s[%d]",prefix.Data(),counter);
133       DumpObjectRecursive(arObject,prefixArr,counterRec);
134       counter++;
135     } 
136     counterRec++;
137     return;
138   }
139
140   TClass * cl = obj->IsA();
141   if (!(cl->GetListOfRealData())) cl->BuildRealData();
142   TRealData* rd = 0;
143   TIter next(cl->GetListOfRealData());  
144   while ((rd = (TRealData*) next())) {
145     counterRec++;
146     TDataMember* dm = rd->GetDataMember();
147     TDataType* dtype = dm->GetDataType();
148     Int_t offset = rd->GetThisOffset();
149     char* pointer = ((char*) obj) + offset;
150     
151     if (dm->IsaPointer()) {
152       // We have a pointer to an object or a pointer to an array of basic types.
153       TClass* clobj = 0;
154       if (!dm->IsBasic()) {
155         clobj = TClass::GetClass(dm->GetTypeName());
156       }
157       if (clobj) {
158         // We have a pointer to an object.
159         //
160         if (!clobj->InheritsFrom(TObject::Class())) {
161           // It must be a TObject object.
162           continue; 
163         }
164         char** apointer = (char**) pointer;
165         TObject* robj = (TObject*) *apointer;
166         //      
167         if(!robj)
168           printf("M:%s%s\n",prefix.Data(),dm->GetName()); // Missing - 0 pointer
169         else{
170           printf("T:%s\t%s%s\n", clobj->GetName(),prefix.Data(), dm->GetName());
171           TString prefixNew=prefix;
172           prefixNew+=dm->GetName();
173           prefixNew+=".";
174           if (robj!=obj) DumpObjectRecursive(robj,prefixNew,counterRec);  // trivial check 
175           if (robj==obj){
176             printf("R:%s\t%s%s\n",clobj->GetName(),prefix.Data(), dm->GetName());
177           }
178         }
179       }
180     } else if (dm->IsBasic()) {
181       //
182       // Basic data type
183       //
184       const char* index = dm->GetArrayIndex();
185       if (dm->GetArrayDim()==0){
186         printf("B:\t%s%s\t%s\n", prefix.Data(),rd->GetName(), dtype->AsString(pointer));
187       }
188       //
189       // Basic array - fixed length
190       //
191       //      if (dm->GetArrayDim()>0 && strlen(index) != 0){
192       if (dm->GetArrayDim()>0 ){
193         printf("A:\t%s%s\t",prefix.Data(),rd->GetName());
194         Int_t counter=0;
195         for  (Int_t idim=0; idim<dm->GetArrayDim(); idim++){
196           //printf("A:%d\t%d\n", dm->GetArrayDim(),dm->GetMaxIndex(idim));
197           for (Int_t j=0; j<dm->GetMaxIndex(idim); j++){
198             printf("%s\t",dtype->AsString(pointer+dm->GetUnitSize()*counter));
199             counter++;
200             if (counter%5==0) printf("\nA:\t%s%s\t",prefix.Data(),rd->GetName());
201           }
202         }
203         printf("\n");
204       }
205       //
206       // Basic array - dynamic length
207       //
208       if (dm->GetArrayDim()>0 && strlen(index) != 0){
209         //
210         // Dump first only for the moment
211         //  
212         printf("B:\t%s%s\t%s\n",prefix.Data(),rd->GetName(), dtype->AsString(pointer));
213       }
214     } else {
215     }
216   }
217 }  
218
219 //
220 // Small checks to test the TRealData and TDataType
221 //
222
223
224
225 void DumpDataSimple(){
226   //
227   // Dump example for elenatr data types 
228   //
229   TObject *obj = new TVectorD(20);
230   TClass * cl = obj->IsA();
231   if (!cl->GetListOfRealData()) cl->BuildRealData();
232   //
233   TRealData* rd = 0;
234   rd = (TRealData*)(cl->GetListOfRealData()->FindObject("fNrows"));
235   TDataMember* dm = rd->GetDataMember();
236   TDataType* dtype = dm->GetDataType();
237   //
238   Int_t offset = rd->GetThisOffset();
239   char* pointer = ((char*) obj) + offset;
240   printf("%s\n",dtype->AsString(pointer));
241   
242 }
243
244 void DumpDataArray(){
245   //
246   // print array example
247   // 
248   TObject *obj = new TVectorD(20);
249   TClass * cl = obj->IsA();
250   if (!cl->GetListOfRealData()) cl->BuildRealData();
251   TRealData* rd = 0;
252   rd = (TRealData*)(cl->GetListOfRealData()->FindObject("*fElements"));
253   TDataMember* dm = rd->GetDataMember();
254   TDataType* dtype = dm->GetDataType();
255   dtype->Print();
256   //
257   Int_t offset = rd->GetThisOffset();
258   char* pointer = ((char*) obj) + offset; 
259   printf("%s\n",dtype->AsString(pointer));
260 }
261
262 void DumpTObjectArray(){
263   //
264   //
265   //
266   TObjArray *array = new TObjArray(10);
267   for (Int_t i=0; i<10; i++) array->AddLast(new TNamed(Form("n%d",i), Form("n%d",i)));  
268   DumpObjectRecursive(array);
269   //
270   //
271   TObject *obj = array;
272   TClass * cl = obj->IsA();
273   if (!cl->GetListOfRealData()) cl->BuildRealData();
274   TRealData* rd = 0;
275   rd = (TRealData*)(cl->GetListOfRealData()->FindObject("*fCont"));
276   TDataMember* dm = rd->GetDataMember();
277   TDataType* dtype = dm->GetDataType();
278   //
279   Int_t offset = rd->GetThisOffset();
280   char* pointer = ((char*) obj) + offset;
281   char** apointer = (char**) pointer;
282   //we have pointer to pointer here
283   TObject** ppobj = (TObject**) *apointer;
284   (*ppobj)->Print();
285   //
286   TIter myiter(array);
287   TObject  *arObject; 
288   dtype->Print();
289   while ((arObject = (TObject*)myiter.Next())) {
290     DumpObjectRecursive(arObject);
291   } 
292
293
294 }