]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Dump object functionality
authormivanov <marian.ivanov@cern.ch>
Wed, 9 Apr 2014 15:10:32 +0000 (17:10 +0200)
committermivanov <marian.ivanov@cern.ch>
Wed, 9 Apr 2014 15:10:32 +0000 (17:10 +0200)
STEER/CDB/AliOCDBtoolkit.cxx
STEER/CDB/AliOCDBtoolkit.h

index eb4aab1fc520c35e888743dbcf3a26f8c2730372..e8b26fc55988c28064620cf81d340bee265b169b 100644 (file)
   Implementation assumption:
   All input formats (a .. f) will  be converted to the TMap storages and TList if AliCDBIds 
 
-
   Example usage:
   AliOCDBtoolkit::MakeDiffExampleUseCase();
 
 
+
+  //=============================================================================
+  // Functionality to dump content of  objects in human readable format
+  //=============================================================================
+  Use case examples 
+  1.) compare oontent of alignent OCDB files for differnt yers
+  2.) compare ClusterParam for different periods
+  
+  
+   
+  =================================================================================================================
+  // 1.)
+  // Compare alignment example:
+  // Compare TPC alignemnt 2013 and 2010
+  //
+  AliOCDBtoolkit::DumpOCDBFile("/cvmfs/alice.gsi.de/alice/data/2013/OCDB/TPC/Align/Data/Run0_999999999_v1_s0.root","TPCalign2013.dump",1);
+  AliOCDBtoolkit::DumpOCDBFile("/cvmfs/alice.gsi.de/alice/data/2010/OCDB/TPC/Align/Data/Run0_999999999_v1_s0.root","TPCalign2010.dump",1);
+  diff  TPCalign2013.dump TPCalign2010.dump > TPCalign2013_TPCalign2010.diff
+  //
+  //    
+  =================================================================================================================
+  //  2.) 
+  // Compare CluterParam OCDB etry
+  //
+  AliOCDBtoolkit::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);
+  AliOCDBtoolkit:: AliOCDBtoolkit::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);
+  AliOCDBtoolkit::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);
+  diff 2010_TPC_Calib_ClusterParam_Run131541_999999999_v2_s0.dump 2010_TPC_Calib_ClusterParam_Run0_999999999_v1_s0.dump
+
 */
 
 
@@ -50,6 +79,11 @@ using namespace std;
 #include <sstream>
 #include <stdexcept>
 #include <functional>
+#include "TRealData.h"
+#include "TDataMember.h"
+#include "TClass.h"
+#include "TROOT.h"
+#include <TVectorD.h>
 //
 #include "TSystem.h"
 #include "TObjArray.h"
@@ -401,3 +435,231 @@ void AliOCDBtoolkit::DumpOCDB(const TMap *cdbMap0, const TList *cdbList0){
     delete file;
   }
 }
+
+
+//====================================================================================================
+//  Dump object part
+//==================================================================================================== 
+
+
+
+
+
+void AliOCDBtoolkit::DumpOCDBFile(const char *finput , const char *foutput, Bool_t dumpMetaData, Bool_t xml){
+  //
+  //  
+  //  DumpOCDBFile("$ALICE_ROOT/OCDB/ITS/Align/Data/Run0_999999999_v0_s0.root", "ITS_Align_Data_Run0_999999999_v0_s0.dump")
+  //
+  if (finput==0) return ;
+  TFile *falignITS  = TFile::Open(finput);
+  AliCDBEntry *entry  = (AliCDBEntry*)falignITS->Get("AliCDBEntry");
+  if (!entry) return; 
+  TObject *obj = ((AliCDBEntry*)falignITS->Get("AliCDBEntry"))->GetObject();  
+
+  //
+  if (!xml){
+    if (dumpMetaData) gROOT->ProcessLine(TString::Format("((TObject*)%p)->Dump(); >%s",entry, foutput).Data());
+    if (!obj) return;
+    gROOT->ProcessLine(TString::Format("DumpObjectRecursive((TObject*)%p); >>%s",obj, foutput).Data());
+  }
+  if (xml){
+    TFile * f = TFile::Open(TString::Format("%s.xml",foutput).Data(),"recreate");
+    if (dumpMetaData) entry->Write("AliCDBEntry");
+    else obj->Write("AliCDBEntry");
+    f->Close();
+  }
+}
+
+
+
+void AliOCDBtoolkit::DumpObjectRecursive(TObject *obj){
+  //
+  //
+  //
+  Int_t counterRec=0;
+  printf("==> Dumping object at: %p, name=%s, class=%s)\n", obj, obj->GetName(), (obj->IsA()->GetName()));
+  DumpObjectRecursive(obj, TString(obj->IsA()->GetName())+".",counterRec);
+}
+//
+//
+//
+void AliOCDBtoolkit::DumpObjectRecursive(TObject *obj, TString prefix, Int_t &counterRec){
+  //
+  // Recursive dump of the TObject
+  // Dump all basic types and follow pointers to the objects
+  // current limitation:
+  //    a.) clases and structures not derived from TObject not followed (to fix)
+  //    b.) dynamic arrays not followed
+  //    c.) std maps,array ....  not followed
+  //    
+  //
+  if (!obj) return;
+  //
+  // Special case of Collection classes
+  //
+  if (obj->IsA()->InheritsFrom(TCollection::Class())) {
+    TIter myiter((TCollection*)obj);
+    TObject  *arObject=0;
+    Int_t counter=0;
+    while ((arObject = (TObject*)myiter.Next())) {
+      TString prefixArr = TString::Format("%s[%d]",prefix.Data(),counter);
+      DumpObjectRecursive(arObject,prefixArr,counterRec);
+      counter++;
+    } 
+    counterRec++;
+    return;
+  }
+
+  TClass * cl = obj->IsA();
+  if (!(cl->GetListOfRealData())) cl->BuildRealData();
+  TRealData* rd = 0;
+  TIter next(cl->GetListOfRealData());  
+  while ((rd = (TRealData*) next())) {
+    counterRec++;
+    TDataMember* dm = rd->GetDataMember();
+    TDataType* dtype = dm->GetDataType();
+    Int_t offset = rd->GetThisOffset();
+    char* pointer = ((char*) obj) + offset;
+    
+    if (dm->IsaPointer()) {
+      // We have a pointer to an object or a pointer to an array of basic types.
+      TClass* clobj = 0;
+      if (!dm->IsBasic()) {
+       clobj = TClass::GetClass(dm->GetTypeName());
+      }
+      if (clobj) {
+       // We have a pointer to an object.
+       //
+       if (!clobj->InheritsFrom(TObject::Class())) {
+         // It must be a TObject object.
+         continue; 
+       }
+       char** apointer = (char**) pointer;
+       TObject* robj = (TObject*) *apointer;
+       //      
+       if(!robj)
+         printf("M:%s%s\n",prefix.Data(),dm->GetName()); // Missing - 0 pointer
+       else{
+         printf("T:%s\t%s%s\n", clobj->GetName(),prefix.Data(), dm->GetName());
+         TString prefixNew=prefix;
+         prefixNew+=dm->GetName();
+         prefixNew+=".";
+         if (robj!=obj) DumpObjectRecursive(robj,prefixNew,counterRec);  // trivial check 
+         if (robj==obj){
+           printf("R:%s\t%s%s\n",clobj->GetName(),prefix.Data(), dm->GetName());
+         }
+       }
+      }
+    } else if (dm->IsBasic()) {
+      //
+      // Basic data type
+      //
+      const char* index = dm->GetArrayIndex();
+      if (dm->GetArrayDim()==0){
+       printf("B:\t%s%s\t%s\n", prefix.Data(),rd->GetName(), dtype->AsString(pointer));
+      }
+      //
+      // Basic array - fixed length
+      //
+      //      if (dm->GetArrayDim()>0 && strlen(index) != 0){
+      if (dm->GetArrayDim()>0 ){
+       printf("A:\t%s%s\t",prefix.Data(),rd->GetName());
+       Int_t counter=0;
+       for  (Int_t idim=0; idim<dm->GetArrayDim(); idim++){
+         //printf("A:%d\t%d\n", dm->GetArrayDim(),dm->GetMaxIndex(idim));
+         for (Int_t j=0; j<dm->GetMaxIndex(idim); j++){
+           printf("%s\t",dtype->AsString(pointer+dm->GetUnitSize()*counter));
+           counter++;
+           if (counter%5==0) printf("\nA:\t%s%s\t",prefix.Data(),rd->GetName());
+         }
+       }
+       printf("\n");
+      }
+      //
+      // Basic array - dynamic length
+      //
+      if (dm->GetArrayDim()>0 && strlen(index) != 0){
+       //
+       // Dump first only for the moment
+       //  
+       printf("B:\t%s%s\t%s\n",prefix.Data(),rd->GetName(), dtype->AsString(pointer));
+      }
+    } else {
+    }
+  }
+}  
+
+//
+// Small checks to test the TRealData and TDataType
+//
+
+
+
+void DumpDataSimple(){
+  //
+  // Dump example for elenatr data types 
+  //
+  TObject *obj = new TVectorD(20);
+  TClass * cl = obj->IsA();
+  if (!cl->GetListOfRealData()) cl->BuildRealData();
+  //
+  TRealData* rd = 0;
+  rd = (TRealData*)(cl->GetListOfRealData()->FindObject("fNrows"));
+  TDataMember* dm = rd->GetDataMember();
+  TDataType* dtype = dm->GetDataType();
+  //
+  Int_t offset = rd->GetThisOffset();
+  char* pointer = ((char*) obj) + offset;
+  printf("%s\n",dtype->AsString(pointer));
+}
+
+void DumpDataArray(){
+  //
+  // print array example
+  // 
+  TObject *obj = new TVectorD(20);
+  TClass * cl = obj->IsA();
+  if (!cl->GetListOfRealData()) cl->BuildRealData();
+  TRealData* rd = 0;
+  rd = (TRealData*)(cl->GetListOfRealData()->FindObject("*fElements"));
+  TDataMember* dm = rd->GetDataMember();
+  TDataType* dtype = dm->GetDataType();
+  dtype->Print();
+  //
+  Int_t offset = rd->GetThisOffset();
+  char* pointer = ((char*) obj) + offset; 
+  printf("%s\n",dtype->AsString(pointer));
+}
+
+void DumpTObjectArray(){
+  //
+  //
+  //
+  TObjArray *array = new TObjArray(10);
+  for (Int_t i=0; i<10; i++) array->AddLast(new TNamed(Form("n%d",i), Form("n%d",i)));  
+   AliOCDBtoolkit::DumpObjectRecursive(array);
+  //
+  //
+  TObject *obj = array;
+  TClass * cl = obj->IsA();
+  if (!cl->GetListOfRealData()) cl->BuildRealData();
+  TRealData* rd = 0;
+  rd = (TRealData*)(cl->GetListOfRealData()->FindObject("*fCont"));
+  TDataMember* dm = rd->GetDataMember();
+  TDataType* dtype = dm->GetDataType();
+  //
+  Int_t offset = rd->GetThisOffset();
+  char* pointer = ((char*) obj) + offset;
+  char** apointer = (char**) pointer;
+  //we have pointer to pointer here
+  TObject** ppobj = (TObject**) *apointer;
+  (*ppobj)->Print();
+  //
+  TIter myiter(array);
+  TObject  *arObject; 
+  dtype->Print();
+  while ((arObject = (TObject*)myiter.Next())) {
+    AliOCDBtoolkit::DumpObjectRecursive(arObject);
+  } 
+}
index d78ced320f8cbeb5e73af532535f652032e19668..b64676defac73cc187ad6e64e53dbe1881e563c5 100644 (file)
@@ -28,7 +28,13 @@ public:
   static void LoadOCDBFromMap(const TMap *cdbMap, const TList *cdbList);
   static void MakeDiff(const TMap *cdbMap0, const TList *cdbList0, const TMap *cdbMap1, const TList *cdbList1, Int_t verbose);
   static void DumpOCDB(const TMap *cdbMap0, const TList *cdbList0);
-  
+  //
+  // dump object functionality
+  //  
+  static void DumpObjectRecursive(TObject *obj);
+  static void DumpObjectRecursive(TObject *obj, TString prefix, Int_t &counterRec);
+  static void DumpOCDBFile(const char *finput , const char *foutput, Bool_t dumpMetaData, Bool_t xml);
+  //
 private:
   AliOCDBtoolkit(const AliOCDBtoolkit& source);
   AliOCDBtoolkit& operator= (const AliOCDBtoolkit& rec);